45#include <TObjString.h>
140 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PRunSingleHistoRRF(): **SEVERE ERROR**: no GLOBAL-block present!";
141 std::cerr << std::endl <<
">> For Single Histo RRF the GLOBAL-block is mandatory! Please fix this first.";
142 std::cerr << std::endl;
148 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PRunSingleHistoRRF(): **SEVERE ERROR**: no RRF-Frequency found!";
149 std::cerr << std::endl;
156 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PRunSingleHistoRRF(): **SEVERE ERROR**: no RRF-Packing found!";
157 std::cerr << std::endl;
170 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PRunSingleHistoRRF(): **SEVERE ERROR**: Couldn't prepare data for fitting!";
171 std::cerr << std::endl <<
">> This is very bad :-(, will quit ...";
172 std::cerr << std::endl;
241 Double_t chisq = 0.0;
245 for (Int_t i=0; i<
fMsrInfo->GetNoOfFuncs(); i++) {
246 UInt_t funcNo =
fMsrInfo->GetFuncNo(i);
264 #pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
267 time =
fData.GetDataTimeStart() +
static_cast<Double_t
>(i)*
fData.GetDataTimeStep();
269 chisq += diff*diff / (
fData.GetError()->at(i)*
fData.GetError()->at(i));
308 Double_t chisq = 0.0;
313 for (Int_t i=0; i<
fMsrInfo->GetNoOfFuncs(); i++) {
314 UInt_t funcNo =
fMsrInfo->GetFuncNo(i);
332 #pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
335 time =
fData.GetDataTimeStart() +
static_cast<Double_t
>(i)*
fData.GetDataTimeStep();
337 diff =
fData.GetValue()->at(i) - theo;
338 chisq += diff*diff / theo;
426 std::vector<Double_t> par;
428 for (UInt_t i=0; i<paramList->size(); i++)
429 par.push_back((*paramList)[i].fValue);
432 for (Int_t i=0; i<
fMsrInfo->GetNoOfFuncs(); i++) {
437 UInt_t size =
fData.GetValue()->size();
438 Double_t start =
fData.GetDataTimeStart();
439 Double_t resolution =
fData.GetDataTimeStep();
441 for (UInt_t i=0; i<size; i++) {
442 time = start +
static_cast<Double_t
>(i)*resolution;
531 TObjArray *tok =
nullptr;
532 TObjString *ostr =
nullptr;
537 tok = fitRange.Tokenize(
" \t");
539 if (tok->GetEntries() == 3) {
541 ostr =
dynamic_cast<TObjString*
>(tok->At(1));
542 str = ostr->GetString();
544 idx = str.First(
"+");
546 str.Remove(0, idx+1);
553 ostr =
dynamic_cast<TObjString*
>(tok->At(2));
554 str = ostr->GetString();
556 idx = str.First(
"-");
558 str.Remove(0, idx+1);
563 }
else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) {
564 Int_t pos = 2*(
fRunNo+1)-1;
566 if (pos + 1 >= tok->GetEntries()) {
567 std::cerr << std::endl <<
">> PRunSingleHistoRRF::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange <<
"'";
568 std::cerr << std::endl <<
">> will ignore it. Sorry ..." << std::endl;
571 ostr =
dynamic_cast<TObjString*
>(tok->At(pos));
572 str = ostr->GetString();
574 idx = str.First(
"+");
576 str.Remove(0, idx+1);
583 ostr =
dynamic_cast<TObjString*
>(tok->At(pos+1));
584 str = ostr->GetString();
586 idx = str.First(
"-");
588 str.Remove(0, idx+1);
595 std::cerr << std::endl <<
">> PRunSingleHistoRRF::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange <<
"'";
596 std::cerr << std::endl <<
">> will ignore it. Sorry ..." << std::endl;
708 Bool_t success =
true;
719 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PrepareData(): **ERROR** Couldn't get run " <<
fRunInfo->GetRunName()->Data() <<
"!";
720 std::cerr << std::endl;
736 for (UInt_t i=0; i<
fRunInfo->GetForwardHistoNoSize(); i++) {
737 histoNo.push_back(
fRunInfo->GetForwardHistoNo(i));
740 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PrepareData(): **PANIC ERROR**:";
741 std::cerr << std::endl <<
">> histoNo found = " << histoNo[i] <<
", which is NOT present in the data file!?!?";
742 std::cerr << std::endl <<
">> Will quit :-(";
743 std::cerr << std::endl;
751 std::cout.precision(10);
752 std::cout << std::endl <<
">> PRunSingleHisto::PrepareData(): time resolution=" << std::fixed << runData->
GetTimeResolution() <<
"(ns)" << std::endl;
760 std::vector<PDoubleVector> forward;
761 forward.resize(histoNo.size());
762 for (UInt_t i=0; i<histoNo.size(); i++) {
763 forward[i].resize(runData->
GetDataBin(histoNo[i])->size());
764 forward[i] = *runData->
GetDataBin(histoNo[i]);
768 if (
fRunInfo->GetRunNameSize() > 1) {
770 std::vector<PDoubleVector> addForward;
771 for (UInt_t i=1; i<
fRunInfo->GetRunNameSize(); i++) {
775 if (addRunData ==
nullptr) {
776 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PrepareData(): **ERROR** Couldn't get addrun " <<
fRunInfo->GetRunName(i)->Data() <<
"!";
777 std::cerr << std::endl;
782 addForward.resize(histoNo.size());
783 for (UInt_t j=0; j<histoNo.size(); j++) {
784 addForward[j].resize(addRunData->
GetDataBin(histoNo[j])->size());
785 addForward[j] = *addRunData->
GetDataBin(histoNo[j]);
791 for (UInt_t k=0; k<histoNo.size(); k++) {
792 addRunSize = addForward[k].size();
793 for (UInt_t j=0; j<addRunSize; j++) {
795 if ((
static_cast<Int_t
>(j)+
static_cast<Int_t
>(
fAddT0s[i-1][k])-
static_cast<Int_t
>(
fT0s[k]) >= 0) &&
796 (j+
static_cast<Int_t
>(
fAddT0s[i-1][k])-
static_cast<Int_t
>(
fT0s[k]) < addRunSize)) {
797 forward[k][j] += addForward[k][j+
static_cast<Int_t
>(
fAddT0s[i-1][k])-
static_cast<Int_t
>(
fT0s[k])];
806 for (UInt_t i=0; i<
fForward.size(); i++) {
811 for (UInt_t i=1; i<histoNo.size(); i++) {
812 for (UInt_t j=0; j<runData->
GetDataBin(histoNo[i])->size(); j++) {
814 if ((
static_cast<Int_t
>(j)+
static_cast<Int_t
>(
fT0s[i])-
static_cast<Int_t
>(
fT0s[0]) >= 0) &&
815 (j+
static_cast<Int_t
>(
fT0s[i])-
static_cast<Int_t
>(
fT0s[0]) < runData->
GetDataBin(histoNo[i])->size())) {
816 fForward[j] += forward[i][j+
static_cast<Int_t
>(
fT0s[i])-
static_cast<Int_t
>(
fT0s[0])];
919 for (UInt_t i=0; i<
fForward.size(); i++) {
923 std::cout <<
"info> freqMax=" << freqMax <<
" (MHz)" << std::endl;
926 Double_t optNoPoints = 8;
929 std::cout <<
"info> optimal packing: " <<
static_cast<Int_t
>(1.0 / (
fTimeResolution*(freqMax -
fMsrInfo->GetMsrGlobal()->GetRRFFreq(
"MHz"))) / optNoPoints);
936 if (
fRunInfo->GetBkgRange(0) >= 0) {
940 fRunInfo->SetBkgRange(
static_cast<Int_t
>(
fT0s[0]*0.1), 0);
941 fRunInfo->SetBkgRange(
static_cast<Int_t
>(
fT0s[0]*0.6), 1);
942 std::cerr << std::endl <<
">> PRunSingleHistoRRF::PrepareFitData(): **WARNING** Neither fix background nor background bins are given!";
943 std::cerr << std::endl <<
">> Will try the following: bkg start = " <<
fRunInfo->GetBkgRange(0) <<
", bkg end = " <<
fRunInfo->GetBkgRange(1);
944 std::cerr << std::endl <<
">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
945 std::cerr << std::endl;
950 for (UInt_t i=0; i<
fForward.size(); i++)
953 for (UInt_t i=0; i<
fForward.size(); i++) {
960 Int_t t0 =
static_cast<Int_t
>(
fT0s[0]);
965 Double_t time_tau=0.0;
966 Double_t exp_t_tau=0.0;
969 exp_t_tau = exp(time_tau);
976 for (UInt_t i=0; i<
fMerr.size(); i++) {
985 Double_t errN0 = 0.0;
996 exp_t_tau = exp(time_tau);
997 fAerr.push_back(exp_t_tau/n0*sqrt(rawNt[i]+pow(((rawNt[i]-
fBackground)/n0)*errN0,2.0)));
1002 Double_t wRRF = globalBlock->
GetRRFFreq(
"Mc");
1003 Double_t phaseRRF = globalBlock->
GetRRFPhase()*TMath::TwoPi()/180.0;
1004 Double_t time = 0.0;
1007 fForward[i] *= 2.0*cos(wRRF * time + phaseRRF);
1018 fData.AppendValue(dval);
1111 Int_t viewPacking =
fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1112 if (viewPacking > 0) {
1114 std::cerr <<
">> PRunSingleHistoRRF::PrepareViewData(): **WARNING** Found View Packing (" << viewPacking <<
") < RRF Packing (" <<
fRRFPacking <<
").";
1115 std::cerr <<
">> Will ignore View Packing." << std::endl;
1126 std::vector<Double_t> par;
1128 for (UInt_t i=0; i<paramList->size(); i++)
1129 par.push_back((*paramList)[i].fValue);
1132 for (Int_t i=0; i<
fMsrInfo->GetNoOfFuncs(); i++) {
1139 fData.SetTheoryTimeStart(
fData.GetDataTimeStart());
1141 fData.SetTheoryTimeStep(
fData.GetDataTimeStep());
1145 fData.SetTheoryTimeStep(
fData.GetDataTimeStep()/(Double_t)factor);
1149 Double_t time = 0.0;
1150 Double_t theoryValue = 0.0;
1151 for (UInt_t i=0; i<size; i++) {
1152 time =
fData.GetTheoryTimeStart() +
static_cast<Double_t
>(i)*
fData.GetTheoryTimeStep();
1154 if (fabs(theoryValue) > 10.0) {
1157 fData.AppendTheoryValue(theoryValue);
1219 fT0s.resize(histoNo.size());
1220 for (UInt_t i=0; i<
fT0s.size(); i++) {
1225 for (UInt_t i=0; i<
fRunInfo->GetT0BinSize(); i++) {
1231 if (
fT0s[i] == -1.0) {
1237 for (UInt_t i=0; i<histoNo.size(); i++) {
1238 if (
fT0s[i] == -1.0) {
1239 if (runData->
GetT0Bin(histoNo[i]) > 0.0) {
1247 for (UInt_t i=0; i<histoNo.size(); i++) {
1248 if (
fT0s[i] == -1.0) {
1252 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1253 std::cerr << std::endl <<
">> run: " <<
fRunInfo->GetRunName()->Data();
1254 std::cerr << std::endl <<
">> will try the estimated one: forward t0 = " << runData->
GetT0BinEstimated(histoNo[i]);
1255 std::cerr << std::endl <<
">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1256 std::cerr << std::endl;
1261 for (UInt_t i=0; i<
fRunInfo->GetForwardHistoNoSize(); i++) {
1262 if ((
fT0s[i] < 0.0) || (
fT0s[i] >
static_cast<Int_t
>(runData->
GetDataBin(histoNo[i])->size()))) {
1263 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperT0(): **ERROR** t0 data bin (" <<
fT0s[i] <<
") doesn't make any sense!";
1264 std::cerr << std::endl;
1270 if (
fRunInfo->GetRunNameSize() > 1) {
1273 for (UInt_t i=1; i<
fRunInfo->GetRunNameSize(); i++) {
1277 if (addRunData ==
nullptr) {
1278 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperT0(): **ERROR** Couldn't get addrun " <<
fRunInfo->GetRunName(i)->Data() <<
"!";
1279 std::cerr << std::endl;
1285 fAddT0s[i-1].resize(histoNo.size());
1286 for (UInt_t j=0; j<
fAddT0s[i-1].size(); j++) {
1291 for (UInt_t j=0; j<
fRunInfo->GetT0BinSize(); j++) {
1296 for (UInt_t j=0; j<histoNo.size(); j++) {
1298 if (addRunData->
GetT0Bin(histoNo[j]) > 0.0) {
1305 for (UInt_t j=0; j<histoNo.size(); j++) {
1306 if (
fAddT0s[i-1][j] == -1.0) {
1310 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1311 std::cerr << std::endl <<
">> run: " <<
fRunInfo->GetRunName(i)->Data();
1312 std::cerr << std::endl <<
">> will try the estimated one: forward t0 = " << addRunData->
GetT0BinEstimated(histoNo[j]);
1313 std::cerr << std::endl <<
">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1314 std::cerr << std::endl;
1319 for (UInt_t j=0; j<
fRunInfo->GetForwardHistoNoSize(); j++) {
1321 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperT0(): **ERROR** addt0 data bin (" <<
fAddT0s[i-1][j] <<
") doesn't make any sense!";
1322 std::cerr << std::endl;
1381 start =
fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1384 end =
fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1390 start =
static_cast<Int_t
>(
fT0s[0])+offset;
1392 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperDataRange(): **WARNING** data range was not provided, will try data range start = t0+" << offset <<
"(=10ns) = " << start <<
".";
1393 std::cerr << std::endl <<
">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1394 std::cerr << std::endl;
1399 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperDataRange(): **WARNING** data range was not provided, will try data range end = " << end <<
".";
1400 std::cerr << std::endl <<
">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1401 std::cerr << std::endl;
1412 if ((start < 0) || (start >
static_cast<Int_t
>(
fForward.size()))) {
1413 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperDataRange(): **ERROR** start data bin (" << start <<
") doesn't make any sense!";
1414 std::cerr << std::endl;
1419 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperDataRange(): **ERROR** end data bin (" << end <<
") doesn't make any sense!";
1420 std::cerr << std::endl;
1423 if (end >
static_cast<Int_t
>(
fForward.size())) {
1424 std::cerr << std::endl <<
">> PRunSingleHistoRRF::GetProperDataRange(): **WARNING** end data bin (" << end <<
") > histo length (" <<
fForward.size() <<
").";
1425 std::cerr << std::endl <<
">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
1426 std::cerr << std::endl;
1427 end =
static_cast<Int_t
>(
fForward.size()-1);
1438 std::cerr << std::endl <<
">> PRunSingleHisto::GetProperDataRange **WARNING** needed to shift forward lgb,";
1439 std::cerr << std::endl <<
">> from " <<
fGoodBins[1] <<
" to " <<
fForward.size()-1 << std::endl;
1518 std::cerr <<
">> PRunSingleHistoRRF::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
1519 std::cerr <<
">> Will set it to fgb/lgb which given in time is: " <<
fFitStartTime <<
"..." <<
fFitEndTime <<
" (usec)" << std::endl;
1565 Double_t freqMax = 0.0;
1572 histo->SetBinContent(i-
fGoodBins[0]+1, data[i]);
1576 std::unique_ptr<PFourier> ft = std::make_unique<PFourier>(histo.get(),
FOURIER_UNIT_FREQ);
1580 TH1F *power = ft->GetPowerFourier();
1581 Double_t maxFreqVal = 0.0;
1582 for (Int_t i=1; i<power->GetNbinsX(); i++) {
1584 if (i<power->GetNbinsX()-1) {
1585 if (power->GetBinContent(i)>power->GetBinContent(i+1))
1589 if (power->GetBinCenter(i) < 10.0)
1592 if (power->GetBinContent(i) > maxFreqVal) {
1593 maxFreqVal = power->GetBinContent(i);
1594 freqMax = power->GetBinCenter(i);
1662 for (Int_t i=0; i<endBin; i++) {
1671 for (Int_t i=0; i<endBin; i++) {
1674 errN0 = sqrt(errN0)/wN;
1676 std::cout <<
"info> PRunSingleHistoRRF::EstimateN0(): N0=" << n0 <<
"(" << errN0 <<
")" << std::endl;
1729 Double_t beamPeriod = 0.0;
1732 if (
fRunInfo->GetInstitute()->Contains(
"psi"))
1734 else if (
fRunInfo->GetInstitute()->Contains(
"ral"))
1736 else if (
fRunInfo->GetInstitute()->Contains(
"triumf"))
1742 UInt_t start =
fRunInfo->GetBkgRange(0);
1743 UInt_t end =
fRunInfo->GetBkgRange(1);
1745 std::cout << std::endl <<
"PRunSingleHistoRRF::EstimatBkg(): end = " << end <<
" > start = " << start <<
"! Will swap them!";
1752 if (beamPeriod != 0.0) {
1754 UInt_t fullCycles =
static_cast<UInt_t
>(timeBkg/beamPeriod);
1756 end = start +
static_cast<UInt_t
>((fullCycles*beamPeriod)/
fTimeResolution);
1757 std::cout << std::endl <<
"PRunSingleHistoRRF::EstimatBkg(): Background " << start <<
", " << end;
1764 std::cerr << std::endl <<
">> PRunSingleHistoRRF::EstimatBkg(): **ERROR** background bin values out of bound!";
1765 std::cerr << std::endl <<
">> histo lengths = " <<
fForward.size();
1766 std::cerr << std::endl <<
">> background start = " << start;
1767 std::cerr << std::endl;
1773 std::cerr << std::endl <<
">> PRunSingleHistoRRF::EstimatBkg(): **ERROR** background bin values out of bound!";
1774 std::cerr << std::endl <<
">> histo lengths = " <<
fForward.size();
1775 std::cerr << std::endl <<
">> background end = " << end;
1776 std::cerr << std::endl;
1784 for (UInt_t i=start; i<end; i++)
1786 bkg /=
static_cast<Double_t
>(end - start + 1);
1791 for (UInt_t i=start; i<end; i++)
1793 fBkgErr = sqrt(bkg/(
static_cast<Double_t
>(end - start)));
1795 std::cout << std::endl <<
"info> fBackground=" <<
fBackground <<
"(" <<
fBkgErr <<
")" << std::endl;
#define F_APODIZATION_STRONG
Strong apodization (heavy roll-off for best frequency resolution)
std::vector< UInt_t > PUIntVector
#define ACCEL_PERIOD_TRIUMF
TRIUMF accelerator cycle: 43.37 ns.
@ kFit
Fitting mode - perform least-squares fit to data.
@ kView
Viewing mode - display data and theory without fitting.
#define FOURIER_UNIT_FREQ
Frequency in MHz.
std::vector< PMsrParamStructure > PMsrParamList
#define ACCEL_PERIOD_PSI
PSI (Paul Scherrer Institute) accelerator cycle: 19.75 ns.
std::vector< Double_t > PDoubleVector
#define ACCEL_PERIOD_RAL
RAL (Rutherford Appleton Lab) - pulsed beam.
virtual Bool_t IsPresent()
virtual void SetFitRange(Double_t dval, UInt_t idx)
virtual Double_t GetFitRange(UInt_t idx)
virtual UInt_t GetT0BinSize()
virtual Bool_t IsFitRangeInBin()
virtual Int_t GetRRFPacking()
virtual TString GetRRFUnit()
virtual Int_t GetFitRangeOffset(UInt_t idx)
virtual Double_t GetT0Bin(UInt_t idx=0)
virtual Double_t GetRRFPhase()
virtual Double_t GetRRFFreq(const char *unit)
MSR file parser and manager for the musrfit framework.
virtual PMsrGlobalBlock * GetMsrGlobal()
Returns pointer to GLOBAL block settings.
virtual const PDoubleVector * GetDataBin(const UInt_t histoNo)
virtual const Double_t GetT0Bin(const UInt_t histoNo)
virtual const Double_t GetTimeResolution()
virtual const Bool_t IsPresent(UInt_t histoNo)
virtual const UInt_t GetNoOfTemperatures()
virtual const Double_t GetField()
virtual const Double_t GetEnergy()
virtual const Double_t GetT0BinEstimated(const UInt_t histoNo)
virtual const PDoublePairVector * GetTemperature() const
Double_t fTimeResolution
Time resolution of raw histogram data in microseconds (μs), e.g., 0.01953125 μs for PSI GPS.
Bool_t fValid
Flag indicating if run object initialized successfully; false if any error occurred.
Double_t fFitEndTime
Fit range end time in microseconds (μs) relative to t0.
PDoubleVector fFuncValues
Cached values of user-defined functions from FUNCTIONS block, evaluated at current parameters.
PMsrHandler * fMsrInfo
Pointer to MSR file handler (owned externally, not deleted here)
virtual void DeadTimeCorrection(std::vector< PDoubleVector > &histos, PUIntVector &histoNo)
carry out dead time correction
PMetaData fMetaData
Experimental metadata extracted from data file header (magnetic field, temperature,...
std::unique_ptr< PTheory > fTheory
Theory function evaluator (smart pointer, automatically deleted)
std::vector< PDoubleVector > fAddT0s
Time-zero bin values for additional runs to be added to main run.
EPMusrHandleTag fHandleTag
Operation mode: kFit (fitting), kView (display only), kEmpty (uninitialized)
PRunData fData
Processed data container: background-corrected, packed, with theory values.
PRunDataHandler * fRawData
Pointer to raw data handler (owned externally, not deleted here)
PDoubleVector fT0s
Time-zero bin values for all histograms in this run (forward, backward, etc.)
PRunBase()
Default constructor.
Int_t fRunNo
Run number (0-based index in MSR file RUN blocks)
PMsrRunBlock * fRunInfo
Pointer to this run's RUN block settings within fMsrInfo.
Double_t fFitStartTime
Fit range start time in microseconds (μs) relative to t0.
Raw data file reader and format converter for μSR data.
PDoubleVector fForward
Forward detector histogram data (progressively transformed during preparation)
Int_t fEndTimeBin
Last bin index in fit range (exclusive: loop as i < fEndTimeBin)
virtual Double_t CalcChiSquare(const std::vector< Double_t > &par)
Calculates χ² between RRF-transformed data and theory.
PDoubleVector fAerr
Asymmetry errors before RRF packing. Used for packed error calculation.
virtual Double_t GetMainFrequency(PDoubleVector &data)
Finds the dominant precession frequency in raw data.
PRunSingleHistoRRF()
Default constructor creating an empty, invalid RRF single histogram run object.
virtual Double_t EstimateN0(Double_t &errN0, Double_t freqMax)
Estimates initial normalization N₀ from lifetime-corrected data.
virtual Bool_t PrepareFitData(PRawRunData *runData, const UInt_t histoNo)
Performs full RRF transformation for fitting.
Double_t fBackground
Estimated or fixed background level in counts/bin (before packing)
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
Calculates expected χ² using theory variance instead of data variance.
virtual ~PRunSingleHistoRRF()
Virtual destructor releasing allocated resources.
virtual Bool_t PrepareData()
Main data preparation orchestrator for RRF single histogram analysis.
Int_t fStartTimeBin
First bin index in fit range (inclusive, 0-based in RRF-packed data)
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
Determines fit time range from MSR file settings.
Bool_t fTheoAsData
Theory resolution mode: true = at data points only, false = 8× finer grid for smooth Fourier transfor...
virtual Bool_t PrepareViewData(PRawRunData *runData, const UInt_t histoNo)
Prepares RRF data for viewing/plotting.
PDoubleVector fM
Lifetime-corrected histogram: M(t) = [N(t) - B] × exp(+t/τ_μ). Used for N₀ estimation.
virtual Bool_t EstimateBkg(UInt_t histoNo)
Estimates background from pre-t0 bins.
Int_t fGoodBins[2]
Good bin range: [0] = first good bin (fgb), [1] = last good bin (lgb). Used for COMMANDS block fit ra...
Double_t fN0EstimateEndTime
End time (μs) for N₀ estimation window. Rounded to integer number of oscillation cycles based on main...
Double_t fBkgErr
Statistical error on background estimate (std dev of background region)
PDoubleVector fMerr
Error on M(t): σ_M = exp(+t/τ_μ) × √(N(t) + σ_B²). Includes background error.
Int_t fRRFPacking
RRF packing factor from GLOBAL block (number of raw bins averaged together)
UInt_t fNoOfFitBins
Number of RRF-packed bins within fit range [fStartTimeBin, fEndTimeBin)
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &histoNo)
Determines and validates t0 values for all histograms.
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
Calculates maximum likelihood (not yet implemented for RRF).
virtual Bool_t GetProperDataRange()
Determines valid data range (first/last good bins).
PDoubleVector fW
Weights for N₀ estimation: W(t) = 1/σ_M². Used in weighted average.
virtual void SetFitRangeBin(const TString fitRange)
Sets fit range using bin-offset syntax from COMMANDS block.
virtual UInt_t GetNoOfFitBins()
Returns the number of bins included in the current fit range.
virtual void CalcTheory()
Evaluates theory function at all data points for viewing/plotting.
virtual void CalcNoOfFitBins()
Calculates start/end bin indices from fit time range.