musrfit  1.9.2
PRunAsymmetryRRF.cpp
Go to the documentation of this file.
1 /***************************************************************************
2 
3  PRunAsymmetryRRF.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 <stdio.h>
39 
40 #include <iostream>
41 #include <fstream>
42 
43 #include <TString.h>
44 #include <TObjArray.h>
45 #include <TObjString.h>
46 
47 #include "PMusr.h"
48 #include "PRunAsymmetryRRF.h"
49 
50 //--------------------------------------------------------------------------
51 // Constructor
52 //--------------------------------------------------------------------------
57 {
58  fNoOfFitBins = 0;
59  fRRFPacking = -1;
60  fTheoAsData = false;
61 
62  // the 2 following variables are need in case fit range is given in bins, and since
63  // the fit range can be changed in the command block, these variables need to be accessible
64  fGoodBins[0] = -1;
65  fGoodBins[1] = -1;
66 }
67 
68 //--------------------------------------------------------------------------
69 // Constructor
70 //--------------------------------------------------------------------------
79 PRunAsymmetryRRF::PRunAsymmetryRRF(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
80  PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
81 {
82  // the 2 following variables are need in case fit range is given in bins, and since
83  // the fit range can be changed in the command block, these variables need to be accessible
84  fGoodBins[0] = -1;
85  fGoodBins[1] = -1;
86 
88  if (fRRFPacking == -1) { // this should NOT happen, somethin is severely wrong
89  std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **SEVERE ERROR**: Couldn't find any RRF packing information!";
90  std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
91  std::cerr << std::endl;
92  fValid = false;
93  return;
94  }
95 
96  // check if alpha and/or beta is fixed --------------------
97 
98  PMsrParamList *param = msrInfo->GetMsrParamList();
99 
100  // check if alpha is given
101  if (fRunInfo->GetAlphaParamNo() == -1) { // no alpha given
102  std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **ERROR** no alpha parameter given! This is needed for an asymmetry fit!";
103  std::cerr << std::endl;
104  fValid = false;
105  return;
106  }
107  // check if alpha parameter is within proper bounds
108  if ((fRunInfo->GetAlphaParamNo() < 0) || (fRunInfo->GetAlphaParamNo() > static_cast<Int_t>(param->size()))) {
109  std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **ERROR** alpha parameter no = " << fRunInfo->GetAlphaParamNo();
110  std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
111  std::cerr << std::endl;
112  fValid = false;
113  return;
114  }
115  // check if alpha is fixed
116  Bool_t alphaFixedToOne = false;
117  if (((*param)[fRunInfo->GetAlphaParamNo()-1].fStep == 0.0) &&
118  ((*param)[fRunInfo->GetAlphaParamNo()-1].fValue == 1.0))
119  alphaFixedToOne = true;
120 
121  // check if beta is given
122  Bool_t betaFixedToOne = false;
123  if (fRunInfo->GetBetaParamNo() == -1) { // no beta given hence assuming beta == 1
124  betaFixedToOne = true;
125  } else if ((fRunInfo->GetBetaParamNo() < 0) || (fRunInfo->GetBetaParamNo() > static_cast<Int_t>(param->size()))) { // check if beta parameter is within proper bounds
126  std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **ERROR** beta parameter no = " << fRunInfo->GetBetaParamNo();
127  std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
128  std::cerr << std::endl;
129  fValid = false;
130  return;
131  } else { // check if beta is fixed
132  if (((*param)[fRunInfo->GetBetaParamNo()-1].fStep == 0.0) &&
133  ((*param)[fRunInfo->GetBetaParamNo()-1].fValue == 1.0))
134  betaFixedToOne = true;
135  }
136 
137  // set fAlphaBetaTag
138  if (alphaFixedToOne && betaFixedToOne) // alpha == 1, beta == 1
139  fAlphaBetaTag = 1;
140  else if (!alphaFixedToOne && betaFixedToOne) // alpha != 1, beta == 1
141  fAlphaBetaTag = 2;
142  else if (alphaFixedToOne && !betaFixedToOne) // alpha == 1, beta != 1
143  fAlphaBetaTag = 3;
144  else
145  fAlphaBetaTag = 4;
146 
147  // calculate fData
148  if (!PrepareData())
149  fValid = false;
150 }
151 
152 //--------------------------------------------------------------------------
153 // Destructor
154 //--------------------------------------------------------------------------
159 {
160  fForward.clear();
161  fForwardErr.clear();
162  fBackward.clear();
163  fBackwardErr.clear();
164 }
165 
166 //--------------------------------------------------------------------------
167 // CalcChiSquare (public)
168 //--------------------------------------------------------------------------
177 Double_t PRunAsymmetryRRF::CalcChiSquare(const std::vector<Double_t>& par)
178 {
179  Double_t chisq = 0.0;
180  Double_t diff = 0.0;
181  Double_t asymFcnValue = 0.0;
182  Double_t a, b, f;
183 
184  // calculate functions
185  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
187  }
188 
189  // calculate chi square
190  Double_t time(1.0);
191  Int_t i;
192 
193  // determine alpha/beta
194  switch (fAlphaBetaTag) {
195  case 1: // alpha == 1, beta == 1
196  a = 1.0;
197  b = 1.0;
198  break;
199  case 2: // alpha != 1, beta == 1
200  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
201  a = par[fRunInfo->GetAlphaParamNo()-1];
202  } else { // alpha is function
203  // get function number
205  // evaluate function
206  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
207  }
208  b = 1.0;
209  break;
210  case 3: // alpha == 1, beta != 1
211  a = 1.0;
212  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
213  b = par[fRunInfo->GetBetaParamNo()-1];
214  } else { // beta is a function
215  // get function number
216  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
217  // evaluate function
218  b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
219  }
220  break;
221  case 4: // alpha != 1, beta != 1
222  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
223  a = par[fRunInfo->GetAlphaParamNo()-1];
224  } else { // alpha is function
225  // get function number
227  // evaluate function
228  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
229  }
230  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
231  b = par[fRunInfo->GetBetaParamNo()-1];
232  } else { // beta is a function
233  // get function number
234  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
235  // evaluate function
236  b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
237  }
238  break;
239  default:
240  a = 1.0;
241  b = 1.0;
242  break;
243  }
244 
245  // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
246  // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
247  // for a given set of parameters---which should be done outside of the parallelized loop.
248  // For all other functions it means a tiny and acceptable overhead.
249  asymFcnValue = fTheory->Func(time, par, fFuncValues);
250 
251  #ifdef HAVE_GOMP
252  Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
253  if (chunk < 10)
254  chunk = 10;
255  #pragma omp parallel for default(shared) private(i,time,diff,asymFcnValue,f) schedule(dynamic,chunk) reduction(+:chisq)
256  #endif
257  for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
258  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
259  f = fTheory->Func(time, par, fFuncValues);
260  asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
261  diff = fData.GetValue()->at(i) - asymFcnValue;
262  chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
263  }
264 
265  return chisq;
266 }
267 
268 //--------------------------------------------------------------------------
269 // CalcChiSquareExpected (public)
270 //--------------------------------------------------------------------------
279 Double_t PRunAsymmetryRRF::CalcChiSquareExpected(const std::vector<Double_t>& par)
280 {
281  return 0.0;
282 }
283 
284 //--------------------------------------------------------------------------
285 // CalcMaxLikelihood (public)
286 //--------------------------------------------------------------------------
292 Double_t PRunAsymmetryRRF::CalcMaxLikelihood(const std::vector<Double_t>& par)
293 {
294  std::cout << std::endl << "PRunAsymmetryRRF::CalcMaxLikelihood(): not implemented yet ..." << std::endl;
295 
296  return 1.0;
297 }
298 
299 //--------------------------------------------------------------------------
300 // GetNoOfFitBins (public)
301 //--------------------------------------------------------------------------
308 {
309  CalcNoOfFitBins();
310 
311  return fNoOfFitBins;
312 }
313 
314 //--------------------------------------------------------------------------
315 // SetFitRangeBin (public)
316 //--------------------------------------------------------------------------
328 void PRunAsymmetryRRF::SetFitRangeBin(const TString fitRange)
329 {
330  TObjArray *tok = nullptr;
331  TObjString *ostr = nullptr;
332  TString str;
333  Ssiz_t idx = -1;
334  Int_t offset = 0;
335 
336  tok = fitRange.Tokenize(" \t");
337 
338  if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
339  // handle fgb+n0 entry
340  ostr = static_cast<TObjString*>(tok->At(1));
341  str = ostr->GetString();
342  // check if there is an offset present
343  idx = str.First("+");
344  if (idx != -1) { // offset present
345  str.Remove(0, idx+1);
346  if (str.IsFloat()) // if str is a valid number, convert is to an integer
347  offset = str.Atoi();
348  }
349  fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
350 
351  // handle lgb-n1 entry
352  ostr = static_cast<TObjString*>(tok->At(2));
353  str = ostr->GetString();
354  // check if there is an offset present
355  idx = str.First("-");
356  if (idx != -1) { // offset present
357  str.Remove(0, idx+1);
358  if (str.IsFloat()) // if str is a valid number, convert is to an integer
359  offset = str.Atoi();
360  }
361  fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
362  } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
363  Int_t pos = 2*(fRunNo+1)-1;
364 
365  if (pos + 1 >= tok->GetEntries()) {
366  std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
367  std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
368  } else {
369  // handle fgb+n0 entry
370  ostr = static_cast<TObjString*>(tok->At(pos));
371  str = ostr->GetString();
372  // check if there is an offset present
373  idx = str.First("+");
374  if (idx != -1) { // offset present
375  str.Remove(0, idx+1);
376  if (str.IsFloat()) // if str is a valid number, convert is to an integer
377  offset = str.Atoi();
378  }
379  fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
380 
381  // handle lgb-n1 entry
382  ostr = static_cast<TObjString*>(tok->At(pos+1));
383  str = ostr->GetString();
384  // check if there is an offset present
385  idx = str.First("-");
386  if (idx != -1) { // offset present
387  str.Remove(0, idx+1);
388  if (str.IsFloat()) // if str is a valid number, convert is to an integer
389  offset = str.Atoi();
390  }
391  fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
392  }
393  } else { // error
394  std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
395  std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
396  }
397 
398  // clean up
399  if (tok) {
400  delete tok;
401  }
402 }
403 
404 //--------------------------------------------------------------------------
405 // CalcNoOfFitBins (public)
406 //--------------------------------------------------------------------------
411 {
412  // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
413  fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
414  if (fStartTimeBin < 0)
415  fStartTimeBin = 0;
416  fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
417  if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
418  fEndTimeBin = fData.GetValue()->size();
419 
422  else
423  fNoOfFitBins = 0;
424 }
425 
426 //--------------------------------------------------------------------------
427 // CalcTheory (protected)
428 //--------------------------------------------------------------------------
433 {
434  // feed the parameter vector
435  std::vector<Double_t> par;
436  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
437  for (UInt_t i=0; i<paramList->size(); i++)
438  par.push_back((*paramList)[i].fValue);
439 
440  // calculate functions
441  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
443  }
444 
445  // calculate asymmetry
446  Double_t asymFcnValue = 0.0;
447  Double_t a, b, f;
448  Double_t time;
449  for (UInt_t i=0; i<fData.GetValue()->size(); i++) {
450  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
451  switch (fAlphaBetaTag) {
452  case 1: // alpha == 1, beta == 1
453  asymFcnValue = fTheory->Func(time, par, fFuncValues);
454  break;
455  case 2: // alpha != 1, beta == 1
456  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
457  a = par[fRunInfo->GetAlphaParamNo()-1];
458  } else { // alpha is function
459  // get function number
461  // evaluate function
462  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
463  }
464  f = fTheory->Func(time, par, fFuncValues);
465  asymFcnValue = (f*(a+1.0)-(a-1.0))/((a+1.0)-f*(a-1.0));
466  break;
467  case 3: // alpha == 1, beta != 1
468  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
469  b = par[fRunInfo->GetBetaParamNo()-1];
470  } else { // beta is a function
471  // get function number
472  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
473  // evaluate function
474  b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
475  }
476  f = fTheory->Func(time, par, fFuncValues);
477  asymFcnValue = f*(b+1.0)/(2.0-f*(b-1.0));
478  break;
479  case 4: // alpha != 1, beta != 1
480  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
481  a = par[fRunInfo->GetAlphaParamNo()-1];
482  } else { // alpha is function
483  // get function number
485  // evaluate function
486  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
487  }
488  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
489  b = par[fRunInfo->GetBetaParamNo()-1];
490  } else { // beta is a function
491  // get function number
492  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
493  // evaluate function
494  b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
495  }
496  f = fTheory->Func(time, par, fFuncValues);
497  asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
498  break;
499  default:
500  asymFcnValue = 0.0;
501  break;
502  }
503  fData.AppendTheoryValue(asymFcnValue);
504  }
505 
506  // clean up
507  par.clear();
508 }
509 
510 //--------------------------------------------------------------------------
511 // PrepareData (protected)
512 //--------------------------------------------------------------------------
533 {
534  if (!fValid)
535  return false;
536 
537  // keep the Global block info
538  PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
539 
540  // get forward/backward histo from PRunDataHandler object ------------------------
541  // get the correct run
543  if (!runData) { // run not found
544  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
545  std::cerr << std::endl;
546  return false;
547  }
548 
549  // keep the field from the meta-data from the data-file
550  fMetaData.fField = runData->GetField();
551 
552  // keep the energy from the meta-data from the data-file
553  fMetaData.fEnergy = runData->GetEnergy();
554 
555  // keep the temperature(s) from the meta-data from the data-file
556  for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
557  fMetaData.fTemp.push_back(runData->GetTemperature(i));
558 
559  // collect histogram numbers
560  PUIntVector forwardHistoNo;
561  PUIntVector backwardHistoNo;
562  for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
563  forwardHistoNo.push_back(fRunInfo->GetForwardHistoNo(i));
564 
565  if (!runData->IsPresent(forwardHistoNo[i])) {
566  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **PANIC ERROR**:";
567  std::cerr << std::endl << ">> forwardHistoNo found = " << forwardHistoNo[i] << ", which is NOT present in the data file!?!?";
568  std::cerr << std::endl << ">> Will quit :-(";
569  std::cerr << std::endl;
570  // clean up
571  forwardHistoNo.clear();
572  backwardHistoNo.clear();
573  return false;
574  }
575  }
576  for (UInt_t i=0; i<fRunInfo->GetBackwardHistoNoSize(); i++) {
577  backwardHistoNo.push_back(fRunInfo->GetBackwardHistoNo(i));
578 
579  if (!runData->IsPresent(backwardHistoNo[i])) {
580  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **PANIC ERROR**:";
581  std::cerr << std::endl << ">> backwardHistoNo found = " << backwardHistoNo[i] << ", which is NOT present in the data file!?!?";
582  std::cerr << std::endl << ">> Will quit :-(";
583  std::cerr << std::endl;
584  // clean up
585  forwardHistoNo.clear();
586  backwardHistoNo.clear();
587  return false;
588  }
589  }
590 
591  // keep the time resolution in (us)
592  fTimeResolution = runData->GetTimeResolution()/1.0e3;
593  std::cout.precision(10);
594  std::cout << std::endl << ">> PRunAsymmetryRRF::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
595 
596  // get all the proper t0's and addt0's for the current RUN block
597  if (!GetProperT0(runData, globalBlock, forwardHistoNo, backwardHistoNo)) {
598  return false;
599  }
600 
601  // keep the histo of each group at this point (addruns handled below)
602  std::vector<PDoubleVector> forward, backward;
603  forward.resize(forwardHistoNo.size()); // resize to number of groups
604  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
605  forward[i].resize(runData->GetDataBin(forwardHistoNo[i])->size());
606  forward[i] = *runData->GetDataBin(forwardHistoNo[i]);
607  }
608  backward.resize(backwardHistoNo.size()); // resize to number of groups
609  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
610  backward[i].resize(runData->GetDataBin(backwardHistoNo[i])->size());
611  backward[i] = *runData->GetDataBin(backwardHistoNo[i]);
612  }
613 
614  // check if addrun's are present, and if yes add data
615  // check if there are runs to be added to the current one
616  if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
617  PRawRunData *addRunData;
618  for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
619  // get run to be added to the main one
620  addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
621  if (addRunData == nullptr) { // couldn't get run
622  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
623  std::cerr << std::endl;
624  return false;
625  }
626 
627  // add forward run
628  UInt_t addRunSize;
629  for (UInt_t k=0; k<forwardHistoNo.size(); k++) { // fill each group
630  addRunSize = addRunData->GetDataBin(forwardHistoNo[k])->size();
631  for (UInt_t j=0; j<addRunData->GetDataBin(forwardHistoNo[k])->size(); j++) { // loop over the bin indices
632  // make sure that the index stays in the proper range
633  if (((Int_t)j+(Int_t)fAddT0s[i-1][2*k]-(Int_t)fT0s[2*k] >= 0) && (j+(Int_t)fAddT0s[i-1][2*k]-(Int_t)fT0s[2*k] < addRunSize)) {
634  forward[k][j] += addRunData->GetDataBin(forwardHistoNo[k])->at(j+(Int_t)fAddT0s[i-1][2*k]-(Int_t)fT0s[2*k]);
635  }
636  }
637  }
638 
639  // add backward run
640  for (UInt_t k=0; k<backwardHistoNo.size(); k++) { // fill each group
641  addRunSize = addRunData->GetDataBin(backwardHistoNo[k])->size();
642  for (UInt_t j=0; j<addRunData->GetDataBin(backwardHistoNo[k])->size(); j++) { // loop over the bin indices
643  // make sure that the index stays in the proper range
644  if (((Int_t)j+(Int_t)fAddT0s[i-1][2*k+1]-(Int_t)fT0s[2*k+1] >= 0) && (j+(Int_t)fAddT0s[i-1][2*k+1]-(Int_t)fT0s[2*k+1] < addRunSize)) {
645  backward[k][j] += addRunData->GetDataBin(backwardHistoNo[k])->at(j+(Int_t)fAddT0s[i-1][2*k+1]-(Int_t)fT0s[2*k+1]);
646  }
647  }
648  }
649  }
650  }
651 
652  // set forward/backward histo data of the first group
653  fForward.resize(forward[0].size());
654  fBackward.resize(backward[0].size());
655  for (UInt_t i=0; i<fForward.size(); i++) {
656  fForward[i] = forward[0][i];
657  fBackward[i] = backward[0][i];
658  }
659 
660  // group histograms, add all the remaining forward histograms of the group
661  for (UInt_t i=1; i<forwardHistoNo.size(); i++) { // loop over the groupings
662  for (UInt_t j=0; j<runData->GetDataBin(forwardHistoNo[i])->size(); j++) { // loop over the bin indices
663  // make sure that the index stays within proper range
664  if (((Int_t)j+fT0s[2*i]-fT0s[0] >= 0) && (j+fT0s[2*i]-fT0s[0] < runData->GetDataBin(forwardHistoNo[i])->size())) {
665  fForward[j] += forward[i][j+(Int_t)fT0s[2*i]-(Int_t)fT0s[0]];
666  }
667  }
668  }
669 
670  // group histograms, add all the remaining backward histograms of the group
671  for (UInt_t i=1; i<backwardHistoNo.size(); i++) { // loop over the groupings
672  for (UInt_t j=0; j<runData->GetDataBin(backwardHistoNo[i])->size(); j++) { // loop over the bin indices
673  // make sure that the index stays within proper range
674  if ((j+fT0s[2*i+1]-fT0s[1] >= 0) && (j+fT0s[2*i+1]-fT0s[1] < runData->GetDataBin(backwardHistoNo[i])->size())) {
675  fBackward[j] += backward[i][j+(Int_t)fT0s[2*i+1]-(Int_t)fT0s[1]];
676  }
677  }
678  }
679 
680  // subtract background from histogramms ------------------------------------------
681  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
682  if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
683  if (!SubtractEstimatedBkg())
684  return false;
685  } else { // no background given to do the job, try to estimate it
686  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
687  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
688  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.1), 2);
689  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.6), 3);
690  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **WARNING** Neither fix background nor background bins are given!";
691  std::cerr << std::endl << ">> Will try the following:";
692  std::cerr << std::endl << ">> forward: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
693  std::cerr << std::endl << ">> backward: bkg start = " << fRunInfo->GetBkgRange(2) << ", bkg end = " << fRunInfo->GetBkgRange(3);
694  std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
695  std::cerr << std::endl;
696  if (!SubtractEstimatedBkg())
697  return false;
698  }
699  } else { // fixed background given
700  if (!SubtractFixBkg())
701  return false;
702  }
703 
704  UInt_t histoNo[2] = {forwardHistoNo[0], backwardHistoNo[0]};
705 
706  // get the data range (fgb/lgb) for the current RUN block
707  if (!GetProperDataRange(runData, histoNo)) {
708  return false;
709  }
710 
711  // get the fit range for the current RUN block
712  GetProperFitRange(globalBlock);
713 
714  // everything looks fine, hence fill data set
715  Bool_t status;
716  switch(fHandleTag) {
717  case kFit:
719  break;
720  case kView:
721  status = PrepareViewData(runData, histoNo);
722  break;
723  default:
724  status = false;
725  break;
726  }
727 
728  // clean up
729  forwardHistoNo.clear();
730  backwardHistoNo.clear();
731 
732  return status;
733 }
734 
735 //--------------------------------------------------------------------------
736 // SubtractFixBkg (private)
737 //--------------------------------------------------------------------------
755 {
756  Double_t dval;
757  for (UInt_t i=0; i<fForward.size(); i++) {
758  // keep the error, and make sure that the bin is NOT empty
759  if (fForward[i] != 0.0)
760  dval = TMath::Sqrt(fForward[i]);
761  else
762  dval = 1.0;
763  fForwardErr.push_back(dval);
764  fForward[i] -= fRunInfo->GetBkgFix(0);
765 
766  // keep the error, and make sure that the bin is NOT empty
767  if (fBackward[i] != 0.0)
768  dval = TMath::Sqrt(fBackward[i]);
769  else
770  dval = 1.0;
771  fBackwardErr.push_back(dval);
772  fBackward[i] -= fRunInfo->GetBkgFix(1);
773  }
774 
775  return true;
776 }
777 
778 //--------------------------------------------------------------------------
779 // SubtractEstimatedBkg (private)
780 //--------------------------------------------------------------------------
799 {
800  Double_t beamPeriod = 0.0;
801 
802  // check if data are from PSI, RAL, or TRIUMF
803  if (fRunInfo->GetInstitute()->Contains("psi"))
804  beamPeriod = ACCEL_PERIOD_PSI;
805  else if (fRunInfo->GetInstitute()->Contains("ral"))
806  beamPeriod = ACCEL_PERIOD_RAL;
807  else if (fRunInfo->GetInstitute()->Contains("triumf"))
808  beamPeriod = ACCEL_PERIOD_TRIUMF;
809  else
810  beamPeriod = 0.0;
811 
812  // check if start and end are in proper order
813  UInt_t start[2] = {static_cast<UInt_t>(fRunInfo->GetBkgRange(0)), static_cast<UInt_t>(fRunInfo->GetBkgRange(2))};
814  UInt_t end[2] = {static_cast<UInt_t>(fRunInfo->GetBkgRange(1)), static_cast<UInt_t>(fRunInfo->GetBkgRange(3))};
815  for (UInt_t i=0; i<2; i++) {
816  if (end[i] < start[i]) {
817  std::cout << std::endl << "PRunAsymmetryRRF::SubtractEstimatedBkg(): end = " << end[i] << " > start = " << start[i] << "! Will swap them!";
818  UInt_t keep = end[i];
819  end[i] = start[i];
820  start[i] = keep;
821  }
822  }
823 
824  // calculate proper background range
825  for (UInt_t i=0; i<2; i++) {
826  if (beamPeriod != 0.0) {
827  Double_t timeBkg = static_cast<Double_t>(end[i]-start[i])*fTimeResolution; // length of the background intervall in time
828  UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
829  // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
830  end[i] = start[i] + static_cast<UInt_t>((fullCycles*beamPeriod)/fTimeResolution);
831  std::cout << "PRunAsymmetryRRF::SubtractEstimatedBkg(): Background " << start[i] << ", " << end[i] << std::endl;
832  if (end[i] == start[i])
833  end[i] = fRunInfo->GetBkgRange(2*i+1);
834  }
835  }
836 
837  // check if start is within histogram bounds
838  if ((start[0] >= fForward.size()) || (start[1] >= fBackward.size())) {
839  std::cerr << std::endl << ">> PRunAsymmetryRRF::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
840  std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
841  std::cerr << std::endl << ">> background start (f/b) = (" << start[0] << "/" << start[1] << ").";
842  return false;
843  }
844 
845  // check if end is within histogram bounds
846  if ((end[0] >= fForward.size()) || (end[1] >= fBackward.size())) {
847  std::cerr << std::endl << ">> PRunAsymmetryRRF::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
848  std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
849  std::cerr << std::endl << ">> background end (f/b) = (" << end[0] << "/" << end[1] << ").";
850  return false;
851  }
852 
853  // calculate background
854  Double_t bkg[2] = {0.0, 0.0};
855  Double_t errBkg[2] = {0.0, 0.0};
856 
857  // forward
858  for (UInt_t i=start[0]; i<=end[0]; i++)
859  bkg[0] += fForward[i];
860  errBkg[0] = TMath::Sqrt(bkg[0])/(end[0] - start[0] + 1);
861  bkg[0] /= static_cast<Double_t>(end[0] - start[0] + 1);
862  std::cout << std::endl << ">> estimated forward histo background: " << bkg[0];
863 
864  // backward
865  for (UInt_t i=start[1]; i<=end[1]; i++)
866  bkg[1] += fBackward[i];
867  errBkg[1] = TMath::Sqrt(bkg[1])/(end[1] - start[1] + 1);
868  bkg[1] /= static_cast<Double_t>(end[1] - start[1] + 1);
869  std::cout << std::endl << ">> estimated backward histo background: " << bkg[1] << std::endl;
870 
871  // correct error for forward, backward
872  Double_t errVal = 0.0;
873  for (UInt_t i=0; i<fForward.size(); i++) {
874  if (fForward[i] > 0.0)
875  errVal = TMath::Sqrt(fForward[i]+errBkg[0]*errBkg[0]);
876  else
877  errVal = 1.0;
878  fForwardErr.push_back(errVal);
879  if (fBackward[i] > 0.0)
880  errVal = TMath::Sqrt(fBackward[i]+errBkg[1]*errBkg[1]);
881  else
882  errVal = 1.0;
883  fBackwardErr.push_back(errVal);
884  }
885 
886  // subtract background from data
887  for (UInt_t i=0; i<fForward.size(); i++) {
888  fForward[i] -= bkg[0];
889  fBackward[i] -= bkg[1];
890  }
891 
892  fRunInfo->SetBkgEstimated(bkg[0], 0);
893  fRunInfo->SetBkgEstimated(bkg[1], 1);
894 
895  return true;
896 }
897 
898 //--------------------------------------------------------------------------
899 // PrepareFitData (protected)
900 //--------------------------------------------------------------------------
906 {
907  // transform raw histo data. At this point, the raw data are already background corrected.
908 
909  // 1st: form the asymmetry of the original data
910 
911  // forward and backward detectors might have different fgb-t0 offset. Take the maximum of both.
912  Int_t fgbOffset = fGoodBins[0]-static_cast<Int_t>(fT0s[0]);
913  if (fgbOffset < fGoodBins[2]-static_cast<Int_t>(fT0s[1]))
914  fgbOffset = fGoodBins[2]-static_cast<Int_t>(fT0s[1]);
915  // last good bin (lgb) is the minimum of forward/backward lgb
916  Int_t lgb_offset = fGoodBins[1]-static_cast<Int_t>(fT0s[0])+fgbOffset;
917  if (lgb_offset < fGoodBins[3]-static_cast<Int_t>(fT0s[1])+fgbOffset)
918  lgb_offset = fGoodBins[3]-static_cast<Int_t>(fT0s[1])+fgbOffset;
919 
920  Int_t fgb = static_cast<Int_t>(fT0s[0])+fgbOffset;
921  Int_t lgb = fgb + lgb_offset;
922  Int_t dt0 = static_cast<Int_t>(fT0s[0])-static_cast<Int_t>(fT0s[1]);
923 
924  PDoubleVector asym;
925  PDoubleVector asymErr;
926  Double_t asymVal, asymValErr;
927  Double_t ff, bb, eff, ebb;
928  for (Int_t i=fgb; i<lgb; i++) {
929  ff = fForward[i];
930  bb = fBackward[i-dt0];
931  if (ff+bb != 0.0)
932  asymVal = (ff-bb)/(ff+bb);
933  else
934  asymVal = 0.0;
935  asym.push_back(asymVal);
936  eff = fForwardErr[i];
937  ebb = fBackwardErr[i-dt0];
938  if ((asymVal != 0.0) && (ff+bb) > 0.0)
939  asymValErr = 2.0/pow((ff+bb),2.0)*sqrt(bb*bb*eff*eff+ff*ff*ebb*ebb);
940  else
941  asymValErr = 1.0;
942  asymErr.push_back(asymValErr);
943  }
944 
945  // 2nd: a_rrf = a * 2*cos(w_rrf*t + phi_rrf)
946  PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
947  Double_t wRRF = globalBlock->GetRRFFreq("Mc");
948  Double_t phaseRRF = globalBlock->GetRRFPhase()*TMath::TwoPi()/180.0;
949 
950  Double_t startTime = fTimeResolution * static_cast<Double_t>(fgbOffset);
951  Double_t time=0.0;
952  for (UInt_t i=0; i<asym.size(); i++) {
953  time = startTime + i*fTimeResolution;
954  asym[i] *= 2.0*cos(wRRF*time+phaseRRF);
955  }
956 
957  // 3rd: rrf packing
958  PDoubleVector asymRRF;
959  asymVal = 0.0;
960  for (UInt_t i=0; i<asym.size(); i++) {
961  if ((i+1) % fRRFPacking == 0) {
962  asymRRF.push_back(asymVal/fRRFPacking);
963  asymVal = 0.0;
964  }
965  asymVal += asym[i];
966  }
967 
968  // 4th: rrf packing error
969  PDoubleVector asymRRFErr;
970  asymValErr = 0.0;
971  for (UInt_t i=0; i<asymErr.size(); i++) {
972  if ((i+1) % fRRFPacking == 0) {
973  asymRRFErr.push_back(sqrt(2.0*asymValErr)/fRRFPacking); // factor of two is needed due to the rescaling
974  asymValErr = 0.0;
975  }
976  asymValErr += asymErr[i]*asymErr[i];
977  }
978 
979  fData.SetDataTimeStart(startTime+fTimeResolution*(static_cast<Double_t>(fRRFPacking-1)/2.0));
980  fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(fRRFPacking));
981 
982  for (UInt_t i=0; i<asymRRF.size(); i++) {
983  fData.AppendValue(asymRRF[i]);
984  fData.AppendErrorValue(asymRRFErr[i]);
985  }
986 
987  CalcNoOfFitBins();
988 
989  // clean up
990  fForward.clear();
991  fForwardErr.clear();
992  fBackward.clear();
993  fBackwardErr.clear();
994  asym.clear();
995  asymErr.clear();
996  asymRRF.clear();
997  asymRRFErr.clear();
998 
999  return true;
1000 }
1001 
1002 //--------------------------------------------------------------------------
1003 // PrepareViewData (protected)
1004 //--------------------------------------------------------------------------
1021 Bool_t PRunAsymmetryRRF::PrepareViewData(PRawRunData* runData, UInt_t histoNo[2])
1022 {
1023  // feed the parameter vector
1024  std::vector<Double_t> par;
1025  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1026  for (UInt_t i=0; i<paramList->size(); i++)
1027  par.push_back((*paramList)[i].fValue);
1028 
1029  // first get start data, end data, and t0
1030  Int_t start[2] = {fGoodBins[0], fGoodBins[2]};
1031  Int_t end[2] = {fGoodBins[1], fGoodBins[3]};
1032  Int_t t0[2] = {static_cast<Int_t>(fT0s[0]), static_cast<Int_t>(fT0s[1])};
1033 
1034  // check if the data ranges and t0's between forward/backward are compatible
1035  Int_t fgb[2];
1036  if (start[0]-t0[0] != start[1]-t0[1]) { // wrong fgb aligning
1037  if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
1038  fgb[0] = start[0];
1039  fgb[1] = t0[1] + start[0]-t0[0];
1040  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **WARNING** needed to shift backward fgb from ";
1041  std::cerr << start[1] << " to " << fgb[1] << std::endl;
1042  } else {
1043  fgb[0] = t0[0] + start[1]-t0[1];
1044  fgb[1] = start[1];
1045  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **WARNING** needed to shift forward fgb from ";
1046  std::cerr << start[0] << " to " << fgb[0] << std::endl;
1047  }
1048  } else { // fgb aligning is correct
1049  fgb[0] = start[0];
1050  fgb[1] = start[1];
1051  }
1052  start[0] = fgb[0];
1053  start[1] = fgb[1];
1054 
1055  // make sure that there are equal number of bins in forward and backward
1056  UInt_t noOfBins0 = runData->GetDataBin(histoNo[0])->size()-start[0];
1057  UInt_t noOfBins1 = runData->GetDataBin(histoNo[1])->size()-start[1];
1058  if (noOfBins0 > noOfBins1)
1059  noOfBins0 = noOfBins1;
1060  end[0] = start[0] + noOfBins0;
1061  end[1] = start[1] + noOfBins0;
1062 
1063  // check if start, end, and t0 make any sense
1064  for (UInt_t i=0; i<2; i++) {
1065  // 1st check if start is within proper bounds
1066  if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1067  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
1068  std::cerr << std::endl;
1069  return false;
1070  }
1071  // 2nd check if end is within proper bounds
1072  if ((end[i] < 0) || (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1073  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
1074  std::cerr << std::endl;
1075  return false;
1076  }
1077  // 3rd check if t0 is within proper bounds
1078  if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1079  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **ERROR** t0 data bin doesn't make any sense!";
1080  std::cerr << std::endl;
1081  return false;
1082  }
1083  }
1084 
1085  // check if forward and backward histo have the same size, otherwise take the minimum size
1086  UInt_t noOfBins = fForward.size();
1087  if (noOfBins > fBackward.size()) {
1088  noOfBins = fBackward.size();
1089  }
1090 
1091  // form asymmetry including error propagation
1092  Double_t asym, error;
1093  Double_t f, b, ef, eb, alpha = 1.0, beta = 1.0;
1094 
1095  // get the proper alpha and beta
1096  switch (fAlphaBetaTag) {
1097  case 1: // alpha == 1, beta == 1
1098  alpha = 1.0;
1099  beta = 1.0;
1100  break;
1101  case 2: // alpha != 1, beta == 1
1102  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1103  alpha = par[fRunInfo->GetAlphaParamNo()-1];
1104  } else { // alpha is function
1105  // get function number
1106  UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1107  // evaluate function
1108  alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1109  }
1110  beta = 1.0;
1111  break;
1112  case 3: // alpha == 1, beta != 1
1113  alpha = 1.0;
1114  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1115  beta = par[fRunInfo->GetBetaParamNo()-1];
1116  } else { // beta is a function
1117  // get function number
1118  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1119  // evaluate function
1120  beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1121  }
1122  break;
1123  case 4: // alpha != 1, beta != 1
1124  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1125  alpha = par[fRunInfo->GetAlphaParamNo()-1];
1126  } else { // alpha is function
1127  // get function number
1128  UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1129  // evaluate function
1130  alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1131  }
1132  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1133  beta = par[fRunInfo->GetBetaParamNo()-1];
1134  } else { // beta is a function
1135  // get function number
1136  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1137  // evaluate function
1138  beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1139  }
1140  break;
1141  default:
1142  break;
1143  }
1144 
1145  PDoubleVector asymVec, asymErr;
1146  Int_t dtBin = start[1]-start[0];
1147  for (Int_t i=start[0]; i<end[0]; i++) {
1148  // to make the formulae more readable
1149  f = fForward[i];
1150  b = fBackward[i+dtBin];
1151  ef = fForwardErr[i];
1152  eb = fBackwardErr[i+dtBin];
1153  // check that there are indeed bins
1154  if (f+b != 0.0)
1155  asym = (alpha*f-b) / (alpha*beta*f+b);
1156  else
1157  asym = 0.0;
1158  asymVec.push_back(asym);
1159  // calculate the error
1160  if (f+b != 0.0)
1161  error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1162  else
1163  error = 1.0;
1164  asymErr.push_back(error);
1165  }
1166 
1167  // RRF transform
1168  // a_rrf = a * 2*cos(w_rrf*t + phi_rrf)
1169  PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
1170  Double_t wRRF = globalBlock->GetRRFFreq("Mc");
1171  Double_t phaseRRF = globalBlock->GetRRFPhase()*TMath::TwoPi()/180.0;
1172  Double_t startTime=fTimeResolution*(static_cast<Double_t>(start[0])-t0[0]);
1173  Double_t time = 0.0;
1174  for (UInt_t i=0; i<asymVec.size(); i++) {
1175  time = startTime + i * fTimeResolution;
1176  asymVec[i] *= 2.0*cos(wRRF*time+phaseRRF); // factor of 2 needed to keep the asymmetry
1177  }
1178 
1179  Double_t dval = 0.0;
1180  PDoubleVector asymRRF;
1181  for (UInt_t i=0; i<asymVec.size(); i++) {
1182  if ((i+1) % fRRFPacking == 0) {
1183  asymRRF.push_back(dval/fRRFPacking);
1184  dval = 0.0;
1185  }
1186  dval += asymVec[i];
1187  }
1188 
1189  // RRF packing error
1190  PDoubleVector asymRRFErr;
1191  dval = 0.0;
1192  for (UInt_t i=0; i<asymErr.size(); i++) {
1193  if ((i+1) % fRRFPacking == 0) {
1194  asymRRFErr.push_back(sqrt(2.0*dval)/fRRFPacking); // factor of two is needed due to the rescaling
1195  dval = 0.0;
1196  }
1197  dval += asymErr[i]*asymErr[i];
1198  }
1199 
1200  fData.SetDataTimeStart(startTime+fTimeResolution*(static_cast<Double_t>(fRRFPacking-1)/2.0));
1201  fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(fRRFPacking));
1202 
1203  for (UInt_t i=0; i<asymRRF.size(); i++) {
1204  fData.AppendValue(asymRRF[i]);
1205  fData.AppendErrorValue(asymRRFErr[i]);
1206  }
1207 
1208  CalcNoOfFitBins();
1209 
1210  // clean up
1211  fForward.clear();
1212  fForwardErr.clear();
1213  fBackward.clear();
1214  fBackwardErr.clear();
1215  asymVec.clear();
1216  asymErr.clear();
1217  asymRRF.clear();
1218  asymRRFErr.clear();
1219 
1220  // fill theory vector for kView
1221  // calculate functions
1222  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1224  }
1225 
1226  // calculate theory
1227  UInt_t size = runData->GetDataBin(histoNo[0])->size();
1228  Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1230  if (fTheoAsData) { // calculate theory only at the data points
1232  } else {
1233  // finer binning for the theory (8 times as many points = factor)
1234  size *= factor;
1235  fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1236  }
1237 
1238  for (UInt_t i=0; i<size; i++) {
1239  time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1240  dval = fTheory->Func(time, par, fFuncValues);
1241  if (fabs(dval) > 10.0) { // dirty hack needs to be fixed!!
1242  dval = 0.0;
1243  }
1244  fData.AppendTheoryValue(dval);
1245  }
1246 
1247  // clean up
1248  par.clear();
1249 
1250  return true;
1251 }
1252 
1253 //--------------------------------------------------------------------------
1254 // GetProperT0 (private)
1255 //--------------------------------------------------------------------------
1274 Bool_t PRunAsymmetryRRF::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHistoNo, PUIntVector &backwardHistoNo)
1275 {
1276  // feed all T0's
1277  // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1278  fT0s.clear();
1279  // this strange fT0 size estimate is needed in case #forw histos != #back histos
1280  size_t size = 2*forwardHistoNo.size();
1281  if (backwardHistoNo.size() > forwardHistoNo.size())
1282  size = 2*backwardHistoNo.size();
1283  fT0s.resize(size);
1284  for (UInt_t i=0; i<fT0s.size(); i++) {
1285  fT0s[i] = -1.0;
1286  }
1287 
1288  // fill in the T0's from the msr-file (if present)
1289  for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1290  fT0s[i] = fRunInfo->GetT0Bin(i);
1291  }
1292 
1293  // fill in the missing T0's from the GLOBAL block section (if present)
1294  for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1295  if (fT0s[i] == -1.0) { // i.e. not given in the RUN block section
1296  fT0s[i] = globalBlock->GetT0Bin(i);
1297  }
1298  }
1299 
1300  // fill in the missing T0's from the data file, if not already present in the msr-file
1301  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1302  if (fT0s[2*i] == -1.0) // i.e. not present in the msr-file, try the data file
1303  if (runData->GetT0Bin(forwardHistoNo[i]) > 0.0) {
1304  fT0s[2*i] = runData->GetT0Bin(forwardHistoNo[i]);
1305  fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1306  }
1307  }
1308  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1309  if (fT0s[2*i+1] == -1.0) // i.e. not present in the msr-file, try the data file
1310  if (runData->GetT0Bin(backwardHistoNo[i]) > 0.0) {
1311  fT0s[2*i+1] = runData->GetT0Bin(backwardHistoNo[i]);
1312  fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1313  }
1314  }
1315 
1316  // fill in the T0's gaps, i.e. in case the T0's are NEITHER in the msr-file and NOR in the data file
1317  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1318  if (fT0s[2*i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1319  fT0s[2*i] = runData->GetT0BinEstimated(forwardHistoNo[i]);
1320  fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1321 
1322  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1323  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1324  std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(forwardHistoNo[i]);
1325  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1326  std::cerr << std::endl;
1327  }
1328  }
1329  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1330  if (fT0s[2*i+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1331  fT0s[2*i+1] = runData->GetT0BinEstimated(backwardHistoNo[i]);
1332  fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1333 
1334  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1335  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1336  std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[i]);
1337  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1338  std::cerr << std::endl;
1339  }
1340  }
1341 
1342  // check if t0 is within proper bounds
1343  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1344  if ((fT0s[2*i] < 0) || (fT0s[2*i] > static_cast<Int_t>(runData->GetDataBin(forwardHistoNo[i])->size()))) {
1345  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **ERROR** t0 data bin (" << fT0s[2*i] << ") doesn't make any sense!";
1346  std::cerr << std::endl << ">> forwardHistoNo " << forwardHistoNo[i];
1347  std::cerr << std::endl;
1348  return false;
1349  }
1350  }
1351  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1352  if ((fT0s[2*i+1] < 0) || (fT0s[2*i+1] > static_cast<Int_t>(runData->GetDataBin(backwardHistoNo[i])->size()))) {
1353  std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **ERROR** t0 data bin (" << fT0s[2*i+1] << ") doesn't make any sense!";
1354  std::cerr << std::endl << ">> backwardHistoNo " << backwardHistoNo[i];
1355  std::cerr << std::endl;
1356  return false;
1357  }
1358  }
1359 
1360  // check if addrun's are present, and if yes add the necessary t0's
1361  if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1362  PRawRunData *addRunData;
1363  fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1364  for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1365  // get run to be added to the main one
1366  addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
1367  if (addRunData == nullptr) { // couldn't get run
1368  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1369  std::cerr << std::endl;
1370  return false;
1371  }
1372 
1373  // feed all T0's
1374  // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1375  fAddT0s[i-1].clear();
1376  fAddT0s[i-1].resize(2*forwardHistoNo.size());
1377  for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1378  fAddT0s[i-1][j] = -1.0;
1379  }
1380 
1381  // fill in the T0's from the msr-file (if present)
1382  for (Int_t j=0; j<fRunInfo->GetAddT0BinSize(i-1); j++) {
1383  fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1, j);
1384  }
1385 
1386  // fill in the T0's from the data file, if not already present in the msr-file
1387  for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1388  if (fAddT0s[i-1][2*j] == -1.0) // i.e. not present in the msr-file, try the data file
1389  if (addRunData->GetT0Bin(forwardHistoNo[j]) > 0.0) {
1390  fAddT0s[i-1][2*j] = addRunData->GetT0Bin(forwardHistoNo[j]);
1391  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1392  }
1393  }
1394  for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1395  if (fAddT0s[i-1][2*j+1] == -1.0) // i.e. not present in the msr-file, try the data file
1396  if (addRunData->GetT0Bin(backwardHistoNo[j]) > 0.0) {
1397  fAddT0s[i-1][2*j+1] = addRunData->GetT0Bin(backwardHistoNo[j]);
1398  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1399  }
1400  }
1401 
1402  // 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
1403  for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1404  if (fAddT0s[i-1][2*j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1405  fAddT0s[i-1][2*j] = addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1406  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1407 
1408  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1409  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1410  std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1411  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1412  std::cerr << std::endl;
1413  }
1414  }
1415  for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1416  if (fAddT0s[i-1][2*j+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1417  fAddT0s[i-1][2*j+1] = addRunData->GetT0BinEstimated(backwardHistoNo[j]);
1418  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1419 
1420  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1421  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1422  std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[j]);
1423  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1424  std::cerr << std::endl;
1425  }
1426  }
1427  }
1428  }
1429 
1430  return true;
1431 }
1432 
1433 //--------------------------------------------------------------------------
1434 // GetProperDataRange (private)
1435 //--------------------------------------------------------------------------
1449 Bool_t PRunAsymmetryRRF::GetProperDataRange(PRawRunData* runData, UInt_t histoNo[2])
1450 {
1451  // first get start/end data
1452  Int_t start[2] = {fRunInfo->GetDataRange(0), fRunInfo->GetDataRange(2)};
1453  Int_t end[2] = {fRunInfo->GetDataRange(1), fRunInfo->GetDataRange(3)};
1454  // check if data range has been provided in the RUN block. If not, try the GLOBAL block
1455  if (start[0] == -1) {
1456  start[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1457  }
1458  if (start[1] == -1) {
1459  start[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(2);
1460  }
1461  if (end[0] == -1) {
1462  end[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1463  }
1464  if (end[1] == -1) {
1465  end[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(3);
1466  }
1467 
1468  Double_t t0[2] = {fT0s[0], fT0s[1]};
1469  Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution); // needed in case first good bin is not given, default = 10ns
1470 
1471  // check if data range has been provided, and if not try to estimate them
1472  if (start[0] < 0) {
1473  start[0] = static_cast<Int_t>(t0[0])+offset;
1474  fRunInfo->SetDataRange(start[0], 0);
1475  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[0] << ".";
1476  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1477  std::cerr << std::endl;
1478  }
1479  if (start[1] < 0) {
1480  start[1] = static_cast<Int_t>(t0[1])+offset;
1481  fRunInfo->SetDataRange(start[1], 2);
1482  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[1] << ".";
1483  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1484  std::cerr << std::endl;
1485  }
1486  if (end[0] < 0) {
1487  end[0] = runData->GetDataBin(histoNo[0])->size();
1488  fRunInfo->SetDataRange(end[0], 1);
1489  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range end = " << end[0] << ".";
1490  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1491  std::cerr << std::endl;
1492  }
1493  if (end[1] < 0) {
1494  end[1] = runData->GetDataBin(histoNo[1])->size();
1495  fRunInfo->SetDataRange(end[1], 3);
1496  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range end = " << end[1] << ".";
1497  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1498  std::cerr << std::endl;
1499  }
1500 
1501  // check if start, end, and t0 make any sense
1502  // 1st check if start and end are in proper order
1503  for (UInt_t i=0; i<2; i++) {
1504  if (end[i] < start[i]) { // need to swap them
1505  Int_t keep = end[i];
1506  end[i] = start[i];
1507  start[i] = keep;
1508  }
1509  // 2nd check if start is within proper bounds
1510  if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1511  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **ERROR** start data bin doesn't make any sense!";
1512  std::cerr << std::endl;
1513  return false;
1514  }
1515  // 3rd check if end is within proper bounds
1516  if (end[i] < 0) {
1517  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **ERROR** end data bin (" << end[i] << ") doesn't make any sense!";
1518  std::cerr << std::endl;
1519  return false;
1520  }
1521  if (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())) {
1522  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** end data bin (" << end[i] << ") > histo length (" << (Int_t)runData->GetDataBin(histoNo[i])->size() << ").";
1523  std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
1524  std::cerr << std::endl;
1525  end[i] = static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())-1;
1526  }
1527  // 4th check if t0 is within proper bounds
1528  if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1529  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **ERROR** t0 data bin doesn't make any sense!";
1530  std::cerr << std::endl;
1531  return false;
1532  }
1533  }
1534 
1535  // check that start-t0 is the same for forward as for backward, otherwise take max(start[i]-t0[i])
1536  if (fabs(static_cast<Double_t>(start[0])-t0[0]) > fabs(static_cast<Double_t>(start[1])-t0[1])){
1537  start[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(start[0]) - t0[0]);
1538  end[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(end[0]) - t0[0]);
1539  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange **WARNING** needed to shift backward data range.";
1540  std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(2) << ", " << fRunInfo->GetDataRange(3);
1541  std::cerr << std::endl << ">> used : " << start[1] << ", " << end[1];
1542  std::cerr << std::endl;
1543  }
1544  if (fabs(static_cast<Double_t>(start[0])-t0[0]) < fabs(static_cast<Double_t>(start[1])-t0[1])){
1545  start[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(start[1]) - t0[1]);
1546  end[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(end[1]) - t0[1]);
1547  std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange **WARNING** needed to shift forward data range.";
1548  std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(0) << ", " << fRunInfo->GetDataRange(1);
1549  std::cerr << std::endl << ">> used : " << start[0] << ", " << end[0];
1550  std::cerr << std::endl;
1551  }
1552 
1553  // keep good bins for potential latter use
1554  fGoodBins[0] = start[0];
1555  fGoodBins[1] = end[0];
1556  fGoodBins[2] = start[1];
1557  fGoodBins[3] = end[1];
1558 
1559  return true;
1560 }
1561 
1562 //--------------------------------------------------------------------------
1563 // GetProperFitRange (private)
1564 //--------------------------------------------------------------------------
1578 {
1579  // set fit start/end time; first check RUN Block
1582  // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1583  if (fRunInfo->IsFitRangeInBin()) {
1584  fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1585  fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1586  // write these times back into the data structure. This way it is available when writting the log-file
1589  }
1590  if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
1591  fFitStartTime = globalBlock->GetFitRange(0);
1592  fFitEndTime = globalBlock->GetFitRange(1);
1593  // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1594  if (globalBlock->IsFitRangeInBin()) {
1595  fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1596  fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1597  // write these times back into the data structure. This way it is available when writting the log-file
1598  globalBlock->SetFitRange(fFitStartTime, 0);
1599  globalBlock->SetFitRange(fFitEndTime, 1);
1600  }
1601  }
1603  fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
1604  fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
1605  std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
1606  std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
1607  }
1608 }
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
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 Bool_t PrepareFitData()
Int_t fRRFPacking
RRF packing for this particular run. Given in the GLOBAL-block.
virtual void SetDataRange(Int_t ival, Int_t idx)
Definition: PMusr.cpp:1672
virtual const PDoubleVector * GetError()
Definition: PMusr.h:249
Definition: PMusr.h:220
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
Int_t fStartTimeBin
bin at which the fit starts
virtual const Double_t GetTimeResolution()
Definition: PMusr.h:429
virtual Bool_t GetProperDataRange(PRawRunData *runData, UInt_t histoNo[2])
PRunDataHandler * fRawData
holds the raw run data
Definition: PRunBase.h:73
virtual Bool_t PrepareViewData(PRawRunData *runData, UInt_t histoNo[2])
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
PMsrHandler * fMsrInfo
msr-file handler
Definition: PRunBase.h:71
virtual Bool_t PrepareData()
virtual Bool_t IsFitRangeInBin()
Definition: PMusr.h:665
virtual const PDoublePairVector * GetTemperature() const
Definition: PMusr.h:422
Double_t fTimeResolution
time resolution in (us)
Definition: PRunBase.h:76
virtual void CalcTheory()
PMetaData fMetaData
keeps the meta data from the data file like field, temperature, energy, ...
Definition: PRunBase.h:77
virtual PRawRunData * GetRunData(const TString &runName)
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
virtual Double_t GetAddT0Bin(UInt_t addRunIdx, UInt_t histoIdx)
Definition: PMusr.cpp:1761
virtual Int_t GetRRFPacking()
Definition: PMusr.h:580
virtual UInt_t GetT0BinSize()
Definition: PMusr.h:660
PDoubleVector fTemp
temperature(s) in (K)
Definition: PMusr.h:229
virtual void AppendErrorValue(Double_t dval)
Definition: PMusr.h:260
PDoubleVector fBackwardErr
backward histo errors
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 Double_t GetRRFFreq(const char *unit)
Definition: PMusr.cpp:762
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 CalcChiSquare(const std::vector< Double_t > &par)
UInt_t fNoOfFitBins
number of bins to be be fitted
virtual Int_t GetBackwardHistoNo(UInt_t idx=0)
Definition: PMusr.cpp:1444
virtual void SetT0Bin(Double_t dval, Int_t idx=-1)
Definition: PMusr.cpp:1712
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHisto, PUIntVector &backwardHistoNo)
virtual void AppendTheoryValue(Double_t dval)
Definition: PMusr.h:262
virtual PMsrGlobalBlock * GetMsrGlobal()
Definition: PMsrHandler.h:62
virtual Int_t GetAddT0BinSize(UInt_t addRunIdx)
Definition: PMusr.cpp:1737
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
virtual const PDoubleVector * GetValue()
Definition: PMusr.h:248
virtual UInt_t GetNoOfFitBins()
virtual const Double_t GetField()
Definition: PMusr.h:420
virtual void SetAddT0Bin(Double_t dval, UInt_t addRunIdx, UInt_t histoNoIdx)
Definition: PMusr.cpp:1788
virtual Double_t GetBkgFix(UInt_t idx)
Definition: PMusr.cpp:1572
virtual UInt_t GetFuncNo(Int_t idx)
Definition: PMsrHandler.h:94
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
Double_t fEnergy
energy in (keV)
Definition: PMusr.h:228
virtual Int_t GetBkgRange(UInt_t idx)
Definition: PMusr.cpp:1612
UInt_t fAlphaBetaTag
; ; ; .
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
Int_t fGoodBins[4]
keep first/last good bins. 0=fgb, 1=lgb (forward); 2=fgb, 3=lgb (backward)
virtual Double_t GetRRFPhase()
Definition: PMusr.h:579
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
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 void GetProperFitRange(PMsrGlobalBlock *globalBlock)
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
PDoubleVector fForward
forward histo data
#define ACCEL_PERIOD_PSI
Definition: PMusr.h:81
virtual Int_t GetAlphaParamNo()
Definition: PMusr.h:644
virtual void SetFitRange(Double_t dval, UInt_t idx)
Definition: PMusr.cpp:1068
Double_t fFitStartTime
fit start time
Definition: PRunBase.h:81
Int_t fEndTimeBin
bin at which the fit ends
virtual Int_t GetBetaParamNo()
Definition: PMusr.h:645
virtual PMsrParamList * GetMsrParamList()
Definition: PMsrHandler.h:59
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
return status
virtual const Double_t GetT0BinEstimated(const UInt_t histoNo)
Definition: PMusr.h:432
PDoubleVector fBackward
backward histo data
std::unique_ptr< PTheory > fTheory
theory needed to calculate chi-square
Definition: PRunBase.h:85
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 Double_t GetDataTimeStep()
Definition: PMusr.h:243
Int_t fRunNo
number of the run within the msr-file
Definition: PRunBase.h:70
PDoubleVector fForwardErr
forward histo errors
virtual void CalcNoOfFitBins()
virtual TString * GetRunName(UInt_t idx=0)
Definition: PMusr.cpp:1232
virtual Double_t GetFitRange(UInt_t idx)
Definition: PMusr.cpp:1811
virtual void SetFitRangeBin(const TString fitRange)
Definition: PMusr.h:220
PDoubleVector fFuncValues
is keeping the values of the functions from the FUNCTIONS block
Definition: PRunBase.h:84
virtual UInt_t GetBackwardHistoNoSize()
Definition: PMusr.h:654
Double_t fField
field in (G)
Definition: PMusr.h:227
#define ACCEL_PERIOD_RAL
Definition: PMusr.h:83