musrfit 1.10.0
PRunSingleHisto.cpp
Go to the documentation of this file.
1/***************************************************************************
2
3 PRunSingleHisto.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
42#include <TString.h>
43#include <TObjArray.h>
44#include <TObjString.h>
45
46#include "PMusr.h"
47#include "PRunSingleHisto.h"
48
49//--------------------------------------------------------------------------
50// Constructor
51//--------------------------------------------------------------------------
66{
67 fScaleN0AndBkg = true;
68 fNoOfFitBins = 0;
69 fBackground = 0;
70 fPacking = -1;
71 fTheoAsData = false;
72
73 // the 2 following variables are need in case fit range is given in bins, and since
74 // the fit range can be changed in the command block, these variables need to be accessible
75 fGoodBins[0] = -1;
76 fGoodBins[1] = -1;
77
78 fStartTimeBin = -1;
79 fEndTimeBin = -1;
80}
81
82//--------------------------------------------------------------------------
83// Constructor
84//--------------------------------------------------------------------------
108PRunSingleHisto::PRunSingleHisto(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
109 PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
110{
112 fNoOfFitBins = 0;
113 fBackground = 0;
114
115 fPacking = fRunInfo->GetPacking();
116 if (fPacking == -1) { // i.e. packing is NOT given in the RUN-block, it must be given in the GLOBAL-block
117 fPacking = fMsrInfo->GetMsrGlobal()->GetPacking();
118 }
119 if (fPacking == -1) { // this should NOT happen, somethin is severely wrong
120 std::cerr << std::endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't find any packing information!";
121 std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
122 std::cerr << std::endl;
123 fValid = false;
124 return;
125 }
126
127 // the 2 following variables are need in case fit range is given in bins, and since
128 // the fit range can be changed in the command block, these variables need to be accessible
129 fGoodBins[0] = -1;
130 fGoodBins[1] = -1;
131
132 fStartTimeBin = -1;
133 fEndTimeBin = -1;
134
135 if (!PrepareData()) {
136 std::cerr << std::endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't prepare data for fitting!";
137 std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
138 std::cerr << std::endl;
139 fValid = false;
140 }
141}
142
143//--------------------------------------------------------------------------
144// Destructor
145//--------------------------------------------------------------------------
157
158//--------------------------------------------------------------------------
159// CalcChiSquare (public)
160//--------------------------------------------------------------------------
208Double_t PRunSingleHisto::CalcChiSquare(const std::vector<Double_t>& par)
209{
210 Double_t chisq = 0.0;
211 Double_t diff = 0.0;
212
213 Double_t N0 = 0.0;
214
215 // check if norm is a parameter or a function
216 if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
217 N0 = par[fRunInfo->GetNormParamNo()-1];
218 } else { // norm is a function
219 // get function number
220 UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
221 // evaluate function
222 N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
223 }
224
225 // get tau
226 Double_t tau;
227 if (fRunInfo->GetLifetimeParamNo() != -1)
228 tau = par[fRunInfo->GetLifetimeParamNo()-1];
229 else
230 tau = PMUON_LIFETIME;
231
232 // get background
233 Double_t bkg;
234 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
235 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
236 bkg = fBackground;
237 } else { // fixed bkg given
238 bkg = fRunInfo->GetBkgFix(0);
239 }
240 } else { // bkg fitted
241 bkg = par[fRunInfo->GetBkgFitParamNo()-1];
242 }
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) -
269 (N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg);
270 chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
271 }
272
273 // the correction factor is need since the data scales like pack*t_res,
274 // whereas the error scales like sqrt(pack*t_res)
275 if (fScaleN0AndBkg)
276 chisq *= fPacking * (fTimeResolution * 1.0e3);
277
278 return chisq;
279}
280
281//--------------------------------------------------------------------------
282// CalcChiSquareExpected (public)
283//--------------------------------------------------------------------------
323Double_t PRunSingleHisto::CalcChiSquareExpected(const std::vector<Double_t>& par)
324{
325 Double_t chisq = 0.0;
326 Double_t diff = 0.0;
327 Double_t theo = 0.0;
328
329 Double_t N0 = 0.0;
330
331 // check if norm is a parameter or a function
332 if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
333 N0 = par[fRunInfo->GetNormParamNo()-1];
334 } else { // norm is a function
335 // get function number
336 UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
337 // evaluate function
338 N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
339 }
340
341 // get tau
342 Double_t tau;
343 if (fRunInfo->GetLifetimeParamNo() != -1)
344 tau = par[fRunInfo->GetLifetimeParamNo()-1];
345 else
346 tau = PMUON_LIFETIME;
347
348 // get background
349 Double_t bkg;
350 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
351 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
352 bkg = fBackground;
353 } else { // fixed bkg given
354 bkg = fRunInfo->GetBkgFix(0);
355 }
356 } else { // bkg fitted
357 bkg = par[fRunInfo->GetBkgFitParamNo()-1];
358 }
359
360 // calculate functions
361 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
362 Int_t funcNo = fMsrInfo->GetFuncNo(i);
363 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
364 }
365
366 // calculate chi square
367 Double_t time(1.0);
368 Int_t i;
369
370 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
371 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
372 // for a given set of parameters---which should be done outside of the parallelized loop.
373 // For all other functions it means a tiny and acceptable overhead.
374 time = fTheory->Func(time, par, fFuncValues);
375
376 #ifdef HAVE_GOMP
377 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
378 if (chunk < 10)
379 chunk = 10;
380 #pragma omp parallel for default(shared) private(i,time,theo,diff) schedule(dynamic,chunk) reduction(+:chisq)
381 #endif
382 for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
383 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
384 theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
385 diff = fData.GetValue()->at(i) - theo;
386 chisq += diff*diff / theo;
387 }
388
389 // the correction factor is need since the data scales like pack*t_res,
390 // whereas the error scales like sqrt(pack*t_res)
391 if (fScaleN0AndBkg)
392 chisq *= fPacking * (fTimeResolution * 1.0e3);
393
394 return chisq;
395}
396
397//--------------------------------------------------------------------------
398// CalcMaxLikelihood (public)
399//--------------------------------------------------------------------------
455Double_t PRunSingleHisto::CalcMaxLikelihood(const std::vector<Double_t>& par)
456{
457 Double_t mllh = 0.0; // maximum log likelihood assuming poisson distribution for the single bin
458
459 Double_t N0;
460
461 // check if norm is a parameter or a function
462 if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
463 N0 = par[fRunInfo->GetNormParamNo()-1];
464 } else { // norm is a function
465 // get function number
466 Int_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
467 // evaluate function
468 N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
469 }
470
471 // get tau
472 Double_t tau;
473 if (fRunInfo->GetLifetimeParamNo() != -1)
474 tau = par[fRunInfo->GetLifetimeParamNo()-1];
475 else
476 tau = PMUON_LIFETIME;
477
478 // get background
479 Double_t bkg;
480 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
481 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
482 bkg = fBackground;
483 } else { // fixed bkg given
484 bkg = fRunInfo->GetBkgFix(0);
485 }
486 } else { // bkg fitted
487 bkg = par[fRunInfo->GetBkgFitParamNo()-1];
488 }
489
490 // calculate functions
491 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
492 UInt_t funcNo = fMsrInfo->GetFuncNo(i);
493 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
494 }
495
496 // calculate maximum log likelihood
497 Double_t theo;
498 Double_t data;
499 Double_t time(1.0);
500 Int_t i;
501
502 // norm is needed since there is no simple scaling like in chisq case to get the correct Max.Log.Likelihood value when normlizing N(t) to 1/ns
503 Double_t normalizer = 1.0;
504
505 if (fScaleN0AndBkg)
506 normalizer = fPacking * (fTimeResolution * 1.0e3);
507
508 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
509 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
510 // for a given set of parameters---which should be done outside of the parallelized loop.
511 // For all other functions it means a tiny and acceptable overhead.
512 time = fTheory->Func(time, par, fFuncValues);
513
514 #ifdef HAVE_GOMP
515 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
516 if (chunk < 10)
517 chunk = 10;
518 #pragma omp parallel for default(shared) private(i,time,theo,data) schedule(dynamic,chunk) reduction(-:mllh)
519 #endif
520 for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
521 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
522 // calculate theory for the given parameter set
523 theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
524
525 data = fData.GetValue()->at(i);
526
527 if (theo <= 0.0) {
528 std::cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << std::endl;
529 continue;
530 }
531
532 if (data > 1.0e-9) {
533 mllh += (theo-data) + data*log(data/theo);
534 } else {
535 mllh += (theo-data);
536 }
537 }
538
539 return normalizer*2.0*mllh;
540}
541
542//--------------------------------------------------------------------------
543// CalcMaxLikelihoodExpected (public)
544//--------------------------------------------------------------------------
590Double_t PRunSingleHisto::CalcMaxLikelihoodExpected(const std::vector<Double_t>& par)
591{
592 Double_t mllh = 0.0; // maximum log likelihood assuming poisson distribution for the single bin
593
594 Double_t N0;
595
596 // check if norm is a parameter or a function
597 if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
598 N0 = par[fRunInfo->GetNormParamNo()-1];
599 } else { // norm is a function
600 // get function number
601 Int_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
602 // evaluate function
603 N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
604 }
605
606 // get tau
607 Double_t tau;
608 if (fRunInfo->GetLifetimeParamNo() != -1)
609 tau = par[fRunInfo->GetLifetimeParamNo()-1];
610 else
611 tau = PMUON_LIFETIME;
612
613 // get background
614 Double_t bkg;
615 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
616 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
617 bkg = fBackground;
618 } else { // fixed bkg given
619 bkg = fRunInfo->GetBkgFix(0);
620 }
621 } else { // bkg fitted
622 bkg = par[fRunInfo->GetBkgFitParamNo()-1];
623 }
624
625 // calculate functions
626 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
627 UInt_t funcNo = fMsrInfo->GetFuncNo(i);
628 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
629 }
630
631 // calculate maximum log likelihood
632 Double_t theo;
633 Double_t data;
634 Double_t time(1.0);
635 Int_t i;
636
637 // norm is needed since there is no simple scaling like in chisq case to get the correct Max.Log.Likelihood value when normlizing N(t) to 1/ns
638 Double_t normalizer = 1.0;
639
640 if (fScaleN0AndBkg)
641 normalizer = fPacking * (fTimeResolution * 1.0e3);
642
643 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
644 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
645 // for a given set of parameters---which should be done outside of the parallelized loop.
646 // For all other functions it means a tiny and acceptable overhead.
647 time = fTheory->Func(time, par, fFuncValues);
648
649 #ifdef HAVE_GOMP
650 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
651 if (chunk < 10)
652 chunk = 10;
653 #pragma omp parallel for default(shared) private(i,time,theo,data) schedule(dynamic,chunk) reduction(-:mllh)
654 #endif
655 for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
656 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
657 // calculate theory for the given parameter set
658 theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
659
660 data = fData.GetValue()->at(i);
661
662 if (theo <= 0.0) {
663 std::cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << std::endl;
664 continue;
665 }
666
667 if (data > 1.0e-9) { // is this correct?? needs to be checked. See G-test
668 mllh += data*log(data/theo);
669 }
670 }
671
672 return normalizer*2.0*mllh;
673}
674
675//--------------------------------------------------------------------------
676// CalcTheory (public)
677//--------------------------------------------------------------------------
715{
716 // feed the parameter vector
717 std::vector<Double_t> par;
718 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
719 for (UInt_t i=0; i<paramList->size(); i++)
720 par.push_back((*paramList)[i].fValue);
721
722 // calculate asymmetry
723 Double_t N0;
724 // check if norm is a parameter or a function
725 if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
726 N0 = par[fRunInfo->GetNormParamNo()-1];
727 } else { // norm is a function
728 // get function number
729 Int_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
730 // evaluate function
731 N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
732 }
733
734 // get tau
735 Double_t tau;
736 if (fRunInfo->GetLifetimeParamNo() != -1)
737 tau = par[fRunInfo->GetLifetimeParamNo()-1];
738 else
739 tau = PMUON_LIFETIME;
740
741 // get background
742 Double_t bkg;
743 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
744 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
745 bkg = fBackground;
746 } else { // fixed bkg given
747 bkg = fRunInfo->GetBkgFix(0);
748 }
749 } else { // bkg fitted
750 bkg = par[fRunInfo->GetBkgFitParamNo()-1];
751 }
752
753 // calculate functions
754 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
755 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
756 }
757
758 // calculate theory
759 UInt_t size = fData.GetValue()->size();
760 Double_t start = fData.GetDataTimeStart();
761 Double_t resolution = fData.GetDataTimeStep();
762 Double_t time;
763 for (UInt_t i=0; i<size; i++) {
764 time = start + static_cast<Double_t>(i)*resolution;
765 fData.AppendTheoryValue(N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg);
766 }
767
768 // clean up
769 par.clear();
770}
771
772//--------------------------------------------------------------------------
773// GetNoOfFitBins (public)
774//--------------------------------------------------------------------------
793{
795
796 return fNoOfFitBins;
797}
798
799//--------------------------------------------------------------------------
800// SetFitRangeBin (public)
801//--------------------------------------------------------------------------
845void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
846{
847 TObjArray *tok = nullptr;
848 TObjString *ostr = nullptr;
849 TString str;
850 Ssiz_t idx = -1;
851 Int_t offset = 0;
852
853 tok = fitRange.Tokenize(" \t");
854
855 if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
856 // handle fgb+n0 entry
857 ostr = dynamic_cast<TObjString*>(tok->At(1));
858 str = ostr->GetString();
859 // check if there is an offset present
860 idx = str.First("+");
861 if (idx != -1) { // offset present
862 str.Remove(0, idx+1);
863 if (str.IsFloat()) // if str is a valid number, convert is to an integer
864 offset = str.Atoi();
865 }
866 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
867
868 // handle lgb-n1 entry
869 ostr = dynamic_cast<TObjString*>(tok->At(2));
870 str = ostr->GetString();
871 // check if there is an offset present
872 idx = str.First("-");
873 if (idx != -1) { // offset present
874 str.Remove(0, idx+1);
875 if (str.IsFloat()) // if str is a valid number, convert is to an integer
876 offset = str.Atoi();
877 }
878 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
879 } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
880 Int_t pos = 2*(fRunNo+1)-1;
881
882 if (pos + 1 >= tok->GetEntries()) {
883 std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
884 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
885 } else {
886 // handle fgb+n0 entry
887 ostr = dynamic_cast<TObjString*>(tok->At(pos));
888 str = ostr->GetString();
889 // check if there is an offset present
890 idx = str.First("+");
891 if (idx != -1) { // offset present
892 str.Remove(0, idx+1);
893 if (str.IsFloat()) // if str is a valid number, convert is to an integer
894 offset = str.Atoi();
895 }
896 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
897
898 // handle lgb-n1 entry
899 ostr = dynamic_cast<TObjString*>(tok->At(pos+1));
900 str = ostr->GetString();
901 // check if there is an offset present
902 idx = str.First("-");
903 if (idx != -1) { // offset present
904 str.Remove(0, idx+1);
905 if (str.IsFloat()) // if str is a valid number, convert is to an integer
906 offset = str.Atoi();
907 }
908 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
909 }
910 } else { // error
911 std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
912 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
913 }
914
915 // clean up
916 if (tok) {
917 delete tok;
918 }
919}
920
921//--------------------------------------------------------------------------
922// CalcNoOfFitBins (public)
923//--------------------------------------------------------------------------
952{
953 // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
954 fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
955 if (fStartTimeBin < 0)
956 fStartTimeBin = 0;
957 fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
958 if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
959 fEndTimeBin = fData.GetValue()->size();
960
963 else
964 fNoOfFitBins = 0;
965}
966
967//--------------------------------------------------------------------------
968// PrepareData (protected)
969//--------------------------------------------------------------------------
1012{
1013 Bool_t success = true;
1014
1015 if (!fValid)
1016 return false;
1017
1018 // keep the Global block info
1019 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
1020
1021 // get the proper run
1022 PRawRunData* runData = fRawData->GetRunData(*fRunInfo->GetRunName());
1023 if (!runData) { // couldn't get run
1024 std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
1025 std::cerr << std::endl;
1026 return false;
1027 }
1028
1029 // keep the field from the meta-data from the data-file
1030 fMetaData.fField = runData->GetField();
1031
1032 // keep the energy from the meta-data from the data-file
1033 fMetaData.fEnergy = runData->GetEnergy();
1034
1035 // keep the temperature(s) from the meta-data from the data-file
1036 for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
1037 fMetaData.fTemp.push_back(runData->GetTemperature(i));
1038
1039 // collect histogram numbers
1040 PUIntVector histoNo; // histoNo = msr-file forward + redGreen_offset - 1
1041 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
1042 histoNo.push_back(fRunInfo->GetForwardHistoNo(i));
1043
1044 if (!runData->IsPresent(histoNo[i])) {
1045 std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **PANIC ERROR**:";
1046 std::cerr << std::endl << ">> histoNo found = " << histoNo[i] << ", which is NOT present in the data file!?!?";
1047 std::cerr << std::endl << ">> Will quit :-(";
1048 std::cerr << std::endl;
1049 histoNo.clear();
1050 return false;
1051 }
1052 }
1053
1054 // keep the time resolution in (us)
1055 fTimeResolution = runData->GetTimeResolution()/1.0e3;
1056 std::cout.precision(10);
1057 std::cout << std::endl << ">> PRunSingleHisto::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
1058
1059 // get all the proper t0's and addt0's for the current RUN block
1060 if (!GetProperT0(runData, globalBlock, histoNo)) {
1061 return false;
1062 }
1063
1064 // keep the histo of each group at this point (addruns handled below)
1065 std::vector<PDoubleVector> forward;
1066 forward.resize(histoNo.size()); // resize to number of groups
1067 for (UInt_t i=0; i<histoNo.size(); i++) {
1068 forward[i].resize(runData->GetDataBin(histoNo[i])->size());
1069 forward[i] = *runData->GetDataBin(histoNo[i]);
1070 }
1071
1072 // check if a dead time correction has to be done
1073 // this will be done automatically in the function itself, which also
1074 // checks in the global and run section
1075 DeadTimeCorrection(forward, histoNo);
1076
1077 // check if there are runs to be added to the current one
1078 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1079 PRawRunData *addRunData;
1080 std::vector<PDoubleVector> addForward;
1081 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) { // loop over all ADDRUN's
1082
1083 // get run to be added to the main one
1084 addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
1085 if (addRunData == nullptr) { // couldn't get run
1086 std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1087 std::cerr << std::endl;
1088 return false;
1089 }
1090
1091 addForward.clear();
1092 addForward.resize(histoNo.size()); // resize to number of groups
1093 for (UInt_t j=0; j<histoNo.size(); j++) {
1094 addForward[j].resize(addRunData->GetDataBin(histoNo[j])->size());
1095 addForward[j] = *addRunData->GetDataBin(histoNo[j]);
1096 }
1097 DeadTimeCorrection(addForward, histoNo);
1098
1099 // add forward run
1100 UInt_t addRunSize;
1101 for (UInt_t k=0; k<histoNo.size(); k++) { // fill each group
1102 addRunSize = addForward[k].size();
1103 for (UInt_t j=0; j<addRunSize; j++) { // loop over the bin indices
1104 // make sure that the index stays in the proper range
1105 if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) >= 0) &&
1106 (j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) < addRunSize)) {
1107 forward[k][j] += addForward[k][j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k])];
1108 }
1109 }
1110 }
1111 }
1112 }
1113
1114 // set forward histo data of the first group
1115 fForward.resize(forward[0].size());
1116 for (UInt_t i=0; i<fForward.size(); i++) {
1117 fForward[i] = forward[0][i];
1118 }
1119
1120 // group histograms, add all the remaining forward histograms of the group
1121 for (UInt_t i=1; i<histoNo.size(); i++) { // loop over the groupings
1122 for (UInt_t j=0; j<runData->GetDataBin(histoNo[i])->size(); j++) { // loop over the bin indices
1123 // make sure that the index stays within proper range
1124 if ((static_cast<Int_t>(j)+fT0s[i]-fT0s[0] >= 0) && (j+fT0s[i]-fT0s[0] < runData->GetDataBin(histoNo[i])->size())) {
1125 fForward[j] += forward[i][j+static_cast<Int_t>(fT0s[i])-static_cast<Int_t>(fT0s[0])];
1126 }
1127 }
1128 }
1129
1130 // get the data range (fgb/lgb) for the current RUN block
1131 if (!GetProperDataRange()) {
1132 return false;
1133 }
1134
1135 // get the fit range for the current RUN block
1136 GetProperFitRange(globalBlock);
1137
1138 // get the lifetimecorrection flag
1139 Bool_t lifetimecorrection = false;
1140 PMsrPlotList *plot = fMsrInfo->GetMsrPlotList();
1141 lifetimecorrection = plot->at(0).fLifeTimeCorrection;
1142
1143 // do the more fit/view specific stuff
1144 if (fHandleTag == kFit)
1145 success = PrepareFitData(runData, histoNo[0]);
1146 else if ((fHandleTag == kView) && !lifetimecorrection)
1147 success = PrepareRawViewData(runData, histoNo[0]);
1148 else if ((fHandleTag == kView) && lifetimecorrection)
1149 success = PrepareViewData(runData, histoNo[0]);
1150 else
1151 success = false;
1152
1153 // cleanup
1154 histoNo.clear();
1155
1156 return success;
1157}
1158
1159//--------------------------------------------------------------------------
1160// PrepareFitData (protected)
1161//--------------------------------------------------------------------------
1213Bool_t PRunSingleHisto::PrepareFitData(PRawRunData* runData, const UInt_t histoNo)
1214{
1215 if (fMsrInfo->EstimateN0()) {
1216 EstimateN0();
1217 }
1218
1219 // transform raw histo data. This is done the following way (for details see the manual):
1220 // for the single histo fit, just the rebinned raw data are copied
1221
1222 // check how the background shall be handled
1223 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg shall **NOT** be fitted
1224 // subtract background from histogramms ------------------------------------------
1225 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
1226 if (fRunInfo->GetBkgRange(0) >= 0) {
1227 if (!EstimateBkg(histoNo))
1228 return false;
1229 } else { // no background given to do the job, try estimate
1230 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
1231 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
1232 std::cerr << std::endl << ">> PRunSingleHisto::PrepareFitData(): **WARNING** Neither fix background nor background bins are given!";
1233 std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
1234 std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
1235 std::cerr << std::endl;
1236 if (!EstimateBkg(histoNo))
1237 return false;
1238 }
1239 } else { // fixed background given
1240 for (UInt_t i=0; i<fForward.size(); i++) {
1241 fForward[i] -= fRunInfo->GetBkgFix(0);
1242 }
1243 }
1244 }
1245
1246 // everything looks fine, hence fill data set
1247 Int_t t0 = static_cast<Int_t>(fT0s[0]);
1248 Double_t value = 0.0;
1249 Double_t normalizer = 1.0;
1250 // in order that after rebinning the fit does not need to be redone (important for plots)
1251 // the value is normalize to per 1 nsec if scaling is whished
1252 if (fScaleN0AndBkg)
1253 normalizer = fPacking * (fTimeResolution * 1.0e3); // fTimeResolution us->ns
1254 // data start at data_start-t0
1255 // time shifted so that packing is included correctly, i.e. t0 == t0 after packing
1256 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(fGoodBins[0])-0.5) + static_cast<Double_t>(fPacking)/2.0 - static_cast<Double_t>(t0)));
1257 fData.SetDataTimeStep(fTimeResolution*fPacking);
1258 for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
1259 if (fPacking == 1) {
1260 value = fForward[i];
1261 value /= normalizer;
1262 fData.AppendValue(value);
1263 if (value == 0.0)
1264 fData.AppendErrorValue(1.0/normalizer);
1265 else
1266 fData.AppendErrorValue(TMath::Sqrt(value));
1267 } else { // packed data, i.e. fPacking > 1
1268 if (((i-fGoodBins[0]) % fPacking == 0) && (i != fGoodBins[0])) { // fill data
1269 value /= normalizer;
1270 fData.AppendValue(value);
1271 if (value == 0.0)
1272 fData.AppendErrorValue(1.0/normalizer);
1273 else
1274 fData.AppendErrorValue(TMath::Sqrt(value));
1275 // reset values
1276 value = 0.0;
1277 }
1278 value += fForward[i];
1279 }
1280 }
1281
1283
1284 return true;
1285}
1286
1287//--------------------------------------------------------------------------
1288// PrepareRawViewData (protected)
1289//--------------------------------------------------------------------------
1306Bool_t PRunSingleHisto::PrepareRawViewData(PRawRunData* runData, const UInt_t histoNo)
1307{
1308 // check if view_packing is wished
1309 Int_t packing = fPacking;
1310 if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
1311 packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1312 }
1313
1314 // calculate necessary norms
1315 Double_t dataNorm = 1.0, theoryNorm = 1.0;
1316 if (fScaleN0AndBkg) {
1317 dataNorm = 1.0/ (packing * (fTimeResolution * 1.0e3)); // fTimeResolution us->ns
1318 } else if (!fScaleN0AndBkg && (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0)) {
1319 theoryNorm = static_cast<Double_t>(fMsrInfo->GetMsrPlotList()->at(0).fViewPacking)/static_cast<Double_t>(fPacking);
1320 }
1321
1322 // raw data, since PMusrCanvas is doing ranging etc.
1323 // start = the first bin which is a multiple of packing backward from first good data bin
1324 Int_t start = fGoodBins[0] - (fGoodBins[0]/packing)*packing;
1325 // end = last bin starting from start which is a multiple of packing and still within the data
1326 Int_t end = start + ((fForward.size()-start)/packing)*packing;
1327 // check if data range has been provided, and if not try to estimate them
1328 if (start < 0) {
1329 Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
1330 start = (static_cast<Int_t>(fT0s[0])+offset) - ((static_cast<Int_t>(fT0s[0])+offset)/packing)*packing;
1331 end = start + ((fForward.size()-start)/packing)*packing;
1332 std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
1333 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1334 std::cerr << std::endl;
1335 }
1336 // check if start, end, and t0 make any sense
1337 // 1st check if start and end are in proper order
1338 if (end < start) { // need to swap them
1339 Int_t keep = end;
1340 end = start;
1341 start = keep;
1342 }
1343 // 2nd check if start is within proper bounds
1344 if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
1345 std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** start data bin doesn't make any sense!";
1346 std::cerr << std::endl;
1347 return false;
1348 }
1349 // 3rd check if end is within proper bounds
1350 if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
1351 std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** end data bin doesn't make any sense!";
1352 std::cerr << std::endl;
1353 return false;
1354 }
1355
1356 // everything looks fine, hence fill data set
1357 Int_t t0 = static_cast<Int_t>(fT0s[0]);
1358 Double_t value = 0.0;
1359 // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1360 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start)-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0)));
1361 fData.SetDataTimeStep(fTimeResolution*packing);
1362
1363 for (Int_t i=start; i<end; i++) {
1364 if (((i-start) % packing == 0) && (i != start)) { // fill data
1365 value *= dataNorm;
1366 fData.AppendValue(value);
1367 if (value == 0.0)
1368 fData.AppendErrorValue(1.0);
1369 else
1370 fData.AppendErrorValue(TMath::Sqrt(value*dataNorm));
1371 // reset values
1372 value = 0.0;
1373 }
1374 value += fForward[i];
1375 }
1376
1378
1379 // fill theory vector for kView
1380 // feed the parameter vector
1381 std::vector<Double_t> par;
1382 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1383 for (UInt_t i=0; i<paramList->size(); i++)
1384 par.push_back((*paramList)[i].fValue);
1385
1386 // calculate asymmetry
1387 Double_t N0;
1388 // check if norm is a parameter or a function
1389 if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
1390 N0 = par[fRunInfo->GetNormParamNo()-1];
1391 } else { // norm is a function
1392 // get function number
1393 UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
1394 // evaluate function
1395 N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1396 }
1397 N0 *= theoryNorm;
1398
1399 // get tau
1400 Double_t tau;
1401 if (fRunInfo->GetLifetimeParamNo() != -1)
1402 tau = par[fRunInfo->GetLifetimeParamNo()-1];
1403 else
1404 tau = PMUON_LIFETIME;
1405
1406 // get background
1407 Double_t bkg;
1408 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
1409 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
1410 if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
1411 if (!EstimateBkg(histoNo))
1412 return false;
1413 } else { // no background given to do the job, try estimate
1414 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
1415 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
1416 std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** Neither fix background nor background bins are given!";
1417 std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
1418 std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
1419 std::cerr << std::endl;
1420 if (!EstimateBkg(histoNo))
1421 return false;
1422 }
1423 bkg = fBackground;
1424 } else { // fixed bkg given
1425 bkg = fRunInfo->GetBkgFix(0);
1426 }
1427 } else { // bkg fitted
1428 bkg = par[fRunInfo->GetBkgFitParamNo()-1];
1429 }
1430 bkg *= theoryNorm;
1431
1432 // calculate functions
1433 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1434 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
1435 }
1436
1437 // calculate theory
1438 UInt_t size = fForward.size();
1439 Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1440 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1441 if (fTheoAsData) { // calculate theory only at the data points
1442 fData.SetTheoryTimeStep(fData.GetDataTimeStep());
1443 } else {
1444 // finer binning for the theory (8 times as many points = factor)
1445 size *= factor;
1446 fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1447 }
1448
1449 Double_t time;
1450 Double_t theoryValue;
1451 for (UInt_t i=0; i<size; i++) {
1452 time = fData.GetTheoryTimeStart() + i*fData.GetTheoryTimeStep();
1453 theoryValue = fTheory->Func(time, par, fFuncValues);
1454 if (fabs(theoryValue) > 1.0e10) { // dirty hack needs to be fixed!!
1455 theoryValue = 0.0;
1456 }
1457 fData.AppendTheoryValue(N0*TMath::Exp(-time/tau)*(1.0+theoryValue)+bkg);
1458 }
1459
1460 // clean up
1461 par.clear();
1462
1463 return true;
1464}
1465
1466//--------------------------------------------------------------------------
1467// PrepareViewData (protected)
1468//--------------------------------------------------------------------------
1495Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histoNo)
1496{
1497 // check if view_packing is wished. This is a global option for all PLOT blocks!
1498 Int_t packing = fPacking;
1499 if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
1500 packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1501 }
1502 // check if rrf_packing is present. This is a global option for all PLOT blocks, since operated on a single set of data.
1503 if (fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking > 0) {
1504 packing = fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking;
1505 }
1506
1507 // calculate necessary norms
1508 Double_t dataNorm = 1.0, theoryNorm = 1.0;
1509 if (fScaleN0AndBkg) {
1510 dataNorm = 1.0/ (packing * (fTimeResolution * 1.0e3)); // fTimeResolution us->ns
1511 } else if (!fScaleN0AndBkg && (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0)) {
1512 theoryNorm = static_cast<Double_t>(fMsrInfo->GetMsrPlotList()->at(0).fViewPacking)/static_cast<Double_t>(fPacking);
1513 }
1514
1515 // transform raw histo data. This is done the following way (for details see the manual):
1516 // for the single histo fit, just the rebinned raw data are copied
1517 // first get start data, end data, and t0
1518 Int_t t0 = static_cast<Int_t>(fT0s[0]);
1519
1520 // start = the first bin which is a multiple of packing backward from first good data bin
1521 Int_t start = fGoodBins[0] - (fGoodBins[0]/packing)*packing;
1522 // end = last bin starting from start which is a multiple of packing and still within the data
1523 Int_t end = start + ((fForward.size()-start)/packing)*packing;
1524
1525 // check if data range has been provided, and if not try to estimate them
1526 if (start < 0) {
1527 Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
1528 start = (static_cast<Int_t>(fT0s[0])+offset) - ((static_cast<Int_t>(fT0s[0])+offset)/packing)*packing;
1529 end = start + ((fForward.size()-start)/packing)*packing;
1530 std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
1531 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1532 std::cerr << std::endl;
1533 }
1534
1535 // check if start, end, and t0 make any sense
1536 // 1st check if start and end are in proper order
1537 if (end < start) { // need to swap them
1538 Int_t keep = end;
1539 end = start;
1540 start = keep;
1541 }
1542 // 2nd check if start is within proper bounds
1543 if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
1544 std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
1545 std::cerr << std::endl;
1546 return false;
1547 }
1548 // 3rd check if end is within proper bounds
1549 if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
1550 std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
1551 std::cerr << std::endl;
1552 return false;
1553 }
1554
1555 // everything looks fine, hence fill data set
1556
1557 // feed the parameter vector
1558 std::vector<Double_t> par;
1559 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1560 for (UInt_t i=0; i<paramList->size(); i++)
1561 par.push_back((*paramList)[i].fValue);
1562
1563 // calculate asymmetry
1564 Double_t N0;
1565 // check if norm is a parameter or a function
1566 if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
1567 N0 = par[fRunInfo->GetNormParamNo()-1];
1568 } else { // norm is a function
1569 // get function number
1570 UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
1571 // evaluate function
1572 N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1573 }
1574 N0 *= theoryNorm;
1575
1576 // get tau
1577 Double_t tau;
1578 if (fRunInfo->GetLifetimeParamNo() != -1)
1579 tau = par[fRunInfo->GetLifetimeParamNo()-1];
1580 else
1581 tau = PMUON_LIFETIME;
1582
1583 // get background
1584 Double_t bkg;
1585 if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
1586 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
1587 if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
1588 if (!EstimateBkg(histoNo))
1589 return false;
1590 } else { // no background given to do the job, try estimate
1591 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
1592 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
1593 std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** Neither fix background nor background bins are given!";
1594 std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
1595 std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
1596 std::cerr << std::endl;
1597 if (!EstimateBkg(histoNo))
1598 return false;
1599 }
1600 bkg = fBackground;
1601 } else { // fixed bkg given
1602 bkg = fRunInfo->GetBkgFix(0);
1603 }
1604 } else { // bkg fitted
1605 bkg = par[fRunInfo->GetBkgFitParamNo()-1];
1606 }
1607 bkg *= theoryNorm;
1608
1609 Double_t value = 0.0;
1610 Double_t expval = 0.0;
1611 Double_t rrf_val = 0.0;
1612 Double_t time = 0.0;
1613
1614 // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1615 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start)-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0)));
1616 fData.SetDataTimeStep(fTimeResolution*packing);
1617
1618 // data is always normalized to (per nsec!!)
1619 Double_t gammaRRF = 0.0, wRRF = 0.0, phaseRRF = 0.0;
1620 if (fMsrInfo->GetMsrPlotList()->at(0).fRRFFreq == 0.0) { // normal Data representation
1621 for (Int_t i=start; i<end; i++) {
1622 if (((i-start) % packing == 0) && (i != start)) { // fill data
1623 value *= dataNorm;
1624 // since the packing counter is already at the end of the bin, the time needs be shifted back by pack*time_resolution
1625 time = (((static_cast<Double_t>(i)-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0)))*fTimeResolution - static_cast<Double_t>(packing)*fTimeResolution;
1626 expval = TMath::Exp(+time/tau)/N0;
1627 fData.AppendValue(-1.0+expval*(value-bkg));
1628 fData.AppendErrorValue(expval*TMath::Sqrt(value*dataNorm));
1629 value = 0.0;
1630 }
1631 value += fForward[i];
1632 }
1633 } else { // RRF representation
1634 // check which units shall be used
1635 switch (fMsrInfo->GetMsrPlotList()->at(0).fRRFUnit) {
1636 case RRF_UNIT_kHz:
1637 gammaRRF = TMath::TwoPi()*1.0e-3;
1638 break;
1639 case RRF_UNIT_MHz:
1640 gammaRRF = TMath::TwoPi();
1641 break;
1642 case RRF_UNIT_Mcs:
1643 gammaRRF = 1.0;
1644 break;
1645 case RRF_UNIT_G:
1646 gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi();
1647 break;
1648 case RRF_UNIT_T:
1649 gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi()*1.0e4;
1650 break;
1651 default:
1652 gammaRRF = TMath::TwoPi();
1653 break;
1654 }
1655 wRRF = gammaRRF * fMsrInfo->GetMsrPlotList()->at(0).fRRFFreq;
1656 phaseRRF = fMsrInfo->GetMsrPlotList()->at(0).fRRFPhase / 180.0 * TMath::Pi();
1657
1658 Double_t error = 0.0;
1659 for (Int_t i=start; i<end; i++) {
1660 if (((i-start) % packing == 0) && (i != start)) { // fill data
1661 fData.AppendValue(2.0*value/packing); // factor 2 needed because cos(a)cos(b) = 1/2(cos(a+b)+cos(a-b))
1662 fData.AppendErrorValue(expval*TMath::Sqrt(error/packing));
1663 value = 0.0;
1664 error = 0.0;
1665 }
1666 time = (static_cast<Double_t>(i)-t0)*fTimeResolution;
1667 expval = TMath::Exp(+time/tau)/N0;
1668 rrf_val = (-1.0+expval*(fForward[i]/(fTimeResolution*1.0e3)-bkg))*TMath::Cos(wRRF * time + phaseRRF);
1669 value += rrf_val;
1670 error += fForward[i]*dataNorm;
1671 }
1672 }
1673
1675
1676 // calculate functions
1677 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1678 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
1679 }
1680
1681 // calculate theory
1682 Double_t theoryValue;
1683 UInt_t size = fForward.size()/packing;
1684 const Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1685 UInt_t rebinRRF = 0;
1686
1687 if (wRRF == 0) { // no RRF
1688 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1689 if (fTheoAsData) { // calculate theory only at the data points
1690 fData.SetTheoryTimeStep(fData.GetDataTimeStep());
1691 } else {
1692 // finer binning for the theory (8 times as many points = factor)
1693 size *= factor;
1694 fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1695 }
1696 } else { // RRF
1697 rebinRRF = static_cast<UInt_t>((TMath::Pi()/2.0/wRRF)/fTimeResolution); // RRF time resolution / data time resolution
1698 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1699 fData.SetTheoryTimeStep(TMath::Pi()/2.0/wRRF/rebinRRF); // = theory time resolution as close as possible to the data time resolution compatible with wRRF
1700 }
1701
1702 for (UInt_t i=0; i<size; i++) {
1703 time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1704 theoryValue = fTheory->Func(time, par, fFuncValues);
1705 if (wRRF != 0.0) {
1706 theoryValue *= 2.0*TMath::Cos(wRRF * time + phaseRRF);
1707 }
1708 if (fabs(theoryValue) > 10.0) { // dirty hack needs to be fixed!!
1709 theoryValue = 0.0;
1710 }
1711 fData.AppendTheoryValue(theoryValue);
1712 }
1713
1714 // if RRF filter the theory with a FIR Kaiser low pass filter
1715 if (wRRF != 0.0) {
1716 // rebin theory to the RRF frequency
1717 if (rebinRRF != 0) {
1718 Double_t dval = 0.0;
1719 PDoubleVector theo;
1720 for (UInt_t i=0; i<fData.GetTheory()->size(); i++) {
1721 if ((i % rebinRRF == 0) && (i != 0)) {
1722 theo.push_back(dval/rebinRRF);
1723 dval = 0.0;
1724 }
1725 dval += fData.GetTheory()->at(i);
1726 }
1727 fData.SetTheoryTimeStart(fData.GetTheoryTimeStart()+static_cast<Double_t>(rebinRRF-1)*fData.GetTheoryTimeStep()/2.0);
1728 fData.SetTheoryTimeStep(rebinRRF*fData.GetTheoryTimeStep());
1729 fData.ReplaceTheory(theo);
1730 theo.clear();
1731 }
1732
1733 // filter theory
1734 CalculateKaiserFilterCoeff(wRRF, 60.0, 0.2); // w_c = wRRF, A = -20 log_10(delta), Delta w / w_c = (w_s - w_p) / (2 w_c)
1735 FilterTheo();
1736 }
1737
1738 // clean up
1739 par.clear();
1740
1741 return true;
1742}
1743
1744//--------------------------------------------------------------------------
1745// GetProperT0 (private)
1746//--------------------------------------------------------------------------
1798{
1799 // feed all T0's
1800 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1801 fT0s.clear();
1802 fT0s.resize(histoNo.size());
1803 for (UInt_t i=0; i<fT0s.size(); i++) {
1804 fT0s[i] = -1.0;
1805 }
1806
1807 // fill in the T0's from the msr-file (if present)
1808 for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1809 fT0s[i] = fRunInfo->GetT0Bin(i);
1810 }
1811
1812 // fill in the T0's from the GLOBAL block section (if present)
1813 for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1814 if (fT0s[i] == -1.0) { // i.e. not given in the RUN block section
1815 fT0s[i] = globalBlock->GetT0Bin(i);
1816 }
1817 }
1818
1819 // fill in the T0's from the data file, if not already present in the msr-file
1820 for (UInt_t i=0; i<histoNo.size(); i++) {
1821 if (fT0s[i] == -1.0) { // i.e. not present in the msr-file, try the data file
1822 if (runData->GetT0Bin(histoNo[i]) > 0.0) {
1823 fT0s[i] = runData->GetT0Bin(histoNo[i]);
1824 fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
1825 }
1826 }
1827 }
1828
1829 // 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
1830 for (UInt_t i=0; i<histoNo.size(); i++) {
1831 if (fT0s[i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1832 fT0s[i] = runData->GetT0BinEstimated(histoNo[i]);
1833 fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
1834
1835 std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1836 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1837 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(histoNo[i]);
1838 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1839 std::cerr << std::endl;
1840 }
1841 }
1842
1843 // check if t0 is within proper bounds
1844 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
1845 if ((fT0s[i] < 0.0) || (fT0s[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1846 std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** t0 data bin (" << fT0s[i] << ") doesn't make any sense!";
1847 std::cerr << std::endl;
1848 return false;
1849 }
1850 }
1851
1852 // check if there are runs to be added to the current one. If yes keep the needed t0's
1853 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1854 PRawRunData *addRunData;
1855 fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1856 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1857
1858 // get run to be added to the main one
1859 addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
1860 if (addRunData == nullptr) { // couldn't get run
1861 std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1862 std::cerr << std::endl;
1863 return false;
1864 }
1865
1866 // feed all T0's
1867 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1868 fAddT0s[i-1].resize(histoNo.size());
1869 for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1870 fAddT0s[i-1][j] = -1.0;
1871 }
1872
1873 // fill in the T0's from the msr-file (if present)
1874 for (UInt_t j=0; j<fRunInfo->GetT0BinSize(); j++) {
1875 fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1,j); // addRunIdx starts at 0
1876 }
1877
1878 // fill in the T0's from the data file, if not already present in the msr-file
1879 for (UInt_t j=0; j<histoNo.size(); j++) {
1880 if (fAddT0s[i-1][j] == -1.0) // i.e. not present in the msr-file, try the data file
1881 if (addRunData->GetT0Bin(histoNo[j]) > 0.0) {
1882 fAddT0s[i-1][j] = addRunData->GetT0Bin(histoNo[j]);
1883 fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1884 }
1885 }
1886
1887 // 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
1888 for (UInt_t j=0; j<histoNo.size(); j++) {
1889 if (fAddT0s[i-1][j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1890 fAddT0s[i-1][j] = addRunData->GetT0BinEstimated(histoNo[j]);
1891 fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1892
1893 std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1894 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1895 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(histoNo[j]);
1896 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1897 std::cerr << std::endl;
1898 }
1899 }
1900
1901 // check if t0 is within proper bounds
1902 for (UInt_t j=0; j<fRunInfo->GetForwardHistoNoSize(); j++) {
1903 if ((fAddT0s[i-1][j] < 0.0) || (fAddT0s[i-1][j] > static_cast<Int_t>(addRunData->GetDataBin(histoNo[j])->size()))) {
1904 std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** addt0 data bin (" << fAddT0s[i-1][j] << ") doesn't make any sense!";
1905 std::cerr << std::endl;
1906 return false;
1907 }
1908 }
1909 }
1910 }
1911
1912 return true;
1913}
1914
1915//--------------------------------------------------------------------------
1916// GetProperDataRange (private)
1917//--------------------------------------------------------------------------
1963{
1964 // get start/end data
1965 Int_t start;
1966 Int_t end;
1967 start = fRunInfo->GetDataRange(0);
1968 end = fRunInfo->GetDataRange(1);
1969
1970 // check if data range has been given in the RUN block, if not try to get it from the GLOBAL block
1971 if (start < 0) {
1972 start = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1973 }
1974 if (end < 0) {
1975 end = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1976 }
1977
1978 // check if data range has been provided, and if not try to estimate them
1979 if (start < 0) {
1980 Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
1981 start = static_cast<Int_t>(fT0s[0])+offset;
1982 fRunInfo->SetDataRange(start, 0);
1983 std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start << ".";
1984 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1985 std::cerr << std::endl;
1986 }
1987 if (end < 0) {
1988 end = fForward.size();
1989 fRunInfo->SetDataRange(end, 1);
1990 std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range end = " << end << ".";
1991 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1992 std::cerr << std::endl;
1993 }
1994
1995 // check if start and end make any sense
1996 // 1st check if start and end are in proper order
1997 if (end < start) { // need to swap them
1998 Int_t keep = end;
1999 end = start;
2000 start = keep;
2001 }
2002 // 2nd check if start is within proper bounds
2003 if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
2004 std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** start data bin (" << start << ") doesn't make any sense!";
2005 std::cerr << std::endl;
2006 return false;
2007 }
2008 // 3rd check if end is within proper bounds
2009 if (end < 0) {
2010 std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** end data bin (" << end << ") doesn't make any sense!";
2011 std::cerr << std::endl;
2012 return false;
2013 }
2014 if (end > static_cast<Int_t>(fForward.size())) {
2015 std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** end data bin (" << end << ") > histo length (" << fForward.size() << ").";
2016 std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
2017 std::cerr << std::endl;
2018 end = static_cast<Int_t>(fForward.size())-1;
2019 }
2020
2021 // keep good bins for potential later use
2022 fGoodBins[0] = start;
2023 fGoodBins[1] = end;
2024
2025 // make sure that fGoodBins are in proper range for fForward
2026 if (fGoodBins[0] < 0)
2027 fGoodBins[0]=0;
2028 if (fGoodBins[1] > fForward.size()) {
2029 std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange **WARNING** needed to shift forward lgb,";
2030 std::cerr << std::endl << ">> from " << fGoodBins[1] << " to " << fForward.size()-1 << std::endl;
2031 fGoodBins[1]=fForward.size()-1;
2032 }
2033
2034 return true;
2035}
2036
2037//--------------------------------------------------------------------------
2038// GetProperFitRange (private)
2039//--------------------------------------------------------------------------
2097{
2098 // set fit start/end time; first check RUN Block
2099 fFitStartTime = fRunInfo->GetFitRange(0);
2100 fFitEndTime = fRunInfo->GetFitRange(1);
2101 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
2102 if (fRunInfo->IsFitRangeInBin()) {
2103 fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
2104 fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
2105 // write these times back into the data structure. This way it is available when writting the log-file
2106 fRunInfo->SetFitRange(fFitStartTime, 0);
2107 fRunInfo->SetFitRange(fFitEndTime, 1);
2108 }
2109 if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
2110 fFitStartTime = globalBlock->GetFitRange(0);
2111 fFitEndTime = globalBlock->GetFitRange(1);
2112 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
2113 if (globalBlock->IsFitRangeInBin()) {
2114 fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
2115 fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
2116 // write these times back into the data structure. This way it is available when writting the log-file
2117 globalBlock->SetFitRange(fFitStartTime, 0);
2118 globalBlock->SetFitRange(fFitEndTime, 1);
2119 }
2120 }
2122 fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
2123 fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
2124 std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
2125 std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
2126 }
2127}
2128
2129//--------------------------------------------------------------------------
2130// EstimateN0 (private)
2131//--------------------------------------------------------------------------
2175{
2176 // check that 'norm' in the msr-file run block is indeed a parameter number.
2177 // in case it is a function, nothing will be done.
2178 UInt_t paramNo = fRunInfo->GetNormParamNo();
2179 if (paramNo > 10000) // i.e. fun or map
2180 return;
2181
2182 // get the parameters
2183 PMsrParamList *param = fMsrInfo->GetMsrParamList();
2184 assert(param);
2185
2186 if (paramNo > param->size()) {
2187 std::cerr << std::endl << ">> PRunSingleHisto::EstimateN0: **ERROR** found parameter number " << paramNo << ", which is larger than the number of parameters = " << param->size() << std::endl;
2188 return;
2189 }
2190
2191 // check if N0 is fixed. If this is the case, do NOT estimate N0
2192 if (param->at(paramNo-1).fStep == 0.0) // N0 parameter fixed
2193 return;
2194
2195
2196 // check that 'backgr.fit' in the msr-file run block is indeed a parameter number.
2197 // in case it is a function, nothing will be done.
2198 Int_t paramNoBkg = fRunInfo->GetBkgFitParamNo();
2199 Bool_t scaleBkg = true;
2200 Double_t bkg=0.0, errBkg=1.0;
2201 if ((paramNoBkg > 10000) || (paramNoBkg == -1)) { // i.e. fun or map
2202 scaleBkg = false;
2203 } else {
2204 if (paramNoBkg-1 < static_cast<Int_t>(param->size())) {
2205 bkg = param->at(paramNoBkg-1).fValue;
2206 errBkg = param->at(paramNoBkg-1).fStep;
2207 }
2208 }
2209
2210 // estimate N0
2211 Double_t dt = fTimeResolution;
2212 Double_t tau = PMUON_LIFETIME;
2213
2214 UInt_t t0 = static_cast<UInt_t>(round(fT0s[0]));
2215 Double_t dval = 0.0;
2216 Double_t nom = 0.0;
2217 Double_t denom = 0.0;
2218 Double_t xx = 0.0;
2219
2220 // calc nominator
2221 for (UInt_t i=t0; i<fForward.size(); i++) {
2222 xx = exp(-dt*static_cast<Double_t>(i-t0)/tau);
2223 nom += xx;
2224 }
2225
2226 // calc denominator
2227 for (UInt_t i=t0; i<fForward.size(); i++) {
2228 xx = exp(-dt*static_cast<Double_t>(i-t0)/tau);
2229 dval = fForward[i];
2230 if (dval > 0)
2231 denom += xx*xx/dval;
2232 }
2233 Double_t N0 = nom/denom;
2234
2235 if (fScaleN0AndBkg) {
2236 N0 /= fTimeResolution*1.0e3;
2237 } else {
2238 N0 *= fPacking;
2239 }
2240
2241 Double_t rescale = 1;
2242 if ((param->at(paramNo-1).fValue != 0.0) && scaleBkg) {
2243 rescale = N0 / param->at(paramNo-1).fValue;
2244 bkg *= rescale;
2245 errBkg *= rescale;
2246 }
2247
2248 std::cout << ">> PRunSingleHisto::EstimateN0: found N0=" << param->at(paramNo-1).fValue << ", will set it to N0=" << N0 << std::endl;
2249 if (scaleBkg)
2250 std::cout << ">> PRunSingleHisto::EstimateN0: found Bkg=" << param->at(paramNoBkg-1).fValue << ", will set it to Bkg=" << bkg << std::endl;
2251 fMsrInfo->SetMsrParamValue(paramNo-1, N0);
2252 fMsrInfo->SetMsrParamStep(paramNo-1, sqrt(fabs(N0)));
2253 if (scaleBkg) {
2254 fMsrInfo->SetMsrParamValue(paramNoBkg-1, bkg);
2255 fMsrInfo->SetMsrParamStep(paramNoBkg-1, errBkg);
2256 }
2257}
2258
2259//--------------------------------------------------------------------------
2260// EstimateBkg (private)
2261//--------------------------------------------------------------------------
2310Bool_t PRunSingleHisto::EstimateBkg(UInt_t histoNo)
2311{
2312 Double_t beamPeriod = 0.0;
2313
2314 // check if data are from PSI, RAL, or TRIUMF
2315 if (fRunInfo->GetInstitute()->Contains("psi"))
2316 beamPeriod = ACCEL_PERIOD_PSI;
2317 else if (fRunInfo->GetInstitute()->Contains("ral"))
2318 beamPeriod = ACCEL_PERIOD_RAL;
2319 else if (fRunInfo->GetInstitute()->Contains("triumf"))
2320 beamPeriod = ACCEL_PERIOD_TRIUMF;
2321 else
2322 beamPeriod = 0.0;
2323
2324 // check if start and end are in proper order
2325 Int_t start = fRunInfo->GetBkgRange(0);
2326 Int_t end = fRunInfo->GetBkgRange(1);
2327 if (end < start) {
2328 std::cout << std::endl << "PRunSingleHisto::EstimatBkg(): end = " << end << " > start = " << start << "! Will swap them!";
2329 Int_t keep = end;
2330 end = start;
2331 start = keep;
2332 }
2333
2334 // calculate proper background range
2335 if (beamPeriod != 0.0) {
2336 Double_t timeBkg = static_cast<Double_t>(end-start)*(fTimeResolution*fPacking); // length of the background intervall in time
2337 UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
2338 // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
2339 end = start + static_cast<UInt_t>((fullCycles*beamPeriod)/(fTimeResolution*fPacking));
2340 std::cout << std::endl << "PRunSingleHisto::EstimatBkg(): Background " << start << ", " << end;
2341 if (end == start)
2342 end = fRunInfo->GetBkgRange(1);
2343 }
2344
2345 // check if start is within histogram bounds
2346 if (start >= fForward.size()) {
2347 std::cerr << std::endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
2348 std::cerr << std::endl << ">> histo lengths = " << fForward.size();
2349 std::cerr << std::endl << ">> background start = " << start;
2350 std::cerr << std::endl;
2351 return false;
2352 }
2353
2354 // check if end is within histogram bounds
2355 if (end >= fForward.size()) {
2356 std::cerr << std::endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
2357 std::cerr << std::endl << ">> histo lengths = " << fForward.size();
2358 std::cerr << std::endl << ">> background end = " << end;
2359 std::cerr << std::endl;
2360 return false;
2361 }
2362
2363 // calculate background
2364 Double_t bkg = 0.0;
2365
2366 // forward
2367 for (UInt_t i=start; i<end; i++)
2368 bkg += fForward[i];
2369 bkg /= static_cast<Double_t>(end - start + 1);
2370
2371 if (fScaleN0AndBkg)
2372 fBackground = bkg / (fTimeResolution * 1e3); // keep background (per 1 nsec) for chisq, max.log.likelihood, fTimeResolution us->ns
2373 else
2374 fBackground = bkg * fPacking; // keep background (per bin)
2375
2376 fRunInfo->SetBkgEstimated(fBackground, 0);
2377
2378 return true;
2379}
2380
2381//--------------------------------------------------------------------------
2382// IsScaleN0AndBkg (private)
2383//--------------------------------------------------------------------------
2427{
2428 Bool_t willScale = true;
2429
2430 PMsrLines *cmd = fMsrInfo->GetMsrCommands();
2431 for (UInt_t i=0; i<cmd->size(); i++) {
2432 if (cmd->at(i).fLine.Contains("SCALE_N0_BKG", TString::kIgnoreCase)) {
2433 TObjArray *tokens = nullptr;
2434 TObjString *ostr = nullptr;
2435 TString str;
2436 tokens = cmd->at(i).fLine.Tokenize(" \t");
2437 if (tokens->GetEntries() != 2) {
2438 std::cerr << std::endl << ">> PRunSingleHisto::IsScaleN0AndBkg(): **WARNING** Found uncorrect 'SCALE_N0_BKG' command, will ignore it.";
2439 std::cerr << std::endl << ">> Allowed commands: SCALE_N0_BKG TRUE | FALSE" << std::endl;
2440 return willScale;
2441 }
2442 ostr = dynamic_cast<TObjString*>(tokens->At(1));
2443 str = ostr->GetString();
2444 if (!str.CompareTo("FALSE", TString::kIgnoreCase)) {
2445 willScale = false;
2446 }
2447 // clean up
2448 if (tokens)
2449 delete tokens;
2450 }
2451 }
2452
2453 return willScale;
2454}
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 RRF_UNIT_MHz
Frequency in MHz (megahertz)
Definition PMusr.h:335
#define PMUSR_UNDEFINED
Definition PMusr.h:172
#define GAMMA_BAR_MUON
Definition PMusr.h:138
#define PMUON_LIFETIME
Definition PMusr.h:119
#define MSR_PARAM_FUN_OFFSET
Offset added to function indices for parameter parsing.
Definition PMusr.h:260
std::vector< PMsrPlotStructure > PMsrPlotList
Definition PMusr.h:1312
std::vector< PMsrLineStructure > PMsrLines
Definition PMusr.h:989
std::vector< PMsrParamStructure > PMsrParamList
Definition PMusr.h:1022
#define RRF_UNIT_Mcs
Angular frequency in Mc/s (Mega-cycles per second)
Definition PMusr.h:337
#define RRF_UNIT_G
Equivalent magnetic field in Gauss (G)
Definition PMusr.h:339
#define RRF_UNIT_kHz
Frequency in kHz (kilohertz)
Definition PMusr.h:333
#define ACCEL_PERIOD_PSI
PSI (Paul Scherrer Institute) accelerator cycle: 19.75 ns.
Definition PMusr.h:149
#define RRF_UNIT_T
Equivalent magnetic field in Tesla (T)
Definition PMusr.h:341
std::vector< Double_t > PDoubleVector
Definition PMusr.h:385
#define ACCEL_PERIOD_RAL
RAL (Rutherford Appleton Lab) - pulsed beam.
Definition PMusr.h:153
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 GetFitRangeOffset(UInt_t idx)
Definition PMusr.cpp:1157
virtual Double_t GetT0Bin(UInt_t idx=0)
Definition PMusr.cpp:1004
MSR file parser and manager for the musrfit framework.
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
virtual void CalculateKaiserFilterCoeff(Double_t wc, Double_t A, Double_t dw)
Calculates Kaiser window FIR filter coefficients for RRF smoothing.
Definition PRunBase.cpp:339
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
virtual void FilterTheo()
Applies Kaiser FIR filter to theory values for RRF fits.
Definition PRunBase.cpp:405
Raw data file reader and format converter for μSR data.
virtual Double_t CalcChiSquare(const std::vector< Double_t > &par)
Calculates χ² between histogram data and theory.
virtual void CalcNoOfFitBins()
Calculates start/end bin indices from fit time range.
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
Calculates expected χ² based on theory predictions.
virtual Bool_t IsScaleN0AndBkg()
Determines if N₀ and background should be scaled to 1/ns.
Double_t fBackground
Background level in counts/bin (estimated from pre-t0 bins or fixed value from RUN block)
UInt_t fNoOfFitBins
Number of bins within fit range (fStartTimeBin to fEndTimeBin)
virtual ~PRunSingleHisto()
Virtual destructor cleaning up allocated resources.
virtual Bool_t PrepareViewData(PRawRunData *runData, const UInt_t histoNo)
Prepares processed histogram data for viewing/plotting.
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
Calculates maximum likelihood for Poisson-distributed histogram counts.
virtual Bool_t EstimateBkg(UInt_t histoNo)
Estimates background from pre-t0 bins.
virtual void CalcTheory()
Evaluates theory function at all data points or high-resolution grid.
virtual UInt_t GetNoOfFitBins()
Returns the number of bins included in the fit range.
virtual void EstimateN0()
Estimates initial normalization N₀ from histogram data.
Int_t fEndTimeBin
Last bin index in fit range (exclusive: loop as i < fEndTimeBin)
PRunSingleHisto()
Default constructor creating an empty, invalid single histogram run object.
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
Determines fit range from MSR file settings.
virtual void SetFitRangeBin(const TString fitRange)
Sets fit range using bin-offset specification (COMMANDS block syntax).
Bool_t fScaleN0AndBkg
Scaling mode: true = scale N₀ and B to 1/ns, false = leave as 1/bin (determined by IsScaleN0AndBkg())
Int_t fStartTimeBin
First bin index in fit range (inclusive, 0-based after packing)
virtual Bool_t GetProperDataRange()
Determines data range (region of valid histogram data).
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &histoNo)
Determines and validates t0 values for histogram.
Bool_t fTheoAsData
Theory mode: true = at data points, false = high-resolution grid for smooth Fourier transforms.
Int_t fGoodBins[2]
Good bin markers for COMMANDS block: [0]=fgb (first good bin/t0), [1]=lgb (last good bin)
Int_t fPacking
Bin packing factor (REQUIRED: from RUN or GLOBAL block)
virtual Bool_t PrepareRawViewData(PRawRunData *runData, const UInt_t histoNo)
Prepares raw histogram data for viewing (minimal processing).
virtual Double_t CalcMaxLikelihoodExpected(const std::vector< Double_t > &par)
Calculates expected maximum likelihood.
virtual Bool_t PrepareData()
Main data preparation orchestrator.
virtual Bool_t PrepareFitData(PRawRunData *runData, const UInt_t histoNo)
Prepares histogram data for fitting.
PDoubleVector fForward
Forward detector histogram (background-corrected, packed)