musrfit  1.9.2
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-2023 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 //--------------------------------------------------------------------------
56 {
57  fScaleN0AndBkg = true;
58  fNoOfFitBins = 0;
59  fBackground = 0;
60  fPacking = -1;
61  fTheoAsData = false;
62 
63  // the 2 following variables are need in case fit range is given in bins, and since
64  // the fit range can be changed in the command block, these variables need to be accessible
65  fGoodBins[0] = -1;
66  fGoodBins[1] = -1;
67 
68  fStartTimeBin = -1;
69  fEndTimeBin = -1;
70 }
71 
72 //--------------------------------------------------------------------------
73 // Constructor
74 //--------------------------------------------------------------------------
83 PRunSingleHisto::PRunSingleHisto(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
84  PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
85 {
87  fNoOfFitBins = 0;
88  fBackground = 0;
89 
91  if (fPacking == -1) { // i.e. packing is NOT given in the RUN-block, it must be given in the GLOBAL-block
93  }
94  if (fPacking == -1) { // this should NOT happen, somethin is severely wrong
95  std::cerr << std::endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't find any packing information!";
96  std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
97  std::cerr << std::endl;
98  fValid = false;
99  return;
100  }
101 
102  // the 2 following variables are need in case fit range is given in bins, and since
103  // the fit range can be changed in the command block, these variables need to be accessible
104  fGoodBins[0] = -1;
105  fGoodBins[1] = -1;
106 
107  fStartTimeBin = -1;
108  fEndTimeBin = -1;
109 
110  if (!PrepareData()) {
111  std::cerr << std::endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't prepare data for fitting!";
112  std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
113  std::cerr << std::endl;
114  fValid = false;
115  }
116 }
117 
118 //--------------------------------------------------------------------------
119 // Destructor
120 //--------------------------------------------------------------------------
125 {
126  fForward.clear();
127 }
128 
129 //--------------------------------------------------------------------------
130 // CalcChiSquare (public)
131 //--------------------------------------------------------------------------
140 Double_t PRunSingleHisto::CalcChiSquare(const std::vector<Double_t>& par)
141 {
142  Double_t chisq = 0.0;
143  Double_t diff = 0.0;
144 
145  Double_t N0 = 0.0;
146 
147  // check if norm is a parameter or a function
148  if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
149  N0 = par[fRunInfo->GetNormParamNo()-1];
150  } else { // norm is a function
151  // get function number
152  UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
153  // evaluate function
154  N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
155  }
156 
157  // get tau
158  Double_t tau;
159  if (fRunInfo->GetLifetimeParamNo() != -1)
160  tau = par[fRunInfo->GetLifetimeParamNo()-1];
161  else
162  tau = PMUON_LIFETIME;
163 
164  // get background
165  Double_t bkg;
166  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
167  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
168  bkg = fBackground;
169  } else { // fixed bkg given
170  bkg = fRunInfo->GetBkgFix(0);
171  }
172  } else { // bkg fitted
173  bkg = par[fRunInfo->GetBkgFitParamNo()-1];
174  }
175 
176  // calculate functions
177  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
178  UInt_t funcNo = fMsrInfo->GetFuncNo(i);
179  fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
180  }
181 
182  // calculate chi square
183  Double_t time(1.0);
184  Int_t i;
185 
186  // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
187  // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
188  // for a given set of parameters---which should be done outside of the parallelized loop.
189  // For all other functions it means a tiny and acceptable overhead.
190  time = fTheory->Func(time, par, fFuncValues);
191 
192  #ifdef HAVE_GOMP
193  Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
194  if (chunk < 10)
195  chunk = 10;
196  #pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
197  #endif
198  for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
199  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
200  diff = fData.GetValue()->at(i) -
201  (N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg);
202  chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
203  }
204 
205  // the correction factor is need since the data scales like pack*t_res,
206  // whereas the error scales like sqrt(pack*t_res)
207  if (fScaleN0AndBkg)
208  chisq *= fPacking * (fTimeResolution * 1.0e3);
209 
210  return chisq;
211 }
212 
213 //--------------------------------------------------------------------------
214 // CalcChiSquareExpected (public)
215 //--------------------------------------------------------------------------
224 Double_t PRunSingleHisto::CalcChiSquareExpected(const std::vector<Double_t>& par)
225 {
226  Double_t chisq = 0.0;
227  Double_t diff = 0.0;
228  Double_t theo = 0.0;
229 
230  Double_t N0 = 0.0;
231 
232  // check if norm is a parameter or a function
233  if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
234  N0 = par[fRunInfo->GetNormParamNo()-1];
235  } else { // norm is a function
236  // get function number
237  UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
238  // evaluate function
239  N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
240  }
241 
242  // get tau
243  Double_t tau;
244  if (fRunInfo->GetLifetimeParamNo() != -1)
245  tau = par[fRunInfo->GetLifetimeParamNo()-1];
246  else
247  tau = PMUON_LIFETIME;
248 
249  // get background
250  Double_t bkg;
251  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
252  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
253  bkg = fBackground;
254  } else { // fixed bkg given
255  bkg = fRunInfo->GetBkgFix(0);
256  }
257  } else { // bkg fitted
258  bkg = par[fRunInfo->GetBkgFitParamNo()-1];
259  }
260 
261  // calculate functions
262  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
263  Int_t funcNo = fMsrInfo->GetFuncNo(i);
264  fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
265  }
266 
267  // calculate chi square
268  Double_t time(1.0);
269  Int_t i;
270 
271  // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
272  // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
273  // for a given set of parameters---which should be done outside of the parallelized loop.
274  // For all other functions it means a tiny and acceptable overhead.
275  time = fTheory->Func(time, par, fFuncValues);
276 
277  #ifdef HAVE_GOMP
278  Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
279  if (chunk < 10)
280  chunk = 10;
281  #pragma omp parallel for default(shared) private(i,time,theo,diff) schedule(dynamic,chunk) reduction(+:chisq)
282  #endif
283  for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
284  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
285  theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
286  diff = fData.GetValue()->at(i) - theo;
287  chisq += diff*diff / theo;
288  }
289 
290  // the correction factor is need since the data scales like pack*t_res,
291  // whereas the error scales like sqrt(pack*t_res)
292  if (fScaleN0AndBkg)
293  chisq *= fPacking * (fTimeResolution * 1.0e3);
294 
295  return chisq;
296 }
297 
298 //--------------------------------------------------------------------------
299 // CalcMaxLikelihood (public)
300 //--------------------------------------------------------------------------
309 Double_t PRunSingleHisto::CalcMaxLikelihood(const std::vector<Double_t>& par)
310 {
311  Double_t mllh = 0.0; // maximum log likelihood assuming poisson distribution for the single bin
312 
313  Double_t N0;
314 
315  // check if norm is a parameter or a function
316  if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
317  N0 = par[fRunInfo->GetNormParamNo()-1];
318  } else { // norm is a function
319  // get function number
321  // evaluate function
322  N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
323  }
324 
325  // get tau
326  Double_t tau;
327  if (fRunInfo->GetLifetimeParamNo() != -1)
328  tau = par[fRunInfo->GetLifetimeParamNo()-1];
329  else
330  tau = PMUON_LIFETIME;
331 
332  // get background
333  Double_t bkg;
334  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
335  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
336  bkg = fBackground;
337  } else { // fixed bkg given
338  bkg = fRunInfo->GetBkgFix(0);
339  }
340  } else { // bkg fitted
341  bkg = par[fRunInfo->GetBkgFitParamNo()-1];
342  }
343 
344  // calculate functions
345  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
346  UInt_t funcNo = fMsrInfo->GetFuncNo(i);
347  fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
348  }
349 
350  // calculate maximum log likelihood
351  Double_t theo;
352  Double_t data;
353  Double_t time(1.0);
354  Int_t i;
355 
356  // 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
357  Double_t normalizer = 1.0;
358 
359  if (fScaleN0AndBkg)
360  normalizer = fPacking * (fTimeResolution * 1.0e3);
361 
362  // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
363  // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
364  // for a given set of parameters---which should be done outside of the parallelized loop.
365  // For all other functions it means a tiny and acceptable overhead.
366  time = fTheory->Func(time, par, fFuncValues);
367 
368  #ifdef HAVE_GOMP
369  Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
370  if (chunk < 10)
371  chunk = 10;
372  #pragma omp parallel for default(shared) private(i,time,theo,data) schedule(dynamic,chunk) reduction(-:mllh)
373  #endif
374  for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
375  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
376  // calculate theory for the given parameter set
377  theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
378 
379  data = fData.GetValue()->at(i);
380 
381  if (theo <= 0.0) {
382  std::cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << std::endl;
383  continue;
384  }
385 
386  if (data > 1.0e-9) {
387  mllh += (theo-data) + data*log(data/theo);
388  } else {
389  mllh += (theo-data);
390  }
391  }
392 
393  return normalizer*2.0*mllh;
394 }
395 
396 //--------------------------------------------------------------------------
397 // CalcMaxLikelihoodExpected (public)
398 //--------------------------------------------------------------------------
407 Double_t PRunSingleHisto::CalcMaxLikelihoodExpected(const std::vector<Double_t>& par)
408 {
409  Double_t mllh = 0.0; // maximum log likelihood assuming poisson distribution for the single bin
410 
411  Double_t N0;
412 
413  // check if norm is a parameter or a function
414  if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
415  N0 = par[fRunInfo->GetNormParamNo()-1];
416  } else { // norm is a function
417  // get function number
419  // evaluate function
420  N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
421  }
422 
423  // get tau
424  Double_t tau;
425  if (fRunInfo->GetLifetimeParamNo() != -1)
426  tau = par[fRunInfo->GetLifetimeParamNo()-1];
427  else
428  tau = PMUON_LIFETIME;
429 
430  // get background
431  Double_t bkg;
432  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
433  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
434  bkg = fBackground;
435  } else { // fixed bkg given
436  bkg = fRunInfo->GetBkgFix(0);
437  }
438  } else { // bkg fitted
439  bkg = par[fRunInfo->GetBkgFitParamNo()-1];
440  }
441 
442  // calculate functions
443  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
444  UInt_t funcNo = fMsrInfo->GetFuncNo(i);
445  fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
446  }
447 
448  // calculate maximum log likelihood
449  Double_t theo;
450  Double_t data;
451  Double_t time(1.0);
452  Int_t i;
453 
454  // 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
455  Double_t normalizer = 1.0;
456 
457  if (fScaleN0AndBkg)
458  normalizer = fPacking * (fTimeResolution * 1.0e3);
459 
460  // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
461  // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
462  // for a given set of parameters---which should be done outside of the parallelized loop.
463  // For all other functions it means a tiny and acceptable overhead.
464  time = fTheory->Func(time, par, fFuncValues);
465 
466  #ifdef HAVE_GOMP
467  Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
468  if (chunk < 10)
469  chunk = 10;
470  #pragma omp parallel for default(shared) private(i,time,theo,data) schedule(dynamic,chunk) reduction(-:mllh)
471  #endif
472  for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
473  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
474  // calculate theory for the given parameter set
475  theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
476 
477  data = fData.GetValue()->at(i);
478 
479  if (theo <= 0.0) {
480  std::cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << std::endl;
481  continue;
482  }
483 
484  if (data > 1.0e-9) { // is this correct?? needs to be checked. See G-test
485  mllh += data*log(data/theo);
486  }
487  }
488 
489  return normalizer*2.0*mllh;
490 }
491 
492 //--------------------------------------------------------------------------
493 // CalcTheory (public)
494 //--------------------------------------------------------------------------
499 {
500  // feed the parameter vector
501  std::vector<Double_t> par;
502  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
503  for (UInt_t i=0; i<paramList->size(); i++)
504  par.push_back((*paramList)[i].fValue);
505 
506  // calculate asymmetry
507  Double_t N0;
508  // check if norm is a parameter or a function
509  if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
510  N0 = par[fRunInfo->GetNormParamNo()-1];
511  } else { // norm is a function
512  // get function number
514  // evaluate function
515  N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
516  }
517 
518  // get tau
519  Double_t tau;
520  if (fRunInfo->GetLifetimeParamNo() != -1)
521  tau = par[fRunInfo->GetLifetimeParamNo()-1];
522  else
523  tau = PMUON_LIFETIME;
524 
525  // get background
526  Double_t bkg;
527  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
528  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
529  bkg = fBackground;
530  } else { // fixed bkg given
531  bkg = fRunInfo->GetBkgFix(0);
532  }
533  } else { // bkg fitted
534  bkg = par[fRunInfo->GetBkgFitParamNo()-1];
535  }
536 
537  // calculate functions
538  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
540  }
541 
542  // calculate theory
543  UInt_t size = fData.GetValue()->size();
544  Double_t start = fData.GetDataTimeStart();
545  Double_t resolution = fData.GetDataTimeStep();
546  Double_t time;
547  for (UInt_t i=0; i<size; i++) {
548  time = start + static_cast<Double_t>(i)*resolution;
549  fData.AppendTheoryValue(N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg);
550  }
551 
552  // clean up
553  par.clear();
554 }
555 
556 //--------------------------------------------------------------------------
557 // GetNoOfFitBins (public)
558 //--------------------------------------------------------------------------
565 {
566  CalcNoOfFitBins();
567 
568  return fNoOfFitBins;
569 }
570 
571 //--------------------------------------------------------------------------
572 // SetFitRangeBin (public)
573 //--------------------------------------------------------------------------
585 void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
586 {
587  TObjArray *tok = nullptr;
588  TObjString *ostr = nullptr;
589  TString str;
590  Ssiz_t idx = -1;
591  Int_t offset = 0;
592 
593  tok = fitRange.Tokenize(" \t");
594 
595  if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
596  // handle fgb+n0 entry
597  ostr = dynamic_cast<TObjString*>(tok->At(1));
598  str = ostr->GetString();
599  // check if there is an offset present
600  idx = str.First("+");
601  if (idx != -1) { // offset present
602  str.Remove(0, idx+1);
603  if (str.IsFloat()) // if str is a valid number, convert is to an integer
604  offset = str.Atoi();
605  }
606  fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
607 
608  // handle lgb-n1 entry
609  ostr = dynamic_cast<TObjString*>(tok->At(2));
610  str = ostr->GetString();
611  // check if there is an offset present
612  idx = str.First("-");
613  if (idx != -1) { // offset present
614  str.Remove(0, idx+1);
615  if (str.IsFloat()) // if str is a valid number, convert is to an integer
616  offset = str.Atoi();
617  }
618  fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
619  } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
620  Int_t pos = 2*(fRunNo+1)-1;
621 
622  if (pos + 1 >= tok->GetEntries()) {
623  std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
624  std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
625  } else {
626  // handle fgb+n0 entry
627  ostr = dynamic_cast<TObjString*>(tok->At(pos));
628  str = ostr->GetString();
629  // check if there is an offset present
630  idx = str.First("+");
631  if (idx != -1) { // offset present
632  str.Remove(0, idx+1);
633  if (str.IsFloat()) // if str is a valid number, convert is to an integer
634  offset = str.Atoi();
635  }
636  fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
637 
638  // handle lgb-n1 entry
639  ostr = dynamic_cast<TObjString*>(tok->At(pos+1));
640  str = ostr->GetString();
641  // check if there is an offset present
642  idx = str.First("-");
643  if (idx != -1) { // offset present
644  str.Remove(0, idx+1);
645  if (str.IsFloat()) // if str is a valid number, convert is to an integer
646  offset = str.Atoi();
647  }
648  fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
649  }
650  } else { // error
651  std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
652  std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
653  }
654 
655  // clean up
656  if (tok) {
657  delete tok;
658  }
659 }
660 
661 //--------------------------------------------------------------------------
662 // CalcNoOfFitBins (public)
663 //--------------------------------------------------------------------------
668 {
669  // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
670  fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
671  if (fStartTimeBin < 0)
672  fStartTimeBin = 0;
673  fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
674  if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
675  fEndTimeBin = fData.GetValue()->size();
676 
679  else
680  fNoOfFitBins = 0;
681 }
682 
683 //--------------------------------------------------------------------------
684 // PrepareData (protected)
685 //--------------------------------------------------------------------------
700 {
701  Bool_t success = true;
702 
703  if (!fValid)
704  return false;
705 
706  // keep the Global block info
707  PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
708 
709  // get the proper run
711  if (!runData) { // couldn't get run
712  std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
713  std::cerr << std::endl;
714  return false;
715  }
716 
717  // keep the field from the meta-data from the data-file
718  fMetaData.fField = runData->GetField();
719 
720  // keep the energy from the meta-data from the data-file
721  fMetaData.fEnergy = runData->GetEnergy();
722 
723  // keep the temperature(s) from the meta-data from the data-file
724  for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
725  fMetaData.fTemp.push_back(runData->GetTemperature(i));
726 
727  // collect histogram numbers
728  PUIntVector histoNo; // histoNo = msr-file forward + redGreen_offset - 1
729  for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
730  histoNo.push_back(fRunInfo->GetForwardHistoNo(i));
731 
732  if (!runData->IsPresent(histoNo[i])) {
733  std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **PANIC ERROR**:";
734  std::cerr << std::endl << ">> histoNo found = " << histoNo[i] << ", which is NOT present in the data file!?!?";
735  std::cerr << std::endl << ">> Will quit :-(";
736  std::cerr << std::endl;
737  histoNo.clear();
738  return false;
739  }
740  }
741 
742  // keep the time resolution in (us)
743  fTimeResolution = runData->GetTimeResolution()/1.0e3;
744  std::cout.precision(10);
745  std::cout << std::endl << ">> PRunSingleHisto::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
746 
747  // get all the proper t0's and addt0's for the current RUN block
748  if (!GetProperT0(runData, globalBlock, histoNo)) {
749  return false;
750  }
751 
752  // keep the histo of each group at this point (addruns handled below)
753  std::vector<PDoubleVector> forward;
754  forward.resize(histoNo.size()); // resize to number of groups
755  for (UInt_t i=0; i<histoNo.size(); i++) {
756  forward[i].resize(runData->GetDataBin(histoNo[i])->size());
757  forward[i] = *runData->GetDataBin(histoNo[i]);
758  }
759 
760  // check if there are runs to be added to the current one
761  if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
762  PRawRunData *addRunData;
763  for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
764 
765  // get run to be added to the main one
766  addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
767  if (addRunData == nullptr) { // couldn't get run
768  std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
769  std::cerr << std::endl;
770  return false;
771  }
772 
773  // add forward run
774  UInt_t addRunSize;
775  for (UInt_t k=0; k<histoNo.size(); k++) { // fill each group
776  addRunSize = addRunData->GetDataBin(histoNo[k])->size();
777  for (UInt_t j=0; j<addRunData->GetDataBin(histoNo[k])->size(); j++) { // loop over the bin indices
778  // make sure that the index stays in the proper range
779  if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) >= 0) &&
780  (j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) < addRunSize)) {
781  forward[k][j] += addRunData->GetDataBin(histoNo[k])->at(j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]));
782  }
783  }
784  }
785  }
786  }
787 
788  // set forward histo data of the first group
789  fForward.resize(forward[0].size());
790  for (UInt_t i=0; i<fForward.size(); i++) {
791  fForward[i] = forward[0][i];
792  }
793 
794  // group histograms, add all the remaining forward histograms of the group
795  for (UInt_t i=1; i<histoNo.size(); i++) { // loop over the groupings
796  for (UInt_t j=0; j<runData->GetDataBin(histoNo[i])->size(); j++) { // loop over the bin indices
797  // make sure that the index stays within proper range
798  if ((static_cast<Int_t>(j)+fT0s[i]-fT0s[0] >= 0) && (j+fT0s[i]-fT0s[0] < runData->GetDataBin(histoNo[i])->size())) {
799  fForward[j] += forward[i][j+static_cast<Int_t>(fT0s[i])-static_cast<Int_t>(fT0s[0])];
800  }
801  }
802  }
803 
804  // get the data range (fgb/lgb) for the current RUN block
805  if (!GetProperDataRange()) {
806  return false;
807  }
808 
809  // get the fit range for the current RUN block
810  GetProperFitRange(globalBlock);
811 
812  // get the lifetimecorrection flag
813  Bool_t lifetimecorrection = false;
815  lifetimecorrection = plot->at(0).fLifeTimeCorrection;
816 
817  // do the more fit/view specific stuff
818  if (fHandleTag == kFit)
819  success = PrepareFitData(runData, histoNo[0]);
820  else if ((fHandleTag == kView) && !lifetimecorrection)
821  success = PrepareRawViewData(runData, histoNo[0]);
822  else if ((fHandleTag == kView) && lifetimecorrection)
823  success = PrepareViewData(runData, histoNo[0]);
824  else
825  success = false;
826 
827  // cleanup
828  histoNo.clear();
829 
830  return success;
831 }
832 
833 //--------------------------------------------------------------------------
834 // PrepareFitData (protected)
835 //--------------------------------------------------------------------------
851 Bool_t PRunSingleHisto::PrepareFitData(PRawRunData* runData, const UInt_t histoNo)
852 {
853  if (fMsrInfo->EstimateN0()) {
854  EstimateN0();
855  }
856 
857  // transform raw histo data. This is done the following way (for details see the manual):
858  // for the single histo fit, just the rebinned raw data are copied
859 
860  // check how the background shall be handled
861  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg shall **NOT** be fitted
862  // subtract background from histogramms ------------------------------------------
863  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
864  if (fRunInfo->GetBkgRange(0) >= 0) {
865  if (!EstimateBkg(histoNo))
866  return false;
867  } else { // no background given to do the job, try estimate
868  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
869  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
870  std::cerr << std::endl << ">> PRunSingleHisto::PrepareFitData(): **WARNING** Neither fix background nor background bins are given!";
871  std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
872  std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
873  std::cerr << std::endl;
874  if (!EstimateBkg(histoNo))
875  return false;
876  }
877  } else { // fixed background given
878  for (UInt_t i=0; i<fForward.size(); i++) {
879  fForward[i] -= fRunInfo->GetBkgFix(0);
880  }
881  }
882  }
883 
884  // everything looks fine, hence fill data set
885  Int_t t0 = static_cast<Int_t>(fT0s[0]);
886  Double_t value = 0.0;
887  Double_t normalizer = 1.0;
888  // in order that after rebinning the fit does not need to be redone (important for plots)
889  // the value is normalize to per 1 nsec if scaling is whished
890  if (fScaleN0AndBkg)
891  normalizer = fPacking * (fTimeResolution * 1.0e3); // fTimeResolution us->ns
892  // data start at data_start-t0
893  // time shifted so that packing is included correctly, i.e. t0 == t0 after packing
894  fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(fGoodBins[0])-0.5) + static_cast<Double_t>(fPacking)/2.0 - static_cast<Double_t>(t0)));
896  for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
897  if (fPacking == 1) {
898  value = fForward[i];
899  value /= normalizer;
900  fData.AppendValue(value);
901  if (value == 0.0)
902  fData.AppendErrorValue(1.0/normalizer);
903  else
904  fData.AppendErrorValue(TMath::Sqrt(value));
905  } else { // packed data, i.e. fPacking > 1
906  if (((i-fGoodBins[0]) % fPacking == 0) && (i != fGoodBins[0])) { // fill data
907  value /= normalizer;
908  fData.AppendValue(value);
909  if (value == 0.0)
910  fData.AppendErrorValue(1.0/normalizer);
911  else
912  fData.AppendErrorValue(TMath::Sqrt(value));
913  // reset values
914  value = 0.0;
915  }
916  value += fForward[i];
917  }
918  }
919 
920  CalcNoOfFitBins();
921 
922  return true;
923 }
924 
925 //--------------------------------------------------------------------------
926 // PrepareRawViewData (protected)
927 //--------------------------------------------------------------------------
944 Bool_t PRunSingleHisto::PrepareRawViewData(PRawRunData* runData, const UInt_t histoNo)
945 {
946  // check if view_packing is wished
947  Int_t packing = fPacking;
948  if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
949  packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
950  }
951 
952  // calculate necessary norms
953  Double_t dataNorm = 1.0, theoryNorm = 1.0;
954  if (fScaleN0AndBkg) {
955  dataNorm = 1.0/ (packing * (fTimeResolution * 1.0e3)); // fTimeResolution us->ns
956  } else if (!fScaleN0AndBkg && (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0)) {
957  theoryNorm = static_cast<Double_t>(fMsrInfo->GetMsrPlotList()->at(0).fViewPacking)/static_cast<Double_t>(fPacking);
958  }
959 
960  // raw data, since PMusrCanvas is doing ranging etc.
961  // start = the first bin which is a multiple of packing backward from first good data bin
962  Int_t start = fGoodBins[0] - (fGoodBins[0]/packing)*packing;
963  // end = last bin starting from start which is a multiple of packing and still within the data
964  Int_t end = start + ((fForward.size()-start)/packing)*packing;
965  // check if data range has been provided, and if not try to estimate them
966  if (start < 0) {
967  Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
968  start = (static_cast<Int_t>(fT0s[0])+offset) - ((static_cast<Int_t>(fT0s[0])+offset)/packing)*packing;
969  end = start + ((fForward.size()-start)/packing)*packing;
970  std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
971  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
972  std::cerr << std::endl;
973  }
974  // check if start, end, and t0 make any sense
975  // 1st check if start and end are in proper order
976  if (end < start) { // need to swap them
977  Int_t keep = end;
978  end = start;
979  start = keep;
980  }
981  // 2nd check if start is within proper bounds
982  if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
983  std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** start data bin doesn't make any sense!";
984  std::cerr << std::endl;
985  return false;
986  }
987  // 3rd check if end is within proper bounds
988  if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
989  std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** end data bin doesn't make any sense!";
990  std::cerr << std::endl;
991  return false;
992  }
993 
994  // everything looks fine, hence fill data set
995  Int_t t0 = static_cast<Int_t>(fT0s[0]);
996  Double_t value = 0.0;
997  // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
998  fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start)-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0)));
1000 
1001  for (Int_t i=start; i<end; i++) {
1002  if (((i-start) % packing == 0) && (i != start)) { // fill data
1003  value *= dataNorm;
1004  fData.AppendValue(value);
1005  if (value == 0.0)
1006  fData.AppendErrorValue(1.0);
1007  else
1008  fData.AppendErrorValue(TMath::Sqrt(value*dataNorm));
1009  // reset values
1010  value = 0.0;
1011  }
1012  value += fForward[i];
1013  }
1014 
1015  CalcNoOfFitBins();
1016 
1017  // fill theory vector for kView
1018  // feed the parameter vector
1019  std::vector<Double_t> par;
1020  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1021  for (UInt_t i=0; i<paramList->size(); i++)
1022  par.push_back((*paramList)[i].fValue);
1023 
1024  // calculate asymmetry
1025  Double_t N0;
1026  // check if norm is a parameter or a function
1027  if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
1028  N0 = par[fRunInfo->GetNormParamNo()-1];
1029  } else { // norm is a function
1030  // get function number
1031  UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
1032  // evaluate function
1033  N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1034  }
1035  N0 *= theoryNorm;
1036 
1037  // get tau
1038  Double_t tau;
1039  if (fRunInfo->GetLifetimeParamNo() != -1)
1040  tau = par[fRunInfo->GetLifetimeParamNo()-1];
1041  else
1042  tau = PMUON_LIFETIME;
1043 
1044  // get background
1045  Double_t bkg;
1046  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
1047  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
1048  if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
1049  if (!EstimateBkg(histoNo))
1050  return false;
1051  } else { // no background given to do the job, try estimate
1052  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
1053  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
1054  std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** Neither fix background nor background bins are given!";
1055  std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
1056  std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
1057  std::cerr << std::endl;
1058  if (!EstimateBkg(histoNo))
1059  return false;
1060  }
1061  bkg = fBackground;
1062  } else { // fixed bkg given
1063  bkg = fRunInfo->GetBkgFix(0);
1064  }
1065  } else { // bkg fitted
1066  bkg = par[fRunInfo->GetBkgFitParamNo()-1];
1067  }
1068  bkg *= theoryNorm;
1069 
1070  // calculate functions
1071  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1073  }
1074 
1075  // calculate theory
1076  UInt_t size = fForward.size();
1077  Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1079  if (fTheoAsData) { // calculate theory only at the data points
1081  } else {
1082  // finer binning for the theory (8 times as many points = factor)
1083  size *= factor;
1084  fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1085  }
1086 
1087  Double_t time;
1088  Double_t theoryValue;
1089  for (UInt_t i=0; i<size; i++) {
1091  theoryValue = fTheory->Func(time, par, fFuncValues);
1092  if (fabs(theoryValue) > 1.0e10) { // dirty hack needs to be fixed!!
1093  theoryValue = 0.0;
1094  }
1095  fData.AppendTheoryValue(N0*TMath::Exp(-time/tau)*(1.0+theoryValue)+bkg);
1096  }
1097 
1098  // clean up
1099  par.clear();
1100 
1101  return true;
1102 }
1103 
1104 //--------------------------------------------------------------------------
1105 // PrepareViewData (protected)
1106 //--------------------------------------------------------------------------
1133 Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histoNo)
1134 {
1135  // check if view_packing is wished. This is a global option for all PLOT blocks!
1136  Int_t packing = fPacking;
1137  if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
1138  packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1139  }
1140  // check if rrf_packing is present. This is a global option for all PLOT blocks, since operated on a single set of data.
1141  if (fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking > 0) {
1142  packing = fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking;
1143  }
1144 
1145  // calculate necessary norms
1146  Double_t dataNorm = 1.0, theoryNorm = 1.0;
1147  if (fScaleN0AndBkg) {
1148  dataNorm = 1.0/ (packing * (fTimeResolution * 1.0e3)); // fTimeResolution us->ns
1149  } else if (!fScaleN0AndBkg && (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0)) {
1150  theoryNorm = static_cast<Double_t>(fMsrInfo->GetMsrPlotList()->at(0).fViewPacking)/static_cast<Double_t>(fPacking);
1151  }
1152 
1153  // transform raw histo data. This is done the following way (for details see the manual):
1154  // for the single histo fit, just the rebinned raw data are copied
1155  // first get start data, end data, and t0
1156  Int_t t0 = static_cast<Int_t>(fT0s[0]);
1157 
1158  // start = the first bin which is a multiple of packing backward from first good data bin
1159  Int_t start = fGoodBins[0] - (fGoodBins[0]/packing)*packing;
1160  // end = last bin starting from start which is a multiple of packing and still within the data
1161  Int_t end = start + ((fForward.size()-start)/packing)*packing;
1162 
1163  // check if data range has been provided, and if not try to estimate them
1164  if (start < 0) {
1165  Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
1166  start = (static_cast<Int_t>(fT0s[0])+offset) - ((static_cast<Int_t>(fT0s[0])+offset)/packing)*packing;
1167  end = start + ((fForward.size()-start)/packing)*packing;
1168  std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
1169  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1170  std::cerr << std::endl;
1171  }
1172 
1173  // check if start, end, and t0 make any sense
1174  // 1st check if start and end are in proper order
1175  if (end < start) { // need to swap them
1176  Int_t keep = end;
1177  end = start;
1178  start = keep;
1179  }
1180  // 2nd check if start is within proper bounds
1181  if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
1182  std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
1183  std::cerr << std::endl;
1184  return false;
1185  }
1186  // 3rd check if end is within proper bounds
1187  if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
1188  std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
1189  std::cerr << std::endl;
1190  return false;
1191  }
1192 
1193  // everything looks fine, hence fill data set
1194 
1195  // feed the parameter vector
1196  std::vector<Double_t> par;
1197  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1198  for (UInt_t i=0; i<paramList->size(); i++)
1199  par.push_back((*paramList)[i].fValue);
1200 
1201  // calculate asymmetry
1202  Double_t N0;
1203  // check if norm is a parameter or a function
1204  if (fRunInfo->GetNormParamNo() < MSR_PARAM_FUN_OFFSET) { // norm is a parameter
1205  N0 = par[fRunInfo->GetNormParamNo()-1];
1206  } else { // norm is a function
1207  // get function number
1208  UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
1209  // evaluate function
1210  N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1211  }
1212  N0 *= theoryNorm;
1213 
1214  // get tau
1215  Double_t tau;
1216  if (fRunInfo->GetLifetimeParamNo() != -1)
1217  tau = par[fRunInfo->GetLifetimeParamNo()-1];
1218  else
1219  tau = PMUON_LIFETIME;
1220 
1221  // get background
1222  Double_t bkg;
1223  if (fRunInfo->GetBkgFitParamNo() == -1) { // bkg not fitted
1224  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given (background interval)
1225  if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
1226  if (!EstimateBkg(histoNo))
1227  return false;
1228  } else { // no background given to do the job, try estimate
1229  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
1230  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
1231  std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** Neither fix background nor background bins are given!";
1232  std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
1233  std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
1234  std::cerr << std::endl;
1235  if (!EstimateBkg(histoNo))
1236  return false;
1237  }
1238  bkg = fBackground;
1239  } else { // fixed bkg given
1240  bkg = fRunInfo->GetBkgFix(0);
1241  }
1242  } else { // bkg fitted
1243  bkg = par[fRunInfo->GetBkgFitParamNo()-1];
1244  }
1245  bkg *= theoryNorm;
1246 
1247  Double_t value = 0.0;
1248  Double_t expval = 0.0;
1249  Double_t rrf_val = 0.0;
1250  Double_t time = 0.0;
1251 
1252  // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1253  fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start)-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0)));
1255 
1256  // data is always normalized to (per nsec!!)
1257  Double_t gammaRRF = 0.0, wRRF = 0.0, phaseRRF = 0.0;
1258  if (fMsrInfo->GetMsrPlotList()->at(0).fRRFFreq == 0.0) { // normal Data representation
1259  for (Int_t i=start; i<end; i++) {
1260  if (((i-start) % packing == 0) && (i != start)) { // fill data
1261  value *= dataNorm;
1262  // since the packing counter is already at the end of the bin, the time needs be shifted back by pack*time_resolution
1263  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;
1264  expval = TMath::Exp(+time/tau)/N0;
1265  fData.AppendValue(-1.0+expval*(value-bkg));
1266  fData.AppendErrorValue(expval*TMath::Sqrt(value*dataNorm));
1267  value = 0.0;
1268  }
1269  value += fForward[i];
1270  }
1271  } else { // RRF representation
1272  // check which units shall be used
1273  switch (fMsrInfo->GetMsrPlotList()->at(0).fRRFUnit) {
1274  case RRF_UNIT_kHz:
1275  gammaRRF = TMath::TwoPi()*1.0e-3;
1276  break;
1277  case RRF_UNIT_MHz:
1278  gammaRRF = TMath::TwoPi();
1279  break;
1280  case RRF_UNIT_Mcs:
1281  gammaRRF = 1.0;
1282  break;
1283  case RRF_UNIT_G:
1284  gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi();
1285  break;
1286  case RRF_UNIT_T:
1287  gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi()*1.0e4;
1288  break;
1289  default:
1290  gammaRRF = TMath::TwoPi();
1291  break;
1292  }
1293  wRRF = gammaRRF * fMsrInfo->GetMsrPlotList()->at(0).fRRFFreq;
1294  phaseRRF = fMsrInfo->GetMsrPlotList()->at(0).fRRFPhase / 180.0 * TMath::Pi();
1295 
1296  Double_t error = 0.0;
1297  for (Int_t i=start; i<end; i++) {
1298  if (((i-start) % packing == 0) && (i != start)) { // fill data
1299  fData.AppendValue(2.0*value/packing); // factor 2 needed because cos(a)cos(b) = 1/2(cos(a+b)+cos(a-b))
1300  fData.AppendErrorValue(expval*TMath::Sqrt(error/packing));
1301  value = 0.0;
1302  error = 0.0;
1303  }
1304  time = (static_cast<Double_t>(i)-t0)*fTimeResolution;
1305  expval = TMath::Exp(+time/tau)/N0;
1306  rrf_val = (-1.0+expval*(fForward[i]/(fTimeResolution*1.0e3)-bkg))*TMath::Cos(wRRF * time + phaseRRF);
1307  value += rrf_val;
1308  error += fForward[i]*dataNorm;
1309  }
1310  }
1311 
1312  CalcNoOfFitBins();
1313 
1314  // calculate functions
1315  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1317  }
1318 
1319  // calculate theory
1320  Double_t theoryValue;
1321  UInt_t size = fForward.size()/packing;
1322  const Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1323  UInt_t rebinRRF = 0;
1324 
1325  if (wRRF == 0) { // no RRF
1327  if (fTheoAsData) { // calculate theory only at the data points
1329  } else {
1330  // finer binning for the theory (8 times as many points = factor)
1331  size *= factor;
1332  fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1333  }
1334  } else { // RRF
1335  rebinRRF = static_cast<UInt_t>((TMath::Pi()/2.0/wRRF)/fTimeResolution); // RRF time resolution / data time resolution
1337  fData.SetTheoryTimeStep(TMath::Pi()/2.0/wRRF/rebinRRF); // = theory time resolution as close as possible to the data time resolution compatible with wRRF
1338  }
1339 
1340  for (UInt_t i=0; i<size; i++) {
1341  time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1342  theoryValue = fTheory->Func(time, par, fFuncValues);
1343  if (wRRF != 0.0) {
1344  theoryValue *= 2.0*TMath::Cos(wRRF * time + phaseRRF);
1345  }
1346  if (fabs(theoryValue) > 10.0) { // dirty hack needs to be fixed!!
1347  theoryValue = 0.0;
1348  }
1349  fData.AppendTheoryValue(theoryValue);
1350  }
1351 
1352  // if RRF filter the theory with a FIR Kaiser low pass filter
1353  if (wRRF != 0.0) {
1354  // rebin theory to the RRF frequency
1355  if (rebinRRF != 0) {
1356  Double_t dval = 0.0;
1357  PDoubleVector theo;
1358  for (UInt_t i=0; i<fData.GetTheory()->size(); i++) {
1359  if ((i % rebinRRF == 0) && (i != 0)) {
1360  theo.push_back(dval/rebinRRF);
1361  dval = 0.0;
1362  }
1363  dval += fData.GetTheory()->at(i);
1364  }
1365  fData.SetTheoryTimeStart(fData.GetTheoryTimeStart()+static_cast<Double_t>(rebinRRF-1)*fData.GetTheoryTimeStep()/2.0);
1367  fData.ReplaceTheory(theo);
1368  theo.clear();
1369  }
1370 
1371  // filter theory
1372  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)
1373  FilterTheo();
1374  }
1375 
1376  // clean up
1377  par.clear();
1378 
1379  return true;
1380 }
1381 
1382 //--------------------------------------------------------------------------
1383 // GetProperT0 (private)
1384 //--------------------------------------------------------------------------
1403 {
1404  // feed all T0's
1405  // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1406  fT0s.clear();
1407  fT0s.resize(histoNo.size());
1408  for (UInt_t i=0; i<fT0s.size(); i++) {
1409  fT0s[i] = -1.0;
1410  }
1411 
1412  // fill in the T0's from the msr-file (if present)
1413  for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1414  fT0s[i] = fRunInfo->GetT0Bin(i);
1415  }
1416 
1417  // fill in the T0's from the GLOBAL block section (if present)
1418  for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1419  if (fT0s[i] == -1.0) { // i.e. not given in the RUN block section
1420  fT0s[i] = globalBlock->GetT0Bin(i);
1421  }
1422  }
1423 
1424  // fill in the T0's from the data file, if not already present in the msr-file
1425  for (UInt_t i=0; i<histoNo.size(); i++) {
1426  if (fT0s[i] == -1.0) { // i.e. not present in the msr-file, try the data file
1427  if (runData->GetT0Bin(histoNo[i]) > 0.0) {
1428  fT0s[i] = runData->GetT0Bin(histoNo[i]);
1429  fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
1430  }
1431  }
1432  }
1433 
1434  // 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
1435  for (UInt_t i=0; i<histoNo.size(); i++) {
1436  if (fT0s[i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1437  fT0s[i] = runData->GetT0BinEstimated(histoNo[i]);
1438  fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
1439 
1440  std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1441  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1442  std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(histoNo[i]);
1443  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1444  std::cerr << std::endl;
1445  }
1446  }
1447 
1448  // check if t0 is within proper bounds
1449  for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
1450  if ((fT0s[i] < 0.0) || (fT0s[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1451  std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** t0 data bin (" << fT0s[i] << ") doesn't make any sense!";
1452  std::cerr << std::endl;
1453  return false;
1454  }
1455  }
1456 
1457  // check if there are runs to be added to the current one. If yes keep the needed t0's
1458  if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1459  PRawRunData *addRunData;
1460  fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1461  for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1462 
1463  // get run to be added to the main one
1464  addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
1465  if (addRunData == nullptr) { // couldn't get run
1466  std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1467  std::cerr << std::endl;
1468  return false;
1469  }
1470 
1471  // feed all T0's
1472  // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1473  fAddT0s[i-1].resize(histoNo.size());
1474  for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1475  fAddT0s[i-1][j] = -1.0;
1476  }
1477 
1478  // fill in the T0's from the msr-file (if present)
1479  for (UInt_t j=0; j<fRunInfo->GetT0BinSize(); j++) {
1480  fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1,j); // addRunIdx starts at 0
1481  }
1482 
1483  // fill in the T0's from the data file, if not already present in the msr-file
1484  for (UInt_t j=0; j<histoNo.size(); j++) {
1485  if (fAddT0s[i-1][j] == -1.0) // i.e. not present in the msr-file, try the data file
1486  if (addRunData->GetT0Bin(histoNo[j]) > 0.0) {
1487  fAddT0s[i-1][j] = addRunData->GetT0Bin(histoNo[j]);
1488  fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1489  }
1490  }
1491 
1492  // 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
1493  for (UInt_t j=0; j<histoNo.size(); j++) {
1494  if (fAddT0s[i-1][j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1495  fAddT0s[i-1][j] = addRunData->GetT0BinEstimated(histoNo[j]);
1496  fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1497 
1498  std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1499  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1500  std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(histoNo[j]);
1501  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1502  std::cerr << std::endl;
1503  }
1504  }
1505 
1506  // check if t0 is within proper bounds
1507  for (UInt_t j=0; j<fRunInfo->GetForwardHistoNoSize(); j++) {
1508  if ((fAddT0s[i-1][j] < 0.0) || (fAddT0s[i-1][j] > static_cast<Int_t>(addRunData->GetDataBin(histoNo[j])->size()))) {
1509  std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** addt0 data bin (" << fAddT0s[i-1][j] << ") doesn't make any sense!";
1510  std::cerr << std::endl;
1511  return false;
1512  }
1513  }
1514  }
1515  }
1516 
1517  return true;
1518 }
1519 
1520 //--------------------------------------------------------------------------
1521 // GetProperDataRange (private)
1522 //--------------------------------------------------------------------------
1534 {
1535  // get start/end data
1536  Int_t start;
1537  Int_t end;
1538  start = fRunInfo->GetDataRange(0);
1539  end = fRunInfo->GetDataRange(1);
1540 
1541  // check if data range has been given in the RUN block, if not try to get it from the GLOBAL block
1542  if (start < 0) {
1543  start = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1544  }
1545  if (end < 0) {
1546  end = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1547  }
1548 
1549  // check if data range has been provided, and if not try to estimate them
1550  if (start < 0) {
1551  Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
1552  start = static_cast<Int_t>(fT0s[0])+offset;
1553  fRunInfo->SetDataRange(start, 0);
1554  std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start << ".";
1555  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1556  std::cerr << std::endl;
1557  }
1558  if (end < 0) {
1559  end = fForward.size();
1560  fRunInfo->SetDataRange(end, 1);
1561  std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range end = " << end << ".";
1562  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1563  std::cerr << std::endl;
1564  }
1565 
1566  // check if start and end make any sense
1567  // 1st check if start and end are in proper order
1568  if (end < start) { // need to swap them
1569  Int_t keep = end;
1570  end = start;
1571  start = keep;
1572  }
1573  // 2nd check if start is within proper bounds
1574  if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
1575  std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** start data bin (" << start << ") doesn't make any sense!";
1576  std::cerr << std::endl;
1577  return false;
1578  }
1579  // 3rd check if end is within proper bounds
1580  if (end < 0) {
1581  std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** end data bin (" << end << ") doesn't make any sense!";
1582  std::cerr << std::endl;
1583  return false;
1584  }
1585  if (end > static_cast<Int_t>(fForward.size())) {
1586  std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** end data bin (" << end << ") > histo length (" << fForward.size() << ").";
1587  std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
1588  std::cerr << std::endl;
1589  end = static_cast<Int_t>(fForward.size())-1;
1590  }
1591 
1592  // keep good bins for potential later use
1593  fGoodBins[0] = start;
1594  fGoodBins[1] = end;
1595 
1596  return true;
1597 }
1598 
1599 //--------------------------------------------------------------------------
1600 // GetProperFitRange (private)
1601 //--------------------------------------------------------------------------
1615 {
1616  // set fit start/end time; first check RUN Block
1619  // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1620  if (fRunInfo->IsFitRangeInBin()) {
1621  fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1622  fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1623  // write these times back into the data structure. This way it is available when writting the log-file
1626  }
1627  if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
1628  fFitStartTime = globalBlock->GetFitRange(0);
1629  fFitEndTime = globalBlock->GetFitRange(1);
1630  // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1631  if (globalBlock->IsFitRangeInBin()) {
1632  fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1633  fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1634  // write these times back into the data structure. This way it is available when writting the log-file
1635  globalBlock->SetFitRange(fFitStartTime, 0);
1636  globalBlock->SetFitRange(fFitEndTime, 1);
1637  }
1638  }
1640  fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
1641  fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
1642  std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
1643  std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
1644  }
1645 }
1646 
1647 //--------------------------------------------------------------------------
1648 // EstimateN0 (private)
1649 //--------------------------------------------------------------------------
1654 {
1655  // check that 'norm' in the msr-file run block is indeed a parameter number.
1656  // in case it is a function, nothing will be done.
1657  UInt_t paramNo = fRunInfo->GetNormParamNo();
1658  if (paramNo > 10000) // i.e. fun or map
1659  return;
1660 
1661  // get the parameters
1663  assert(param);
1664 
1665  if (paramNo > param->size()) {
1666  std::cerr << std::endl << ">> PRunSingleHisto::EstimateN0: **ERROR** found parameter number " << paramNo << ", which is larger than the number of parameters = " << param->size() << std::endl;
1667  return;
1668  }
1669 
1670  // check if N0 is fixed. If this is the case, do NOT estimate N0
1671  if (param->at(paramNo-1).fStep == 0.0) // N0 parameter fixed
1672  return;
1673 
1674 
1675  // check that 'backgr.fit' in the msr-file run block is indeed a parameter number.
1676  // in case it is a function, nothing will be done.
1677  Int_t paramNoBkg = fRunInfo->GetBkgFitParamNo();
1678  Bool_t scaleBkg = true;
1679  Double_t bkg=0.0, errBkg=1.0;
1680  if ((paramNoBkg > 10000) || (paramNoBkg == -1)) { // i.e. fun or map
1681  scaleBkg = false;
1682  } else {
1683  if (paramNoBkg-1 < static_cast<Int_t>(param->size())) {
1684  bkg = param->at(paramNoBkg-1).fValue;
1685  errBkg = param->at(paramNoBkg-1).fStep;
1686  }
1687  }
1688 
1689  // estimate N0
1690  Double_t dt = fTimeResolution;
1691  Double_t tau = PMUON_LIFETIME;
1692 
1693  UInt_t t0 = static_cast<UInt_t>(round(fT0s[0]));
1694  Double_t dval = 0.0;
1695  Double_t nom = 0.0;
1696  Double_t denom = 0.0;
1697  Double_t xx = 0.0;
1698 
1699  // calc nominator
1700  for (UInt_t i=t0; i<fForward.size(); i++) {
1701  xx = exp(-dt*static_cast<Double_t>(i-t0)/tau);
1702  nom += xx;
1703  }
1704 
1705  // calc denominator
1706  for (UInt_t i=t0; i<fForward.size(); i++) {
1707  xx = exp(-dt*static_cast<Double_t>(i-t0)/tau);
1708  dval = fForward[i];
1709  if (dval > 0)
1710  denom += xx*xx/dval;
1711  }
1712  Double_t N0 = nom/denom;
1713 
1714  if (fScaleN0AndBkg) {
1715  N0 /= fTimeResolution*1.0e3;
1716  } else {
1717  N0 *= fPacking;
1718  }
1719 
1720  Double_t rescale = 1;
1721  if ((param->at(paramNo-1).fValue != 0.0) && scaleBkg) {
1722  rescale = N0 / param->at(paramNo-1).fValue;
1723  bkg *= rescale;
1724  errBkg *= rescale;
1725  }
1726 
1727  std::cout << ">> PRunSingleHisto::EstimateN0: found N0=" << param->at(paramNo-1).fValue << ", will set it to N0=" << N0 << std::endl;
1728  if (scaleBkg)
1729  std::cout << ">> PRunSingleHisto::EstimateN0: found Bkg=" << param->at(paramNoBkg-1).fValue << ", will set it to Bkg=" << bkg << std::endl;
1730  fMsrInfo->SetMsrParamValue(paramNo-1, N0);
1731  fMsrInfo->SetMsrParamStep(paramNo-1, sqrt(fabs(N0)));
1732  if (scaleBkg) {
1733  fMsrInfo->SetMsrParamValue(paramNoBkg-1, bkg);
1734  fMsrInfo->SetMsrParamStep(paramNoBkg-1, errBkg);
1735  }
1736 }
1737 
1738 //--------------------------------------------------------------------------
1739 // EstimatBkg (private)
1740 //--------------------------------------------------------------------------
1750 Bool_t PRunSingleHisto::EstimateBkg(UInt_t histoNo)
1751 {
1752  Double_t beamPeriod = 0.0;
1753 
1754  // check if data are from PSI, RAL, or TRIUMF
1755  if (fRunInfo->GetInstitute()->Contains("psi"))
1756  beamPeriod = ACCEL_PERIOD_PSI;
1757  else if (fRunInfo->GetInstitute()->Contains("ral"))
1758  beamPeriod = ACCEL_PERIOD_RAL;
1759  else if (fRunInfo->GetInstitute()->Contains("triumf"))
1760  beamPeriod = ACCEL_PERIOD_TRIUMF;
1761  else
1762  beamPeriod = 0.0;
1763 
1764  // check if start and end are in proper order
1765  Int_t start = fRunInfo->GetBkgRange(0);
1766  Int_t end = fRunInfo->GetBkgRange(1);
1767  if (end < start) {
1768  std::cout << std::endl << "PRunSingleHisto::EstimatBkg(): end = " << end << " > start = " << start << "! Will swap them!";
1769  Int_t keep = end;
1770  end = start;
1771  start = keep;
1772  }
1773 
1774  // calculate proper background range
1775  if (beamPeriod != 0.0) {
1776  Double_t timeBkg = static_cast<Double_t>(end-start)*(fTimeResolution*fPacking); // length of the background intervall in time
1777  UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
1778  // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
1779  end = start + static_cast<UInt_t>((fullCycles*beamPeriod)/(fTimeResolution*fPacking));
1780  std::cout << std::endl << "PRunSingleHisto::EstimatBkg(): Background " << start << ", " << end;
1781  if (end == start)
1782  end = fRunInfo->GetBkgRange(1);
1783  }
1784 
1785  // check if start is within histogram bounds
1786  if (start >= fForward.size()) {
1787  std::cerr << std::endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
1788  std::cerr << std::endl << ">> histo lengths = " << fForward.size();
1789  std::cerr << std::endl << ">> background start = " << start;
1790  std::cerr << std::endl;
1791  return false;
1792  }
1793 
1794  // check if end is within histogram bounds
1795  if (end >= fForward.size()) {
1796  std::cerr << std::endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
1797  std::cerr << std::endl << ">> histo lengths = " << fForward.size();
1798  std::cerr << std::endl << ">> background end = " << end;
1799  std::cerr << std::endl;
1800  return false;
1801  }
1802 
1803  // calculate background
1804  Double_t bkg = 0.0;
1805 
1806  // forward
1807  for (UInt_t i=start; i<end; i++)
1808  bkg += fForward[i];
1809  bkg /= static_cast<Double_t>(end - start + 1);
1810 
1811  if (fScaleN0AndBkg)
1812  fBackground = bkg / (fTimeResolution * 1e3); // keep background (per 1 nsec) for chisq, max.log.likelihood, fTimeResolution us->ns
1813  else
1814  fBackground = bkg * fPacking; // keep background (per bin)
1815 
1817 
1818  return true;
1819 }
1820 
1821 //--------------------------------------------------------------------------
1822 // IsScaleN0AndBkg (private)
1823 //--------------------------------------------------------------------------
1835 {
1836  Bool_t willScale = true;
1837 
1838  PMsrLines *cmd = fMsrInfo->GetMsrCommands();
1839  for (UInt_t i=0; i<cmd->size(); i++) {
1840  if (cmd->at(i).fLine.Contains("SCALE_N0_BKG", TString::kIgnoreCase)) {
1841  TObjArray *tokens = nullptr;
1842  TObjString *ostr = nullptr;
1843  TString str;
1844  tokens = cmd->at(i).fLine.Tokenize(" \t");
1845  if (tokens->GetEntries() != 2) {
1846  std::cerr << std::endl << ">> PRunSingleHisto::IsScaleN0AndBkg(): **WARNING** Found uncorrect 'SCALE_N0_BKG' command, will ignore it.";
1847  std::cerr << std::endl << ">> Allowed commands: SCALE_N0_BKG TRUE | FALSE" << std::endl;
1848  return willScale;
1849  }
1850  ostr = dynamic_cast<TObjString*>(tokens->At(1));
1851  str = ostr->GetString();
1852  if (!str.CompareTo("FALSE", TString::kIgnoreCase)) {
1853  willScale = false;
1854  }
1855  // clean up
1856  if (tokens)
1857  delete tokens;
1858  }
1859  }
1860 
1861  return willScale;
1862 }
virtual UInt_t GetNoOfFitBins()
virtual UInt_t GetT0BinSize()
Definition: PMusr.h:583
Bool_t fValid
flag showing if the state of the class is valid
Definition: PRunBase.h:66
virtual void SetDataRange(Int_t ival, Int_t idx)
Definition: PMusr.cpp:1672
virtual const PDoubleVector * GetError()
Definition: PMusr.h:249
virtual Bool_t EstimateBkg(UInt_t histoNo)
Definition: PMusr.h:220
virtual Double_t CalcMaxLikelihoodExpected(const std::vector< Double_t > &par)
virtual Int_t GetNormParamNo()
Definition: PMusr.h:646
virtual Int_t GetNoOfFuncs()
Definition: PMsrHandler.h:93
virtual Int_t GetDataRange(UInt_t idx)
Definition: PMusr.cpp:895
virtual Double_t GetT0Bin(UInt_t idx=0)
Definition: PMusr.cpp:935
#define PMUON_LIFETIME
Definition: PMusr.h:68
virtual const Double_t GetTimeResolution()
Definition: PMusr.h:429
Int_t fPacking
packing for this particular run. Either given in the RUN- or GLOBAL-block.
virtual void CalcNoOfFitBins()
virtual Int_t GetPacking()
Definition: PMusr.h:668
virtual Bool_t PrepareViewData(PRawRunData *runData, const UInt_t histoNo)
virtual Bool_t GetProperDataRange()
PRunDataHandler * fRawData
holds the raw run data
Definition: PRunBase.h:73
std::vector< PDoubleVector > fAddT0s
all t0 bins of all addrun&#39;s of a run! The derived classes will handle it.
Definition: PRunBase.h:79
virtual void SetBkgEstimated(Double_t dval, Int_t idx)
Definition: PMusr.cpp:1549
virtual Int_t GetFitRangeOffset(UInt_t idx)
Definition: PMusr.cpp:1088
virtual void SetTheoryTimeStart(Double_t dval)
Definition: PMusr.h:255
virtual Int_t GetBkgFitParamNo()
Definition: PMusr.h:647
PMsrHandler * fMsrInfo
msr-file handler
Definition: PRunBase.h:71
virtual Bool_t IsFitRangeInBin()
Definition: PMusr.h:665
virtual const PDoublePairVector * GetTemperature() const
Definition: PMusr.h:422
UInt_t fNoOfFitBins
number of bins to be fitted
Double_t fTimeResolution
time resolution in (us)
Definition: PRunBase.h:76
virtual const PDoubleVector * GetTheory()
Definition: PMusr.h:251
PMetaData fMetaData
keeps the meta data from the data file like field, temperature, energy, ...
Definition: PRunBase.h:77
virtual PRawRunData * GetRunData(const TString &runName)
#define RRF_UNIT_MHz
Definition: PMusr.h:155
#define RRF_UNIT_kHz
Definition: PMusr.h:154
virtual Double_t GetAddT0Bin(UInt_t addRunIdx, UInt_t histoIdx)
Definition: PMusr.cpp:1761
virtual UInt_t GetT0BinSize()
Definition: PMusr.h:660
virtual void EstimateN0()
PDoubleVector fTemp
temperature(s) in (K)
Definition: PMusr.h:229
std::vector< PMsrPlotStructure > PMsrPlotList
Definition: PMusr.h:802
virtual void AppendErrorValue(Double_t dval)
Definition: PMusr.h:260
virtual Bool_t IsFitRangeInBin()
Definition: PMusr.h:588
virtual Double_t GetTheoryTimeStep()
Definition: PMusr.h:245
std::vector< UInt_t > PUIntVector
Definition: PMusr.h:172
PRunData fData
data to be fitted, viewed, i.e. binned data
Definition: PRunBase.h:75
virtual void SetDataTimeStep(Double_t dval)
Definition: PMusr.h:254
virtual Double_t GetTheoryTimeStart()
Definition: PMusr.h:244
virtual Double_t GetT0Bin(UInt_t idx=0)
Definition: PMusr.cpp:1695
std::vector< PMsrParamStructure > PMsrParamList
Definition: PMusr.h:564
std::vector< Double_t > PDoubleVector
Definition: PMusr.h:196
virtual const Bool_t IsPresent(UInt_t histoNo)
Definition: PMusr.h:430
#define MSR_PARAM_FUN_OFFSET
Definition: PMusr.h:127
#define ACCEL_PERIOD_TRIUMF
Definition: PMusr.h:82
virtual const Double_t GetT0Bin(const UInt_t histoNo)
Definition: PMusr.h:431
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
Int_t fGoodBins[2]
keep first/last good bins. 0=fgb, 1=lgb
virtual Bool_t IsScaleN0AndBkg()
virtual void SetT0Bin(Double_t dval, Int_t idx=-1)
Definition: PMusr.cpp:1712
#define RRF_UNIT_T
Definition: PMusr.h:158
virtual void AppendTheoryValue(Double_t dval)
Definition: PMusr.h:262
virtual PMsrGlobalBlock * GetMsrGlobal()
Definition: PMsrHandler.h:62
virtual Int_t GetLifetimeParamNo()
Definition: PMusr.h:648
std::vector< PMsrLineStructure > PMsrLines
Definition: PMusr.h:539
virtual void SetDataTimeStart(Double_t dval)
Definition: PMusr.h:253
virtual Int_t GetFitRangeOffset(UInt_t idx)
Definition: PMusr.cpp:1848
virtual void AppendValue(Double_t dval)
Definition: PMusr.h:259
#define RRF_UNIT_G
Definition: PMusr.h:157
virtual const PDoubleVector * GetValue()
Definition: PMusr.h:248
#define RRF_UNIT_Mcs
Definition: PMusr.h:156
virtual const Double_t GetField()
Definition: PMusr.h:420
virtual Bool_t PrepareData()
virtual void SetAddT0Bin(Double_t dval, UInt_t addRunIdx, UInt_t histoNoIdx)
Definition: PMusr.cpp:1788
virtual Bool_t SetMsrParamStep(UInt_t i, Double_t value)
virtual void CalculateKaiserFilterCoeff(Double_t wc, Double_t A, Double_t dw)
Definition: PRunBase.cpp:200
virtual Double_t GetBkgFix(UInt_t idx)
Definition: PMusr.cpp:1572
virtual UInt_t GetFuncNo(Int_t idx)
Definition: PMsrHandler.h:94
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &histoNo)
virtual ~PRunSingleHisto()
virtual const UInt_t GetNoOfTemperatures()
Definition: PMusr.h:421
#define PMUSR_UNDEFINED
Definition: PMusr.h:89
virtual Int_t GetForwardHistoNo(UInt_t idx=0)
Definition: PMusr.cpp:1400
PDoubleVector fT0s
all t0 bins of a run! The derived classes will handle it.
Definition: PRunBase.h:78
virtual PMsrPlotList * GetMsrPlotList()
Definition: PMsrHandler.h:66
Double_t fEnergy
energy in (keV)
Definition: PMusr.h:228
virtual Int_t GetBkgRange(UInt_t idx)
Definition: PMusr.cpp:1612
EPMusrHandleTag fHandleTag
tag telling whether this is used for fit, view, ...
Definition: PRunBase.h:68
virtual UInt_t GetForwardHistoNoSize()
Definition: PMusr.h:652
EPMusrHandleTag
Definition: PMusr.h:220
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
Int_t fStartTimeBin
bin at which the fit starts
virtual Bool_t PrepareRawViewData(PRawRunData *runData, const UInt_t histoNo)
PMsrRunBlock * fRunInfo
run info used to filter out needed infos of a run
Definition: PRunBase.h:72
virtual const Double_t GetEnergy()
Definition: PMusr.h:425
Int_t fEndTimeBin
bin at which the fit ends
virtual Bool_t SetMsrParamValue(UInt_t i, Double_t value)
virtual void SetBkgRange(Int_t ival, Int_t idx)
Definition: PMusr.cpp:1630
virtual void SetFitRange(Double_t dval, UInt_t idx)
Definition: PMusr.cpp:1828
virtual Double_t CalcChiSquare(const std::vector< Double_t > &par)
Double_t fFitEndTime
fit end time
Definition: PRunBase.h:82
virtual void SetTheoryTimeStep(Double_t dval)
Definition: PMusr.h:256
virtual TString * GetInstitute(UInt_t idx=0)
Definition: PMusr.cpp:1316
Double_t fBackground
needed if background range is given (units: 1/bin)
#define ACCEL_PERIOD_PSI
Definition: PMusr.h:81
virtual void SetFitRange(Double_t dval, UInt_t idx)
Definition: PMusr.cpp:1068
Double_t fFitStartTime
fit start time
Definition: PRunBase.h:81
virtual PMsrParamList * GetMsrParamList()
Definition: PMsrHandler.h:59
Bool_t fScaleN0AndBkg
true=scale N0 and background to 1/ns, otherwise 1/bin
virtual PIntVector * GetMap()
Definition: PMusr.h:650
virtual Double_t EvalFunc(UInt_t i, std::vector< Int_t > map, std::vector< Double_t > param, PMetaData metaData)
Definition: PMsrHandler.h:98
virtual const PDoubleVector * GetDataBin(const UInt_t histoNo)
Definition: PMusr.h:438
virtual Double_t GetDataTimeStart()
Definition: PMusr.h:242
virtual UInt_t GetRunNameSize()
Definition: PMusr.h:635
virtual Int_t GetDataRange(UInt_t idx)
Definition: PMusr.cpp:1654
virtual PMsrLines * GetMsrCommands()
Definition: PMsrHandler.h:64
virtual const Double_t GetT0BinEstimated(const UInt_t histoNo)
Definition: PMusr.h:432
virtual void FilterTheo()
Definition: PRunBase.cpp:240
PDoubleVector fForward
forward histo data
std::unique_ptr< PTheory > fTheory
theory needed to calculate chi-square
Definition: PRunBase.h:85
virtual Bool_t EstimateN0()
virtual void SetFitRangeBin(const TString fitRange)
virtual Int_t GetPacking()
Definition: PMusr.h:591
virtual Double_t GetFitRange(UInt_t idx)
Definition: PMusr.cpp:1051
Bool_t fTheoAsData
true=only calculate the theory points at the data points, false=calculate more points for the theory ...
virtual void ReplaceTheory(const PDoubleVector &theo)
Definition: PMusr.cpp:103
virtual Double_t GetDataTimeStep()
Definition: PMusr.h:243
Int_t fRunNo
number of the run within the msr-file
Definition: PRunBase.h:70
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
#define GAMMA_BAR_MUON
Definition: PMusr.h:78
virtual TString * GetRunName(UInt_t idx=0)
Definition: PMusr.cpp:1232
virtual Double_t GetFitRange(UInt_t idx)
Definition: PMusr.cpp:1811
Definition: PMusr.h:220
virtual Bool_t PrepareFitData(PRawRunData *runData, const UInt_t histoNo)
PDoubleVector fFuncValues
is keeping the values of the functions from the FUNCTIONS block
Definition: PRunBase.h:84
virtual void CalcTheory()
Double_t fField
field in (G)
Definition: PMusr.h:227
#define ACCEL_PERIOD_RAL
Definition: PMusr.h:83