musrfit 1.10.0
PRunSingleHistoRRF.cpp
Go to the documentation of this file.
1/***************************************************************************
2
3 PRunSingleHistoRRF.cpp
4
5 Author: Andreas Suter
6 e-mail: andreas.suter@psi.ch
7
8***************************************************************************/
9
10/***************************************************************************
11 * Copyright (C) 2007-2026 by Andreas Suter *
12 * andreas.suter@psi.ch *
13 * *
14 * This program is free software; you can redistribute it and/or modify *
15 * it under the terms of the GNU General Public License as published by *
16 * the Free Software Foundation; either version 2 of the License, or *
17 * (at your option) any later version. *
18 * *
19 * This program is distributed in the hope that it will be useful, *
20 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
21 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
22 * GNU General Public License for more details. *
23 * *
24 * You should have received a copy of the GNU General Public License *
25 * along with this program; if not, write to the *
26 * Free Software Foundation, Inc., *
27 * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
28 ***************************************************************************/
29
30#ifdef HAVE_CONFIG_H
31#include "config.h"
32#endif
33
34#ifdef HAVE_GOMP
35#include <omp.h>
36#endif
37
38#include <cmath>
39#include <iostream>
40#include <fstream>
41#include <memory>
42
43#include <TString.h>
44#include <TObjArray.h>
45#include <TObjString.h>
46#include <TH1F.h>
47
48#include "PMusr.h"
49#include "PFourier.h"
50#include "PRunSingleHistoRRF.h"
51
52//--------------------------------------------------------------------------
53// Constructor
54//--------------------------------------------------------------------------
74{
75 fNoOfFitBins = 0;
76 fBackground = 0.0;
77 fBkgErr = 1.0;
78 fRRFPacking = -1;
79 fTheoAsData = false;
80
81 // the 2 following variables are need in case fit range is given in bins, and since
82 // the fit range can be changed in the command block, these variables need to be accessible
83 fGoodBins[0] = -1;
84 fGoodBins[1] = -1;
85
86 fN0EstimateEndTime = 1.0; // end time in (us) over which N0 is estimated.
87}
88
89//--------------------------------------------------------------------------
90// Constructor
91//--------------------------------------------------------------------------
132PRunSingleHistoRRF::PRunSingleHistoRRF(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
133 PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
134{
135 fNoOfFitBins = 0;
136
137 PMsrGlobalBlock *global = msrInfo->GetMsrGlobal();
138
139 if (!global->IsPresent()) {
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;
143 fValid = false;
144 return;
145 }
146
147 if (!global->GetRRFUnit().CompareTo("??")) {
148 std::cerr << std::endl << ">> PRunSingleHistoRRF::PRunSingleHistoRRF(): **SEVERE ERROR**: no RRF-Frequency found!";
149 std::cerr << std::endl;
150 fValid = false;
151 return;
152 }
153
154 fRRFPacking = global->GetRRFPacking();
155 if (fRRFPacking == -1) {
156 std::cerr << std::endl << ">> PRunSingleHistoRRF::PRunSingleHistoRRF(): **SEVERE ERROR**: no RRF-Packing found!";
157 std::cerr << std::endl;
158 fValid = false;
159 return;
160 }
161
162 // the 2 following variables are need in case fit range is given in bins, and since
163 // the fit range can be changed in the command block, these variables need to be accessible
164 fGoodBins[0] = -1;
165 fGoodBins[1] = -1;
166
167 fN0EstimateEndTime = 1.0; // end time in (us) over which N0 is estimated.
168
169 if (!PrepareData()) {
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;
173 fValid = false;
174 }
175}
176
177//--------------------------------------------------------------------------
178// Destructor
179//--------------------------------------------------------------------------
197
198//--------------------------------------------------------------------------
199// CalcChiSquare (public)
200//--------------------------------------------------------------------------
239Double_t PRunSingleHistoRRF::CalcChiSquare(const std::vector<Double_t>& par)
240{
241 Double_t chisq = 0.0;
242 Double_t diff = 0.0;
243
244 // calculate functions
245 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
246 UInt_t funcNo = fMsrInfo->GetFuncNo(i);
247 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
248 }
249
250 // calculate chi square
251 Double_t time(1.0);
252 Int_t i;
253
254 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
255 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
256 // for a given set of parameters---which should be done outside of the parallelized loop.
257 // For all other functions it means a tiny and acceptable overhead.
258 time = fTheory->Func(time, par, fFuncValues);
259
260 #ifdef HAVE_GOMP
261 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
262 if (chunk < 10)
263 chunk = 10;
264 #pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
265 #endif
266 for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
267 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
268 diff = fData.GetValue()->at(i) - fTheory->Func(time, par, fFuncValues);
269 chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
270 }
271
272 return chisq;
273}
274
275//--------------------------------------------------------------------------
276// CalcChiSquareExpected (public)
277//--------------------------------------------------------------------------
306Double_t PRunSingleHistoRRF::CalcChiSquareExpected(const std::vector<Double_t>& par)
307{
308 Double_t chisq = 0.0;
309 Double_t diff = 0.0;
310 Double_t theo = 0.0;
311
312 // calculate functions
313 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
314 UInt_t funcNo = fMsrInfo->GetFuncNo(i);
315 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
316 }
317
318 // calculate chi square
319 Double_t time(1.0);
320 Int_t i;
321
322 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
323 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
324 // for a given set of parameters---which should be done outside of the parallelized loop.
325 // For all other functions it means a tiny and acceptable overhead.
326 time = fTheory->Func(time, par, fFuncValues);
327
328 #ifdef HAVE_GOMP
329 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
330 if (chunk < 10)
331 chunk = 10;
332 #pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
333 #endif
334 for (i=fStartTimeBin; i < fEndTimeBin; ++i) {
335 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
336 theo = fTheory->Func(time, par, fFuncValues);
337 diff = fData.GetValue()->at(i) - theo;
338 chisq += diff*diff / theo;
339 }
340
341 return chisq;
342}
343
344//--------------------------------------------------------------------------
345// CalcMaxLikelihood (public)
346//--------------------------------------------------------------------------
378Double_t PRunSingleHistoRRF::CalcMaxLikelihood(const std::vector<Double_t>& par)
379{
380 // not yet implemented
381
382 return 0.0;
383}
384
385//--------------------------------------------------------------------------
386// CalcTheory (public)
387//--------------------------------------------------------------------------
424{
425 // feed the parameter vector
426 std::vector<Double_t> par;
427 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
428 for (UInt_t i=0; i<paramList->size(); i++)
429 par.push_back((*paramList)[i].fValue);
430
431 // calculate functions
432 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
433 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
434 }
435
436 // calculate theory
437 UInt_t size = fData.GetValue()->size();
438 Double_t start = fData.GetDataTimeStart();
439 Double_t resolution = fData.GetDataTimeStep();
440 Double_t time;
441 for (UInt_t i=0; i<size; i++) {
442 time = start + static_cast<Double_t>(i)*resolution;
443 fData.AppendTheoryValue(fTheory->Func(time, par, fFuncValues));
444 }
445
446 // clean up
447 par.clear();
448}
449
450//--------------------------------------------------------------------------
451// GetNoOfFitBins (public)
452//--------------------------------------------------------------------------
473{
475
476 return fNoOfFitBins;
477}
478
479//--------------------------------------------------------------------------
480// SetFitRangeBin (public)
481//--------------------------------------------------------------------------
529void PRunSingleHistoRRF::SetFitRangeBin(const TString fitRange)
530{
531 TObjArray *tok = nullptr;
532 TObjString *ostr = nullptr;
533 TString str;
534 Ssiz_t idx = -1;
535 Int_t offset = 0;
536
537 tok = fitRange.Tokenize(" \t");
538
539 if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
540 // handle fgb+n0 entry
541 ostr = dynamic_cast<TObjString*>(tok->At(1));
542 str = ostr->GetString();
543 // check if there is an offset present
544 idx = str.First("+");
545 if (idx != -1) { // offset present
546 str.Remove(0, idx+1);
547 if (str.IsFloat()) // if str is a valid number, convert is to an integer
548 offset = str.Atoi();
549 }
550 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
551
552 // handle lgb-n1 entry
553 ostr = dynamic_cast<TObjString*>(tok->At(2));
554 str = ostr->GetString();
555 // check if there is an offset present
556 idx = str.First("-");
557 if (idx != -1) { // offset present
558 str.Remove(0, idx+1);
559 if (str.IsFloat()) // if str is a valid number, convert is to an integer
560 offset = str.Atoi();
561 }
562 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
563 } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
564 Int_t pos = 2*(fRunNo+1)-1;
565
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;
569 } else {
570 // handle fgb+n0 entry
571 ostr = dynamic_cast<TObjString*>(tok->At(pos));
572 str = ostr->GetString();
573 // check if there is an offset present
574 idx = str.First("+");
575 if (idx != -1) { // offset present
576 str.Remove(0, idx+1);
577 if (str.IsFloat()) // if str is a valid number, convert is to an integer
578 offset = str.Atoi();
579 }
580 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
581
582 // handle lgb-n1 entry
583 ostr = dynamic_cast<TObjString*>(tok->At(pos+1));
584 str = ostr->GetString();
585 // check if there is an offset present
586 idx = str.First("-");
587 if (idx != -1) { // offset present
588 str.Remove(0, idx+1);
589 if (str.IsFloat()) // if str is a valid number, convert is to an integer
590 offset = str.Atoi();
591 }
592 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
593 }
594 } else { // error
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;
597 }
598
599 // clean up
600 if (tok) {
601 delete tok;
602 }
603}
604
605//--------------------------------------------------------------------------
606// CalcNoOfFitBins (public)
607//--------------------------------------------------------------------------
642{
643 // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
644 fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
645 if (fStartTimeBin < 0)
646 fStartTimeBin = 0;
647 fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
648 if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
649 fEndTimeBin = fData.GetValue()->size();
650
653 else
654 fNoOfFitBins = 0;
655}
656
657//--------------------------------------------------------------------------
658// PrepareData (protected)
659//--------------------------------------------------------------------------
707{
708 Bool_t success = true;
709
710 if (!fValid)
711 return false;
712
713 // keep the Global block info
714 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
715
716 // get the proper run
717 PRawRunData* runData = fRawData->GetRunData(*fRunInfo->GetRunName());
718 if (!runData) { // couldn't get run
719 std::cerr << std::endl << ">> PRunSingleHistoRRF::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
720 std::cerr << std::endl;
721 return false;
722 }
723
724 // keep the field from the meta-data from the data-file
725 fMetaData.fField = runData->GetField();
726
727 // keep the energy from the meta-data from the data-file
728 fMetaData.fEnergy = runData->GetEnergy();
729
730 // keep the temperature(s) from the meta-data from the data-file
731 for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
732 fMetaData.fTemp.push_back(runData->GetTemperature(i));
733
734 // collect histogram numbers
735 PUIntVector histoNo; // histoNo = msr-file forward + redGreen_offset - 1
736 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
737 histoNo.push_back(fRunInfo->GetForwardHistoNo(i));
738
739 if (!runData->IsPresent(histoNo[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;
744 histoNo.clear();
745 return false;
746 }
747 }
748
749 // keep the time resolution in (us)
750 fTimeResolution = runData->GetTimeResolution()/1.0e3;
751 std::cout.precision(10);
752 std::cout << std::endl << ">> PRunSingleHisto::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
753
754 // get all the proper t0's and addt0's for the current RUN block
755 if (!GetProperT0(runData, globalBlock, histoNo)) {
756 return false;
757 }
758
759 // keep the histo of each group at this point (addruns handled below)
760 std::vector<PDoubleVector> forward;
761 forward.resize(histoNo.size()); // resize to number of groups
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]);
765 }
766
767 // check if there are runs to be added to the current one
768 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
769 PRawRunData *addRunData;
770 std::vector<PDoubleVector> addForward;
771 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
772
773 // get run to be added to the main one
774 addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
775 if (addRunData == nullptr) { // couldn't get run
776 std::cerr << std::endl << ">> PRunSingleHistoRRF::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
777 std::cerr << std::endl;
778 return false;
779 }
780
781 addForward.clear();
782 addForward.resize(histoNo.size()); // resize to number of groups
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]);
786 }
787 DeadTimeCorrection(addForward, histoNo);
788
789 // add forward run
790 UInt_t addRunSize;
791 for (UInt_t k=0; k<histoNo.size(); k++) { // fill each group
792 addRunSize = addForward[k].size();
793 for (UInt_t j=0; j<addRunSize; j++) { // loop over the bin indices
794 // make sure that the index stays in the proper range
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])];
798 }
799 }
800 }
801 }
802 }
803
804 // set forward histo data of the first group
805 fForward.resize(forward[0].size());
806 for (UInt_t i=0; i<fForward.size(); i++) {
807 fForward[i] = forward[0][i];
808 }
809
810 // group histograms, add all the remaining forward histograms of the group
811 for (UInt_t i=1; i<histoNo.size(); i++) { // loop over the groupings
812 for (UInt_t j=0; j<runData->GetDataBin(histoNo[i])->size(); j++) { // loop over the bin indices
813 // make sure that the index stays within proper range
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])];
817 }
818 }
819 }
820
821 // get the data range (fgb/lgb) for the current RUN block
822 if (!GetProperDataRange()) {
823 return false;
824 }
825
826 // get the fit range for the current RUN block
827 GetProperFitRange(globalBlock);
828
829 // do the more fit/view specific stuff
830 if (fHandleTag == kFit)
831 success = PrepareFitData(runData, histoNo[0]);
832 else if (fHandleTag == kView)
833 success = PrepareViewData(runData, histoNo[0]);
834 else
835 success = false;
836
837 // cleanup
838 histoNo.clear();
839
840 return success;
841}
842
843//--------------------------------------------------------------------------
844// PrepareFitData (protected)
845//--------------------------------------------------------------------------
915Bool_t PRunSingleHistoRRF::PrepareFitData(PRawRunData* runData, const UInt_t histoNo)
916{
917 // keep the raw data for the RRF asymmetry error estimate for later
918 PDoubleVector rawNt;
919 for (UInt_t i=0; i<fForward.size(); i++) {
920 rawNt.push_back(fForward[i]); // N(t) without any corrections
921 }
922 Double_t freqMax = GetMainFrequency(rawNt);
923 std::cout << "info> freqMax=" << freqMax << " (MHz)" << std::endl;
924
925 // "optimal packing"
926 Double_t optNoPoints = 8;
927 if (freqMax < 271.0) // < 271 MHz, i.e ~ 2T
928 optNoPoints = 5;
929 std::cout << "info> optimal packing: " << static_cast<Int_t>(1.0 / (fTimeResolution*(freqMax - fMsrInfo->GetMsrGlobal()->GetRRFFreq("MHz"))) / optNoPoints);
930
931 // initially fForward is the "raw data set" (i.e. grouped histo and raw runs already added) to be fitted. This means fForward = N(t) at this point.
932
933 // 1) check how the background shall be handled
934 // subtract background from histogramms ------------------------------------------
935 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
936 if (fRunInfo->GetBkgRange(0) >= 0) {
937 if (!EstimateBkg(histoNo))
938 return false;
939 } else { // no background given to do the job, try estimate
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;
946 if (!EstimateBkg(histoNo))
947 return false;
948 }
949 // subtract background from fForward
950 for (UInt_t i=0; i<fForward.size(); i++)
951 fForward[i] -= fBackground;
952 } else { // fixed background given
953 for (UInt_t i=0; i<fForward.size(); i++) {
954 fForward[i] -= fRunInfo->GetBkgFix(0);
955 }
956 fBackground = fRunInfo->GetBkgFix(0);
957 }
958 // here fForward = N(t) - Nbkg
959
960 Int_t t0 = static_cast<Int_t>(fT0s[0]);
961
962 // 2) N(t) - Nbkg -> exp(+t/tau) [N(t)-Nbkg]
963 Double_t startTime = fTimeResolution * (static_cast<Double_t>(fGoodBins[0]) - static_cast<Double_t>(t0));
964
965 Double_t time_tau=0.0;
966 Double_t exp_t_tau=0.0;
967 for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
968 time_tau = (startTime + fTimeResolution * (i - fGoodBins[0])) / PMUON_LIFETIME;
969 exp_t_tau = exp(time_tau);
970 fForward[i] *= exp_t_tau;
971 fM.push_back(fForward[i]); // i.e. M(t) = [N(t)-Nbkg] exp(+t/tau); needed to estimate N0 later on
972 fMerr.push_back(exp_t_tau*sqrt(rawNt[i]+fBkgErr*fBkgErr));
973 }
974
975 // calculate weights
976 for (UInt_t i=0; i<fMerr.size(); i++) {
977 if (fMerr[i] > 0.0)
978 fW.push_back(1.0/(fMerr[i]*fMerr[i]));
979 else
980 fW.push_back(1.0);
981 }
982 // now fForward = exp(+t/tau) [N(t)-Nbkg] = M(t)
983
984 // 3) estimate N0
985 Double_t errN0 = 0.0;
986 Double_t n0 = EstimateN0(errN0, freqMax);
987
988 // 4a) A(t) = exp(+t/tau) [N(t)-Nbkg] / N0 - 1.0
989 for (Int_t i=fGoodBins[0]; i<=fGoodBins[1]; i++) {
990 fForward[i] = fForward[i] / n0 - 1.0;
991 }
992
993 // 4b) error estimate of A(t): errA(t) = exp(+t/tau)/N0 sqrt( N(t) + ([N(t)-N_bkg]/N0)^2 errN0^2 )
994 for (Int_t i=fGoodBins[0]; i<=fGoodBins[1]; i++) {
995 time_tau = (startTime + fTimeResolution * (i - fGoodBins[0])) / PMUON_LIFETIME;
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)));
998 }
999
1000 // 5) rotate A(t): A(t) -> 2* A(t) * cos(wRRF t + phiRRF), the factor 2.0 is needed since the high frequency part is suppressed.
1001 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
1002 Double_t wRRF = globalBlock->GetRRFFreq("Mc");
1003 Double_t phaseRRF = globalBlock->GetRRFPhase()*TMath::TwoPi()/180.0;
1004 Double_t time = 0.0;
1005 for (Int_t i=fGoodBins[0]; i<=fGoodBins[1]; i++) {
1006 time = startTime + fTimeResolution * (static_cast<Double_t>(i) - static_cast<Double_t>(fGoodBins[0]));
1007 fForward[i] *= 2.0*cos(wRRF * time + phaseRRF);
1008 }
1009
1010 // 6) RRF packing
1011 Double_t dval=0.0;
1012 for (Int_t i=fGoodBins[0]; i<=fGoodBins[1]; i++) {
1013 if (fRRFPacking == 1) {
1014 fData.AppendValue(fForward[i]);
1015 } else { // RRF packing > 1
1016 if (((i-fGoodBins[0]) % fRRFPacking == 0) && (i != fGoodBins[0])) { // fill data
1017 dval /= fRRFPacking;
1018 fData.AppendValue(dval);
1019 // reset dval
1020 dval = 0.0;
1021 }
1022 dval += fForward[i];
1023 }
1024 }
1025
1026 // 7) estimate packed RRF errors (see log-book p.204)
1027 // the error estimate of the unpacked RRF asymmetry is: errA_RRF(t) \simeq exp(t/tau)/N0 sqrt( [N(t) + ((N(t)-N_bkg)/N0)^2 errN0^2] )
1028 dval = 0.0;
1029 // the packed RRF asymmetry error
1030 for (Int_t i=fGoodBins[0]; i<=fGoodBins[1]; i++) {
1031 if (((i-fGoodBins[0]) % fRRFPacking == 0) && (i != fGoodBins[0])) { // fill data
1032 fData.AppendErrorValue(sqrt(2.0*dval)/fRRFPacking); // the factor 2.0 is needed since the high frequency part is suppressed.
1033 dval = 0.0;
1034 }
1035 dval += fAerr[i-fGoodBins[0]]*fAerr[i-fGoodBins[0]];
1036 }
1037
1038 // set start time and time step
1039 fData.SetDataTimeStart(fTimeResolution*(static_cast<Double_t>(fGoodBins[0])-static_cast<Double_t>(t0)+static_cast<Double_t>(fRRFPacking-1)/2.0));
1040 fData.SetDataTimeStep(fTimeResolution*fRRFPacking);
1041
1043
1044 return true;
1045}
1046
1047//--------------------------------------------------------------------------
1048// PrepareViewData (protected)
1049//--------------------------------------------------------------------------
1101Bool_t PRunSingleHistoRRF::PrepareViewData(PRawRunData* runData, const UInt_t histoNo)
1102{
1103 // --------------
1104 // prepare data
1105 // --------------
1106
1107 // prepare RRF single histo
1108 PrepareFitData(runData, histoNo);
1109
1110 // check for view packing
1111 Int_t viewPacking = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1112 if (viewPacking > 0) {
1113 if (viewPacking < fRRFPacking) {
1114 std::cerr << ">> PRunSingleHistoRRF::PrepareViewData(): **WARNING** Found View Packing (" << viewPacking << ") < RRF Packing (" << fRRFPacking << ").";
1115 std::cerr << ">> Will ignore View Packing." << std::endl;
1116 } else {
1117 // STILL MISSING
1118 }
1119 }
1120
1121 // --------------
1122 // prepare theory
1123 // --------------
1124
1125 // feed the parameter vector
1126 std::vector<Double_t> par;
1127 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1128 for (UInt_t i=0; i<paramList->size(); i++)
1129 par.push_back((*paramList)[i].fValue);
1130
1131 // calculate functions
1132 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1133 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
1134 }
1135
1136 // check if a finer binning for the theory is needed
1137 UInt_t size = fForward.size();
1138 Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1139 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1140 if (fTheoAsData) { // calculate theory only at the data points
1141 fData.SetTheoryTimeStep(fData.GetDataTimeStep());
1142 } else {
1143 // finer binning for the theory (8 times as many points = factor)
1144 size *= factor;
1145 fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1146 }
1147
1148 // calculate theory
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();
1153 theoryValue = fTheory->Func(time, par, fFuncValues);
1154 if (fabs(theoryValue) > 10.0) { // dirty hack needs to be fixed!!
1155 theoryValue = 0.0;
1156 }
1157 fData.AppendTheoryValue(theoryValue);
1158 }
1159
1160 return true;
1161}
1162
1163//--------------------------------------------------------------------------
1164// GetProperT0 (private)
1165//--------------------------------------------------------------------------
1215{
1216 // feed all T0's
1217 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1218 fT0s.clear();
1219 fT0s.resize(histoNo.size());
1220 for (UInt_t i=0; i<fT0s.size(); i++) {
1221 fT0s[i] = -1.0;
1222 }
1223
1224 // fill in the T0's from the msr-file (if present)
1225 for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1226 fT0s[i] = fRunInfo->GetT0Bin(i);
1227 }
1228
1229 // fill in the T0's from the GLOBAL block section (if present)
1230 for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1231 if (fT0s[i] == -1.0) { // i.e. not given in the RUN block section
1232 fT0s[i] = globalBlock->GetT0Bin(i);
1233 }
1234 }
1235
1236 // fill in the T0's from the data file, if not already present in the msr-file
1237 for (UInt_t i=0; i<histoNo.size(); i++) {
1238 if (fT0s[i] == -1.0) { // i.e. not present in the msr-file, try the data file
1239 if (runData->GetT0Bin(histoNo[i]) > 0.0) {
1240 fT0s[i] = runData->GetT0Bin(histoNo[i]);
1241 fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
1242 }
1243 }
1244 }
1245
1246 // fill in the T0's gaps, i.e. in case the T0's are NOT in the msr-file and NOT in the data file
1247 for (UInt_t i=0; i<histoNo.size(); i++) {
1248 if (fT0s[i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1249 fT0s[i] = runData->GetT0BinEstimated(histoNo[i]);
1250 fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
1251
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;
1257 }
1258 }
1259
1260 // check if t0 is within proper bounds
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;
1265 return false;
1266 }
1267 }
1268
1269 // check if there are runs to be added to the current one. If yes keep the needed t0's
1270 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1271 PRawRunData *addRunData;
1272 fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1273 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1274
1275 // get run to be added to the main one
1276 addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
1277 if (addRunData == nullptr) { // couldn't get run
1278 std::cerr << std::endl << ">> PRunSingleHistoRRF::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1279 std::cerr << std::endl;
1280 return false;
1281 }
1282
1283 // feed all T0's
1284 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1285 fAddT0s[i-1].resize(histoNo.size());
1286 for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1287 fAddT0s[i-1][j] = -1.0;
1288 }
1289
1290 // fill in the T0's from the msr-file (if present)
1291 for (UInt_t j=0; j<fRunInfo->GetT0BinSize(); j++) {
1292 fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1,j); // addRunIdx starts at 0
1293 }
1294
1295 // fill in the T0's from the data file, if not already present in the msr-file
1296 for (UInt_t j=0; j<histoNo.size(); j++) {
1297 if (fAddT0s[i-1][j] == -1.0) // i.e. not present in the msr-file, try the data file
1298 if (addRunData->GetT0Bin(histoNo[j]) > 0.0) {
1299 fAddT0s[i-1][j] = addRunData->GetT0Bin(histoNo[j]);
1300 fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1301 }
1302 }
1303
1304 // fill in the T0's gaps, i.e. in case the T0's are NOT in the msr-file and NOT in the data file
1305 for (UInt_t j=0; j<histoNo.size(); j++) {
1306 if (fAddT0s[i-1][j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1307 fAddT0s[i-1][j] = addRunData->GetT0BinEstimated(histoNo[j]);
1308 fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1309
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;
1315 }
1316 }
1317
1318 // check if t0 is within proper bounds
1319 for (UInt_t j=0; j<fRunInfo->GetForwardHistoNoSize(); j++) {
1320 if ((fAddT0s[i-1][j] < 0) || (fAddT0s[i-1][j] > static_cast<Int_t>(addRunData->GetDataBin(histoNo[j])->size()))) {
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;
1323 return false;
1324 }
1325 }
1326 }
1327 }
1328
1329 return true;
1330}
1331
1332//--------------------------------------------------------------------------
1333// GetProperDataRange (private)
1334//--------------------------------------------------------------------------
1372{
1373 // get start/end data
1374 Int_t start;
1375 Int_t end;
1376 start = fRunInfo->GetDataRange(0);
1377 end = fRunInfo->GetDataRange(1);
1378
1379 // check if data range has been given in the RUN block, if not try to get it from the GLOBAL block
1380 if (start < 0) {
1381 start = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1382 }
1383 if (end < 0) {
1384 end = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1385 }
1386
1387 // check if data range has been provided, and if not try to estimate them
1388 if (start < 0) {
1389 Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
1390 start = static_cast<Int_t>(fT0s[0])+offset;
1391 fRunInfo->SetDataRange(start, 0);
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;
1395 }
1396 if (end < 0) {
1397 end = fForward.size();
1398 fRunInfo->SetDataRange(end, 1);
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;
1402 }
1403
1404 // check if start and end make any sense
1405 // 1st check if start and end are in proper order
1406 if (end < start) { // need to swap them
1407 Int_t keep = end;
1408 end = start;
1409 start = keep;
1410 }
1411 // 2nd check if start is within proper bounds
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;
1415 return false;
1416 }
1417 // 3rd check if end is within proper bounds
1418 if (end < 0) {
1419 std::cerr << std::endl << ">> PRunSingleHistoRRF::GetProperDataRange(): **ERROR** end data bin (" << end << ") doesn't make any sense!";
1420 std::cerr << std::endl;
1421 return false;
1422 }
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);
1428 }
1429
1430 // keep good bins for potential later use
1431 fGoodBins[0] = start;
1432 fGoodBins[1] = end;
1433
1434 // make sure that fGoodBins are in proper range for fForward
1435 if (fGoodBins[0] < 0)
1436 fGoodBins[0]=0;
1437 if (fGoodBins[1] > fForward.size()) {
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;
1440 fGoodBins[1]=fForward.size()-1;
1441 }
1442
1443 return true;
1444}
1445
1446//--------------------------------------------------------------------------
1447// GetProperFitRange (private)
1448//--------------------------------------------------------------------------
1491{
1492 // set fit start/end time; first check RUN Block
1493 fFitStartTime = fRunInfo->GetFitRange(0);
1494 fFitEndTime = fRunInfo->GetFitRange(1);
1495 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1496 if (fRunInfo->IsFitRangeInBin()) {
1497 fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1498 fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1499 // write these times back into the data structure. This way it is available when writting the log-file
1500 fRunInfo->SetFitRange(fFitStartTime, 0);
1501 fRunInfo->SetFitRange(fFitEndTime, 1);
1502 }
1503 if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
1504 fFitStartTime = globalBlock->GetFitRange(0);
1505 fFitEndTime = globalBlock->GetFitRange(1);
1506 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1507 if (globalBlock->IsFitRangeInBin()) {
1508 fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1509 fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1510 // write these times back into the data structure. This way it is available when writting the log-file
1511 globalBlock->SetFitRange(fFitStartTime, 0);
1512 globalBlock->SetFitRange(fFitEndTime, 1);
1513 }
1514 }
1516 fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
1517 fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
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;
1520 }
1521}
1522
1523//--------------------------------------------------------------------------
1524// GetMainFrequency (private)
1525//--------------------------------------------------------------------------
1564{
1565 Double_t freqMax = 0.0;
1566
1567 // create histo
1568 Double_t startTime = (fGoodBins[0]-fT0s[0]) * fTimeResolution;
1569 Int_t noOfBins = fGoodBins[1]-fGoodBins[0]+1;
1570 std::unique_ptr<TH1F> histo = std::make_unique<TH1F>("data", "data", noOfBins, startTime-fTimeResolution/2.0, startTime+fTimeResolution/2.0+noOfBins*fTimeResolution);
1571 for (Int_t i=fGoodBins[0]; i<=fGoodBins[1]; i++) {
1572 histo->SetBinContent(i-fGoodBins[0]+1, data[i]);
1573 }
1574
1575 // Fourier transform
1576 std::unique_ptr<PFourier> ft = std::make_unique<PFourier>(histo.get(), FOURIER_UNIT_FREQ);
1577 ft->Transform(F_APODIZATION_STRONG);
1578
1579 // find frequency maximum
1580 TH1F *power = ft->GetPowerFourier();
1581 Double_t maxFreqVal = 0.0;
1582 for (Int_t i=1; i<power->GetNbinsX(); i++) {
1583 // ignore dc part at 0 frequency
1584 if (i<power->GetNbinsX()-1) {
1585 if (power->GetBinContent(i)>power->GetBinContent(i+1))
1586 continue;
1587 }
1588 // ignore everything below 10 MHz
1589 if (power->GetBinCenter(i) < 10.0)
1590 continue;
1591 // check for maximum
1592 if (power->GetBinContent(i) > maxFreqVal) {
1593 maxFreqVal = power->GetBinContent(i);
1594 freqMax = power->GetBinCenter(i);
1595 }
1596 }
1597
1598 // clean up
1599 if (power)
1600 delete power;
1601
1602 return freqMax;
1603}
1604
1605//--------------------------------------------------------------------------
1606// EstimateN0 (private)
1607//--------------------------------------------------------------------------
1654Double_t PRunSingleHistoRRF::EstimateN0(Double_t &errN0, Double_t freqMax)
1655{
1656 // endBin is estimated such that the number of full cycles (according to the maximum frequency of the data)
1657 // is approximately the time fN0EstimateEndTime.
1658 Int_t endBin = static_cast<Int_t>(round(ceil(fN0EstimateEndTime*freqMax/TMath::TwoPi()) * (TMath::TwoPi()/freqMax) / fTimeResolution));
1659
1660 Double_t n0 = 0.0;
1661 Double_t wN = 0.0;
1662 for (Int_t i=0; i<endBin; i++) {
1663// n0 += fW[i]*fM[i];
1664 n0 += fM[i];
1665 wN += fW[i];
1666 }
1667// n0 /= wN;
1668 n0 /= endBin;
1669
1670 errN0 = 0.0;
1671 for (Int_t i=0; i<endBin; i++) {
1672 errN0 += fW[i]*fW[i]*fMerr[i]*fMerr[i];
1673 }
1674 errN0 = sqrt(errN0)/wN;
1675
1676 std::cout << "info> PRunSingleHistoRRF::EstimateN0(): N0=" << n0 << "(" << errN0 << ")" << std::endl;
1677
1678 return n0;
1679}
1680
1681//--------------------------------------------------------------------------
1682// EstimatBkg (private)
1683//--------------------------------------------------------------------------
1728{
1729 Double_t beamPeriod = 0.0;
1730
1731 // check if data are from PSI, RAL, or TRIUMF
1732 if (fRunInfo->GetInstitute()->Contains("psi"))
1733 beamPeriod = ACCEL_PERIOD_PSI;
1734 else if (fRunInfo->GetInstitute()->Contains("ral"))
1735 beamPeriod = ACCEL_PERIOD_RAL;
1736 else if (fRunInfo->GetInstitute()->Contains("triumf"))
1737 beamPeriod = ACCEL_PERIOD_TRIUMF;
1738 else
1739 beamPeriod = 0.0;
1740
1741 // check if start and end are in proper order
1742 UInt_t start = fRunInfo->GetBkgRange(0);
1743 UInt_t end = fRunInfo->GetBkgRange(1);
1744 if (end < start) {
1745 std::cout << std::endl << "PRunSingleHistoRRF::EstimatBkg(): end = " << end << " > start = " << start << "! Will swap them!";
1746 UInt_t keep = end;
1747 end = start;
1748 start = keep;
1749 }
1750
1751 // calculate proper background range
1752 if (beamPeriod != 0.0) {
1753 Double_t timeBkg = static_cast<Double_t>(end-start)*fTimeResolution; // length of the background intervall in time
1754 UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
1755 // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
1756 end = start + static_cast<UInt_t>((fullCycles*beamPeriod)/fTimeResolution);
1757 std::cout << std::endl << "PRunSingleHistoRRF::EstimatBkg(): Background " << start << ", " << end;
1758 if (end == start)
1759 end = fRunInfo->GetBkgRange(1);
1760 }
1761
1762 // check if start is within histogram bounds
1763 if (start >= fForward.size()) {
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;
1768 return false;
1769 }
1770
1771 // check if end is within histogram bounds
1772 if (end >= fForward.size()) {
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;
1777 return false;
1778 }
1779
1780 // calculate background
1781 Double_t bkg = 0.0;
1782
1783 // forward
1784 for (UInt_t i=start; i<end; i++)
1785 bkg += fForward[i];
1786 bkg /= static_cast<Double_t>(end - start + 1);
1787
1788 fBackground = bkg; // keep background (per bin)
1789
1790 bkg = 0.0;
1791 for (UInt_t i=start; i<end; i++)
1792 bkg += pow(fForward[i]-fBackground, 2.0);
1793 fBkgErr = sqrt(bkg/(static_cast<Double_t>(end - start)));
1794
1795 std::cout << std::endl << "info> fBackground=" << fBackground << "(" << fBkgErr << ")" << std::endl;
1796
1797 fRunInfo->SetBkgEstimated(fBackground, 0);
1798
1799 return true;
1800}
#define F_APODIZATION_STRONG
Strong apodization (heavy roll-off for best frequency resolution)
Definition PFourier.h:59
std::vector< UInt_t > PUIntVector
Definition PMusr.h:361
#define ACCEL_PERIOD_TRIUMF
TRIUMF accelerator cycle: 43.37 ns.
Definition PMusr.h:151
EPMusrHandleTag
Definition PMusr.h:413
@ kFit
Fitting mode - perform least-squares fit to data.
Definition PMusr.h:415
@ kView
Viewing mode - display data and theory without fitting.
Definition PMusr.h:416
#define FOURIER_UNIT_FREQ
Frequency in MHz.
Definition PMusr.h:276
#define PMUSR_UNDEFINED
Definition PMusr.h:172
#define PMUON_LIFETIME
Definition PMusr.h:119
std::vector< PMsrParamStructure > PMsrParamList
Definition PMusr.h:1022
#define ACCEL_PERIOD_PSI
PSI (Paul Scherrer Institute) accelerator cycle: 19.75 ns.
Definition PMusr.h:149
std::vector< Double_t > PDoubleVector
Definition PMusr.h:385
#define ACCEL_PERIOD_RAL
RAL (Rutherford Appleton Lab) - pulsed beam.
Definition PMusr.h:153
if(xmlFile.is_open())
virtual Bool_t IsPresent()
Definition PMusr.h:1042
virtual void SetFitRange(Double_t dval, UInt_t idx)
Definition PMusr.cpp:1137
virtual Double_t GetFitRange(UInt_t idx)
Definition PMusr.cpp:1120
virtual UInt_t GetT0BinSize()
Definition PMusr.h:1050
virtual Bool_t IsFitRangeInBin()
Definition PMusr.h:1055
virtual Int_t GetRRFPacking()
Definition PMusr.h:1047
virtual TString GetRRFUnit()
Definition PMusr.cpp:903
virtual Int_t GetFitRangeOffset(UInt_t idx)
Definition PMusr.cpp:1157
virtual Double_t GetT0Bin(UInt_t idx=0)
Definition PMusr.cpp:1004
virtual Double_t GetRRFPhase()
Definition PMusr.h:1046
virtual Double_t GetRRFFreq(const char *unit)
Definition PMusr.cpp:831
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)
Definition PMusr.h:880
virtual const Double_t GetT0Bin(const UInt_t histoNo)
Definition PMusr.h:870
virtual const Double_t GetTimeResolution()
Definition PMusr.h:868
virtual const Bool_t IsPresent(UInt_t histoNo)
Definition PMusr.h:869
virtual const UInt_t GetNoOfTemperatures()
Definition PMusr.h:860
virtual const Double_t GetField()
Definition PMusr.h:859
virtual const Double_t GetEnergy()
Definition PMusr.h:864
virtual const Double_t GetT0BinEstimated(const UInt_t histoNo)
Definition PMusr.h:871
virtual const PDoublePairVector * GetTemperature() const
Definition PMusr.h:861
Double_t fTimeResolution
Time resolution of raw histogram data in microseconds (μs), e.g., 0.01953125 μs for PSI GPS.
Definition PRunBase.h:276
Bool_t fValid
Flag indicating if run object initialized successfully; false if any error occurred.
Definition PRunBase.h:266
Double_t fFitEndTime
Fit range end time in microseconds (μs) relative to t0.
Definition PRunBase.h:282
PDoubleVector fFuncValues
Cached values of user-defined functions from FUNCTIONS block, evaluated at current parameters.
Definition PRunBase.h:284
PMsrHandler * fMsrInfo
Pointer to MSR file handler (owned externally, not deleted here)
Definition PRunBase.h:271
virtual void DeadTimeCorrection(std::vector< PDoubleVector > &histos, PUIntVector &histoNo)
carry out dead time correction
Definition PRunBase.cpp:171
PMetaData fMetaData
Experimental metadata extracted from data file header (magnetic field, temperature,...
Definition PRunBase.h:277
std::unique_ptr< PTheory > fTheory
Theory function evaluator (smart pointer, automatically deleted)
Definition PRunBase.h:285
std::vector< PDoubleVector > fAddT0s
Time-zero bin values for additional runs to be added to main run.
Definition PRunBase.h:279
EPMusrHandleTag fHandleTag
Operation mode: kFit (fitting), kView (display only), kEmpty (uninitialized)
Definition PRunBase.h:268
PRunData fData
Processed data container: background-corrected, packed, with theory values.
Definition PRunBase.h:275
PRunDataHandler * fRawData
Pointer to raw data handler (owned externally, not deleted here)
Definition PRunBase.h:273
PDoubleVector fT0s
Time-zero bin values for all histograms in this run (forward, backward, etc.)
Definition PRunBase.h:278
PRunBase()
Default constructor.
Definition PRunBase.cpp:54
Int_t fRunNo
Run number (0-based index in MSR file RUN blocks)
Definition PRunBase.h:270
PMsrRunBlock * fRunInfo
Pointer to this run's RUN block settings within fMsrInfo.
Definition PRunBase.h:272
Double_t fFitStartTime
Fit range start time in microseconds (μs) relative to t0.
Definition PRunBase.h:281
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.