musrfit  1.9.2
PRunAsymmetry.cpp
Go to the documentation of this file.
1 /***************************************************************************
2 
3  PRunAsymmetry.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 
42 #include <TString.h>
43 #include <TObjArray.h>
44 #include <TObjString.h>
45 
46 #include "PMusr.h"
47 #include "PRunAsymmetry.h"
48 
49 //--------------------------------------------------------------------------
50 // Constructor
51 //--------------------------------------------------------------------------
56 {
57  fNoOfFitBins = 0;
58  fPacking = -1;
59  fTheoAsData = false;
60 
61  // the 2 following variables are need in case fit range is given in bins, and since
62  // the fit range can be changed in the command block, these variables need to be accessible
63  fGoodBins[0] = -1;
64  fGoodBins[1] = -1;
65 
66  fStartTimeBin = -1;
67  fEndTimeBin = -1;
68 }
69 
70 //--------------------------------------------------------------------------
71 // Constructor
72 //--------------------------------------------------------------------------
81 PRunAsymmetry::PRunAsymmetry(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
82  PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
83 {
84  // the 2 following variables are need in case fit range is given in bins, and since
85  // the fit range can be changed in the command block, these variables need to be accessible
86  fGoodBins[0] = -1;
87  fGoodBins[1] = -1;
88 
89  fStartTimeBin = -1;
90  fEndTimeBin = -1;
91 
93  if (fPacking == -1) { // i.e. packing is NOT given in the RUN-block, it must be given in the GLOBAL-block
95  }
96  if (fPacking == -1) { // this should NOT happen, somethin is severely wrong
97  std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **SEVERE ERROR**: Couldn't find any packing information!";
98  std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
99  std::cerr << std::endl;
100  fValid = false;
101  return;
102  }
103 
104  // check if alpha and/or beta is fixed --------------------
105 
106  PMsrParamList *param = msrInfo->GetMsrParamList();
107 
108  // check if alpha is given
109  if (fRunInfo->GetAlphaParamNo() == -1) { // no alpha given
110  std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** no alpha parameter given! This is needed for an asymmetry fit!";
111  std::cerr << std::endl;
112  fValid = false;
113  return;
114  }
115  // check if alpha parameter is within proper bounds
116  if (fRunInfo->GetAlphaParamNo() >= MSR_PARAM_FUN_OFFSET) { // alpha seems to be a function
119  std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** alpha parameter is a function with no = " << fRunInfo->GetAlphaParamNo();
120  std::cerr << std::endl << ">> This is out of bound, since there are only " << msrInfo->GetNoOfFuncs() << " functions.";
121  std::cerr << std::endl;
122  fValid = false;
123  return;
124  }
125  } else if ((fRunInfo->GetAlphaParamNo() < 0) || (fRunInfo->GetAlphaParamNo() > static_cast<Int_t>(param->size()))) {
126  std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** alpha parameter no = " << fRunInfo->GetAlphaParamNo();
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  }
132  // check if alpha is fixed
133  Bool_t alphaFixedToOne = false;
134  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is parameter
135  if (((*param)[fRunInfo->GetAlphaParamNo()-1].fStep == 0.0) &&
136  ((*param)[fRunInfo->GetAlphaParamNo()-1].fValue == 1.0))
137  alphaFixedToOne = true;
138  }
139  // check if beta is given
140  Bool_t betaFixedToOne = false;
141  if (fRunInfo->GetBetaParamNo() == -1) { // no beta given hence assuming beta == 1
142  betaFixedToOne = true;
143  } else if ((fRunInfo->GetBetaParamNo() < 0) || (fRunInfo->GetBetaParamNo() > static_cast<Int_t>(param->size()))) { // check if beta parameter is within proper bounds
144  std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** beta parameter no = " << fRunInfo->GetBetaParamNo();
145  std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
146  std::cerr << std::endl;
147  fValid = false;
148  return;
149  } else { // check if beta is fixed
150  if (((*param)[fRunInfo->GetBetaParamNo()-1].fStep == 0.0) &&
151  ((*param)[fRunInfo->GetBetaParamNo()-1].fValue == 1.0))
152  betaFixedToOne = true;
153  }
154 
155  // set fAlphaBetaTag
156  if (alphaFixedToOne && betaFixedToOne) // alpha == 1, beta == 1
157  fAlphaBetaTag = 1;
158  else if (!alphaFixedToOne && betaFixedToOne) // alpha != 1, beta == 1
159  fAlphaBetaTag = 2;
160  else if (alphaFixedToOne && !betaFixedToOne) // alpha == 1, beta != 1
161  fAlphaBetaTag = 3;
162  else
163  fAlphaBetaTag = 4;
164 
165  // calculate fData
166  if (!PrepareData())
167  fValid = false;
168 }
169 
170 //--------------------------------------------------------------------------
171 // Destructor
172 //--------------------------------------------------------------------------
177 {
178  fForward.clear();
179  fForwardErr.clear();
180  fBackward.clear();
181  fBackwardErr.clear();
182 }
183 
184 //--------------------------------------------------------------------------
185 // CalcChiSquare (public)
186 //--------------------------------------------------------------------------
195 Double_t PRunAsymmetry::CalcChiSquare(const std::vector<Double_t>& par)
196 {
197  Double_t chisq = 0.0;
198  Double_t diff = 0.0;
199  Double_t asymFcnValue = 0.0;
200  Double_t a, b, f;
201 
202  // calculate functions
203  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
205  }
206 
207  // calculate chi square
208  Double_t time(1.0);
209  Int_t i;
210 
211  // determine alpha/beta
212  switch (fAlphaBetaTag) {
213  case 1: // alpha == 1, beta == 1
214  a = 1.0;
215  b = 1.0;
216  break;
217  case 2: // alpha != 1, beta == 1
218  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
219  a = par[fRunInfo->GetAlphaParamNo()-1];
220  } else { // alpha is function
221  // get function number
223  // evaluate function
224  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
225  }
226  b = 1.0;
227  break;
228  case 3: // alpha == 1, beta != 1
229  a = 1.0;
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  case 4: // alpha != 1, beta != 1
240  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
241  a = par[fRunInfo->GetAlphaParamNo()-1];
242  } else { // alpha is function
243  // get function number
245  // evaluate function
246  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
247  }
248  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
249  b = par[fRunInfo->GetBetaParamNo()-1];
250  } else { // beta is a function
251  // get function number
252  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
253  // evaluate function
254  b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
255  }
256  break;
257  default:
258  a = 1.0;
259  b = 1.0;
260  break;
261  }
262 
263  // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
264  // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
265  // for a given set of parameters---which should be done outside of the parallelized loop.
266  // For all other functions it means a tiny and acceptable overhead.
267  asymFcnValue = fTheory->Func(time, par, fFuncValues);
268 
269  #ifdef HAVE_GOMP
270  Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
271  if (chunk < 10)
272  chunk = 10;
273  #pragma omp parallel for default(shared) private(i,time,diff,asymFcnValue,f) schedule(dynamic,chunk) reduction(+:chisq)
274  #endif
275  for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
276  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
277  f = fTheory->Func(time, par, fFuncValues);
278  asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
279  diff = fData.GetValue()->at(i) - asymFcnValue;
280  chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
281  }
282 
283  return chisq;
284 }
285 
286 //--------------------------------------------------------------------------
287 // CalcChiSquareExpected (public)
288 //--------------------------------------------------------------------------
297 Double_t PRunAsymmetry::CalcChiSquareExpected(const std::vector<Double_t>& par)
298 {
299  return 0.0;
300 }
301 
302 //--------------------------------------------------------------------------
303 // CalcMaxLikelihood (public)
304 //--------------------------------------------------------------------------
310 Double_t PRunAsymmetry::CalcMaxLikelihood(const std::vector<Double_t>& par)
311 {
312  std::cout << std::endl << "PRunAsymmetry::CalcMaxLikelihood(): not implemented yet ..." << std::endl;
313 
314  return 1.0;
315 }
316 
317 //--------------------------------------------------------------------------
318 // GetNoOfFitBins (public)
319 //--------------------------------------------------------------------------
326 {
327  CalcNoOfFitBins();
328 
329  return fNoOfFitBins;
330 }
331 
332 //--------------------------------------------------------------------------
333 // SetFitRangeBin (public)
334 //--------------------------------------------------------------------------
346 void PRunAsymmetry::SetFitRangeBin(const TString fitRange)
347 {
348  TObjArray *tok = nullptr;
349  TObjString *ostr = nullptr;
350  TString str;
351  Ssiz_t idx = -1;
352  Int_t offset = 0;
353 
354  tok = fitRange.Tokenize(" \t");
355 
356  if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
357  // handle fgb+n0 entry
358  ostr = dynamic_cast<TObjString*>(tok->At(1));
359  str = ostr->GetString();
360  // check if there is an offset present
361  idx = str.First("+");
362  if (idx != -1) { // offset present
363  str.Remove(0, idx+1);
364  if (str.IsFloat()) // if str is a valid number, convert is to an integer
365  offset = str.Atoi();
366  }
367  fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
368 
369  // handle lgb-n1 entry
370  ostr = dynamic_cast<TObjString*>(tok->At(2));
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  fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
380  } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
381  Int_t pos = 2*(fRunNo+1)-1;
382 
383  if (pos + 1 >= tok->GetEntries()) {
384  std::cerr << std::endl << ">> PRunAsymmetry::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
385  std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
386  } else {
387  // handle fgb+n0 entry
388  ostr = static_cast<TObjString*>(tok->At(pos));
389  str = ostr->GetString();
390  // check if there is an offset present
391  idx = str.First("+");
392  if (idx != -1) { // offset present
393  str.Remove(0, idx+1);
394  if (str.IsFloat()) // if str is a valid number, convert is to an integer
395  offset = str.Atoi();
396  }
397  fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
398 
399  // handle lgb-n1 entry
400  ostr = static_cast<TObjString*>(tok->At(pos+1));
401  str = ostr->GetString();
402  // check if there is an offset present
403  idx = str.First("-");
404  if (idx != -1) { // offset present
405  str.Remove(0, idx+1);
406  if (str.IsFloat()) // if str is a valid number, convert is to an integer
407  offset = str.Atoi();
408  }
409  fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
410  }
411  } else { // error
412  std::cerr << std::endl << ">> PRunAsymmetry::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
413  std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
414  }
415 
416  // clean up
417  if (tok) {
418  delete tok;
419  }
420 }
421 
422 //--------------------------------------------------------------------------
423 // CalcNoOfFitBins (public)
424 //--------------------------------------------------------------------------
429 {
430  // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
431  fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
432  if (fStartTimeBin < 0)
433  fStartTimeBin = 0;
434  fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
435  if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
436  fEndTimeBin = fData.GetValue()->size();
437 
439  fNoOfFitBins = static_cast<UInt_t>(fEndTimeBin - fStartTimeBin);
440  else
441  fNoOfFitBins = 0;
442 }
443 
444 //--------------------------------------------------------------------------
445 // CalcTheory (protected)
446 //--------------------------------------------------------------------------
451 {
452  // feed the parameter vector
453  std::vector<Double_t> par;
454  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
455  for (UInt_t i=0; i<paramList->size(); i++)
456  par.push_back((*paramList)[i].fValue);
457 
458  // calculate functions
459  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
461  }
462 
463  // calculate asymmetry
464  Double_t asymFcnValue = 0.0;
465  Double_t a, b, f;
466  Double_t time;
467  for (UInt_t i=0; i<fData.GetValue()->size(); i++) {
468  time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
469  switch (fAlphaBetaTag) {
470  case 1: // alpha == 1, beta == 1
471  asymFcnValue = fTheory->Func(time, par, fFuncValues);
472  break;
473  case 2: // alpha != 1, beta == 1
474  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
475  a = par[fRunInfo->GetAlphaParamNo()-1];
476  } else { // alpha is function
477  // get function number
479  // evaluate function
480  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
481  }
482  f = fTheory->Func(time, par, fFuncValues);
483  asymFcnValue = (f*(a+1.0)-(a-1.0))/((a+1.0)-f*(a-1.0));
484  break;
485  case 3: // alpha == 1, beta != 1
486  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
487  b = par[fRunInfo->GetBetaParamNo()-1];
488  } else { // beta is a function
489  // get function number
490  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
491  // evaluate function
492  b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
493  }
494  f = fTheory->Func(time, par, fFuncValues);
495  asymFcnValue = f*(b+1.0)/(2.0-f*(b-1.0));
496  break;
497  case 4: // alpha != 1, beta != 1
498  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
499  a = par[fRunInfo->GetAlphaParamNo()-1];
500  } else { // alpha is function
501  // get function number
503  // evaluate function
504  a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
505  }
506  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
507  b = par[fRunInfo->GetBetaParamNo()-1];
508  } else { // beta is a function
509  // get function number
510  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
511  // evaluate function
512  b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
513  }
514  f = fTheory->Func(time, par, fFuncValues);
515  asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
516  break;
517  default:
518  asymFcnValue = 0.0;
519  break;
520  }
521  fData.AppendTheoryValue(asymFcnValue);
522  }
523 
524  // clean up
525  par.clear();
526 }
527 
528 //--------------------------------------------------------------------------
529 // PrepareData (protected)
530 //--------------------------------------------------------------------------
551 {
552  if (!fValid)
553  return false;
554 
555  // keep the Global block info
556  PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
557 
558  // get forward/backward histo from PRunDataHandler object ------------------------
559  // get the correct run
561  if (!runData) { // run not found
562  std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
563  std::cerr << std::endl;
564  return false;
565  }
566 
567  // keep the field from the meta-data from the data-file
568  fMetaData.fField = runData->GetField();
569 
570  // keep the energy from the meta-data from the data-file
571  fMetaData.fEnergy = runData->GetEnergy();
572 
573  // keep the temperature(s) from the meta-data from the data-file
574  for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
575  fMetaData.fTemp.push_back(runData->GetTemperature(i));
576 
577  // collect histogram numbers
578  PUIntVector forwardHistoNo;
579  PUIntVector backwardHistoNo;
580  for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
581  forwardHistoNo.push_back(fRunInfo->GetForwardHistoNo(i));
582 
583  if (!runData->IsPresent(forwardHistoNo[i])) {
584  std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **PANIC ERROR**:";
585  std::cerr << std::endl << ">> forwardHistoNo found = " << forwardHistoNo[i] << ", which is NOT present in the data file!?!?";
586  std::cerr << std::endl << ">> Will quit :-(";
587  std::cerr << std::endl;
588  // clean up
589  forwardHistoNo.clear();
590  backwardHistoNo.clear();
591  return false;
592  }
593  }
594  for (UInt_t i=0; i<fRunInfo->GetBackwardHistoNoSize(); i++) {
595  backwardHistoNo.push_back(fRunInfo->GetBackwardHistoNo(i));
596 
597  if (!runData->IsPresent(backwardHistoNo[i])) {
598  std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **PANIC ERROR**:";
599  std::cerr << std::endl << ">> backwardHistoNo found = " << backwardHistoNo[i] << ", which is NOT present in the data file!?!?";
600  std::cerr << std::endl << ">> Will quit :-(";
601  std::cerr << std::endl;
602  // clean up
603  forwardHistoNo.clear();
604  backwardHistoNo.clear();
605  return false;
606  }
607  }
608 
609  // keep the time resolution in (us)
610  fTimeResolution = runData->GetTimeResolution()/1.0e3;
611  std::cout.precision(10);
612  std::cout << std::endl << ">> PRunAsymmetry::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
613 
614  // get all the proper t0's and addt0's for the current RUN block
615  if (!GetProperT0(runData, globalBlock, forwardHistoNo, backwardHistoNo)) {
616  return false;
617  }
618 
619  // keep the histo of each group at this point (addruns handled below)
620  std::vector<PDoubleVector> forward, backward;
621  forward.resize(forwardHistoNo.size()); // resize to number of groups
622  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
623  forward[i].resize(runData->GetDataBin(forwardHistoNo[i])->size());
624  forward[i] = *runData->GetDataBin(forwardHistoNo[i]);
625  }
626  backward.resize(backwardHistoNo.size()); // resize to number of groups
627  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
628  backward[i].resize(runData->GetDataBin(backwardHistoNo[i])->size());
629  backward[i] = *runData->GetDataBin(backwardHistoNo[i]);
630  }
631 
632  // check if addrun's are present, and if yes add data
633  // check if there are runs to be added to the current one
634  if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
635  PRawRunData *addRunData;
636  for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
637  // get run to be added to the main one
638  addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
639  if (addRunData == nullptr) { // couldn't get run
640  std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
641  std::cerr << std::endl;
642  return false;
643  }
644 
645  // add forward run
646  UInt_t addRunSize;
647  for (UInt_t k=0; k<forwardHistoNo.size(); k++) { // fill each group
648  addRunSize = addRunData->GetDataBin(forwardHistoNo[k])->size();
649  for (UInt_t j=0; j<addRunData->GetDataBin(forwardHistoNo[k])->size(); j++) { // loop over the bin indices
650  // make sure that the index stays in the proper range
651  if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]) >= 0) &&
652  (j+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]) < addRunSize)) {
653  forward[k][j] += addRunData->GetDataBin(forwardHistoNo[k])->at(j+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]));
654  }
655  }
656  }
657 
658  // add backward run
659  for (UInt_t k=0; k<backwardHistoNo.size(); k++) { // fill each group
660  addRunSize = addRunData->GetDataBin(backwardHistoNo[k])->size();
661  for (UInt_t j=0; j<addRunData->GetDataBin(backwardHistoNo[k])->size(); j++) { // loop over the bin indices
662  // make sure that the index stays in the proper range
663  if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][2*k+1])-static_cast<Int_t>(fT0s[2*k+1]) >= 0) &&
664  (j+static_cast<Int_t>(fAddT0s[i-1][2*k+1])-static_cast<Int_t>(fT0s[2*k+1]) < addRunSize)) {
665  backward[k][j] += addRunData->GetDataBin(backwardHistoNo[k])->at(j+static_cast<Int_t>(fAddT0s[i-1][2*k+1])-static_cast<Int_t>(fT0s[2*k+1]));
666  }
667  }
668  }
669  }
670  }
671 
672  // set forward/backward histo data of the first group
673  fForward.resize(forward[0].size());
674  for (UInt_t i=0; i<fForward.size(); i++) {
675  fForward[i] = forward[0][i];
676  }
677  fBackward.resize(backward[0].size());
678  for (UInt_t i=0; i<fBackward.size(); i++) {
679  fBackward[i] = backward[0][i];
680  }
681 
682  // group histograms, add all the remaining forward histograms of the group
683  for (UInt_t i=1; i<forwardHistoNo.size(); i++) { // loop over the groupings
684  for (UInt_t j=0; j<runData->GetDataBin(forwardHistoNo[i])->size(); j++) { // loop over the bin indices
685  // make sure that the index stays within proper range
686  if ((j+fT0s[2*i]-fT0s[0] >= 0) && (j+fT0s[2*i]-fT0s[0] < runData->GetDataBin(forwardHistoNo[i])->size())) {
687  fForward[j] += forward[i][j+static_cast<Int_t>(fT0s[2*i])-static_cast<Int_t>(fT0s[0])];
688  }
689  }
690  }
691 
692  // group histograms, add all the remaining backward histograms of the group
693  for (UInt_t i=1; i<backwardHistoNo.size(); i++) { // loop over the groupings
694  for (UInt_t j=0; j<runData->GetDataBin(backwardHistoNo[i])->size(); j++) { // loop over the bin indices
695  // make sure that the index stays within proper range
696  if ((j+fT0s[2*i+1]-fT0s[1] >= 0) && (j+fT0s[2*i+1]-fT0s[1] < runData->GetDataBin(backwardHistoNo[i])->size())) {
697  fBackward[j] += backward[i][j+static_cast<Int_t>(fT0s[2*i+1])-static_cast<Int_t>(fT0s[1])];
698  }
699  }
700  }
701 
702  // subtract background from histogramms ------------------------------------------
703  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
704  if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
705  if (!SubtractEstimatedBkg())
706  return false;
707  } else { // no background given to do the job, try to estimate it
708  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
709  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
710  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.1), 2);
711  fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.6), 3);
712  std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **WARNING** Neither fix background nor background bins are given!";
713  std::cerr << std::endl << ">> Will try the following:";
714  std::cerr << std::endl << ">> forward: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
715  std::cerr << std::endl << ">> backward: bkg start = " << fRunInfo->GetBkgRange(2) << ", bkg end = " << fRunInfo->GetBkgRange(3);
716  std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
717  std::cerr << std::endl;
718  if (!SubtractEstimatedBkg())
719  return false;
720  }
721  } else { // fixed background given
722  if (!SubtractFixBkg())
723  return false;
724  }
725 
726  UInt_t histoNo[2] = {forwardHistoNo[0], backwardHistoNo[0]};
727 
728  // get the data range (fgb/lgb) for the current RUN block
729  if (!GetProperDataRange(runData, histoNo)) {
730  return false;
731  }
732 
733  // get the fit range for the current RUN block
734  GetProperFitRange(globalBlock);
735 
736  // everything looks fine, hence fill data set
737  Bool_t status;
738  switch(fHandleTag) {
739  case kFit:
741  break;
742  case kView:
743  if (fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking == 0)
744  status = PrepareViewData(runData, histoNo);
745  else
746  status = PrepareRRFViewData(runData, histoNo);
747  break;
748  default:
749  status = false;
750  break;
751  }
752 
753  // clean up
754  forwardHistoNo.clear();
755  backwardHistoNo.clear();
756 
757  return status;
758 }
759 
760 //--------------------------------------------------------------------------
761 // SubtractFixBkg (private)
762 //--------------------------------------------------------------------------
780 {
781  if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) {
782  std::cerr << "PRunAsymmetry::SubtractFixBkg(): **ERROR** no fixed bkg for forward set. Will do nothing here." << std::endl;
783  return false;
784  }
785  if (fRunInfo->GetBkgFix(1) == PMUSR_UNDEFINED) {
786  std::cerr << "PRunAsymmetry::SubtractFixBkg(): **ERROR** no fixed bkg for backward set. Will do nothing here." << std::endl;
787  std::cerr << " you need an entry like:" << std::endl;
788  std::cerr << " backgr.fix 2 3" << std::endl;
789  std::cerr << " i.e. two entries for forward and backward." << std::endl;
790  return false;
791  }
792 
793 
794  Double_t dval;
795 
796  for (UInt_t i=0; i<fForward.size(); i++) {
797  // keep the error, and make sure that the bin is NOT empty
798  if (fForward[i] != 0.0)
799  dval = TMath::Sqrt(fForward[i]);
800  else
801  dval = 1.0;
802  fForwardErr.push_back(dval);
803  fForward[i] -= fRunInfo->GetBkgFix(0);
804 
805  // keep the error, and make sure that the bin is NOT empty
806  if (fBackward[i] != 0.0)
807  dval = TMath::Sqrt(fBackward[i]);
808  else
809  dval = 1.0;
810  fBackwardErr.push_back(dval);
811  fBackward[i] -= fRunInfo->GetBkgFix(1);
812  }
813 
814  return true;
815 }
816 
817 //--------------------------------------------------------------------------
818 // SubtractEstimatedBkg (private)
819 //--------------------------------------------------------------------------
838 {
839  Double_t beamPeriod = 0.0;
840 
841  // check if data are from PSI, RAL, or TRIUMF
842  if (fRunInfo->GetInstitute()->Contains("psi"))
843  beamPeriod = ACCEL_PERIOD_PSI;
844  else if (fRunInfo->GetInstitute()->Contains("ral"))
845  beamPeriod = ACCEL_PERIOD_RAL;
846  else if (fRunInfo->GetInstitute()->Contains("triumf"))
847  beamPeriod = ACCEL_PERIOD_TRIUMF;
848  else
849  beamPeriod = 0.0;
850 
851  // check if start and end are in proper order
852  Int_t start[2] = {fRunInfo->GetBkgRange(0), fRunInfo->GetBkgRange(2)};
853  Int_t end[2] = {fRunInfo->GetBkgRange(1), fRunInfo->GetBkgRange(3)};
854  for (UInt_t i=0; i<2; i++) {
855  if (end[i] < start[i]) {
856  std::cout << std::endl << "PRunAsymmetry::SubtractEstimatedBkg(): end = " << end[i] << " > start = " << start[i] << "! Will swap them!";
857  UInt_t keep = end[i];
858  end[i] = start[i];
859  start[i] = keep;
860  }
861  }
862 
863  // calculate proper background range
864  for (UInt_t i=0; i<2; i++) {
865  if (beamPeriod != 0.0) {
866  Double_t timeBkg = static_cast<Double_t>(end[i]-start[i])*(fTimeResolution*fPacking); // length of the background intervall in time
867  UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
868  // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
869  end[i] = start[i] + static_cast<UInt_t>((fullCycles*beamPeriod)/(fTimeResolution*fPacking));
870  std::cout << "PRunAsymmetry::SubtractEstimatedBkg(): Background " << start[i] << ", " << end[i] << std::endl;
871  if (end[i] == start[i])
872  end[i] = fRunInfo->GetBkgRange(2*i+1);
873  }
874  }
875 
876  // check if start is within histogram bounds
877  if ((start[0] < 0) || (start[0] >= fForward.size()) ||
878  (start[1] < 0) || (start[1] >= fBackward.size())) {
879  std::cerr << std::endl << ">> PRunAsymmetry::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
880  std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
881  std::cerr << std::endl << ">> background start (f/b) = (" << start[0] << "/" << start[1] << ").";
882  return false;
883  }
884 
885  // check if end is within histogram bounds
886  if ((end[0] < 0) || (end[0] >= fForward.size()) ||
887  (end[1] < 0) || (end[1] >= fBackward.size())) {
888  std::cerr << std::endl << ">> PRunAsymmetry::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
889  std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
890  std::cerr << std::endl << ">> background end (f/b) = (" << end[0] << "/" << end[1] << ").";
891  return false;
892  }
893 
894  // calculate background
895  Double_t bkg[2] = {0.0, 0.0};
896  Double_t errBkg[2] = {0.0, 0.0};
897 
898  // forward
899  for (UInt_t i=start[0]; i<=end[0]; i++)
900  bkg[0] += fForward[i];
901  errBkg[0] = TMath::Sqrt(bkg[0])/(end[0] - start[0] + 1);
902  bkg[0] /= static_cast<Double_t>(end[0] - start[0] + 1);
903  std::cout << std::endl << ">> estimated forward histo background: " << bkg[0];
904 
905  // backward
906  for (UInt_t i=start[1]; i<=end[1]; i++)
907  bkg[1] += fBackward[i];
908  errBkg[1] = TMath::Sqrt(bkg[1])/(end[1] - start[1] + 1);
909  bkg[1] /= static_cast<Double_t>(end[1] - start[1] + 1);
910  std::cout << std::endl << ">> estimated backward histo background: " << bkg[1] << std::endl;
911 
912  // correct error for forward, backward
913  Double_t errVal = 0.0;
914  for (UInt_t i=0; i<fForward.size(); i++) {
915  if (fForward[i] > 0.0)
916  errVal = TMath::Sqrt(fForward[i]+errBkg[0]*errBkg[0]);
917  else
918  errVal = 1.0;
919  fForwardErr.push_back(errVal);
920  if (fBackward[i] > 0.0)
921  errVal = TMath::Sqrt(fBackward[i]+errBkg[1]*errBkg[1]);
922  else
923  errVal = 1.0;
924  fBackwardErr.push_back(errVal);
925  }
926 
927  // subtract background from data
928  for (UInt_t i=0; i<fForward.size(); i++) {
929  fForward[i] -= bkg[0];
930  fBackward[i] -= bkg[1];
931  }
932 
933  fRunInfo->SetBkgEstimated(bkg[0], 0);
934  fRunInfo->SetBkgEstimated(bkg[1], 1);
935 
936  return true;
937 }
938 
939 //--------------------------------------------------------------------------
940 // PrepareFitData (protected)
941 //--------------------------------------------------------------------------
954 {
955  // transform raw histo data. This is done the following way (for details see the manual):
956  // first rebin the data, than calculate the asymmetry
957 
958  // everything looks fine, hence fill packed forward and backward histo
959  PRunData forwardPacked;
960  PRunData backwardPacked;
961  Double_t value = 0.0;
962  Double_t error = 0.0;
963  // forward
964  for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
965  if (fPacking == 1) {
966  forwardPacked.AppendValue(fForward[i]);
967  forwardPacked.AppendErrorValue(fForwardErr[i]);
968  } else { // packed data, i.e. fPacking > 1
969  if (((i-fGoodBins[0]) % fPacking == 0) && (i != fGoodBins[0])) { // fill data
970  // in order that after rebinning the fit does not need to be redone (important for plots)
971  // the value is normalize to per bin
972  value /= fPacking;
973  forwardPacked.AppendValue(value);
974  if (value == 0.0)
975  forwardPacked.AppendErrorValue(1.0);
976  else
977  forwardPacked.AppendErrorValue(TMath::Sqrt(error)/fPacking);
978  value = 0.0;
979  error = 0.0;
980  }
981  value += fForward[i];
982  error += fForwardErr[i]*fForwardErr[i];
983  }
984  }
985  // backward
986  for (Int_t i=fGoodBins[2]; i<fGoodBins[3]; i++) {
987  if (fPacking == 1) {
988  backwardPacked.AppendValue(fBackward[i]);
989  backwardPacked.AppendErrorValue(fBackwardErr[i]);
990  } else { // packed data, i.e. fPacking > 1
991  if (((i-fGoodBins[2]) % fPacking == 0) && (i != fGoodBins[2])) { // fill data
992  // in order that after rebinning the fit does not need to be redone (important for plots)
993  // the value is normalize to per bin
994  value /= fPacking;
995  backwardPacked.AppendValue(value);
996  if (value == 0.0)
997  backwardPacked.AppendErrorValue(1.0);
998  else
999  backwardPacked.AppendErrorValue(TMath::Sqrt(error)/fPacking);
1000  value = 0.0;
1001  error = 0.0;
1002  }
1003  value += fBackward[i];
1004  error += fBackwardErr[i]*fBackwardErr[i];
1005  }
1006  }
1007 
1008  // check if packed forward and backward hist have the same size, otherwise take the minimum size
1009  UInt_t noOfBins = forwardPacked.GetValue()->size();
1010  if (forwardPacked.GetValue()->size() != backwardPacked.GetValue()->size()) {
1011  if (forwardPacked.GetValue()->size() > backwardPacked.GetValue()->size())
1012  noOfBins = backwardPacked.GetValue()->size();
1013  }
1014 
1015  // form asymmetry including error propagation
1016  Double_t asym;
1017  Double_t f, b, ef, eb;
1018  // fill data time start, and step
1019  // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1020  fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(fGoodBins[0])-0.5) + static_cast<Double_t>(fPacking)/2.0 - static_cast<Double_t>(fT0s[0])));
1021  fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(fPacking));
1022  for (UInt_t i=0; i<noOfBins; i++) {
1023  // to make the formulae more readable
1024  f = forwardPacked.GetValue()->at(i);
1025  b = backwardPacked.GetValue()->at(i);
1026  ef = forwardPacked.GetError()->at(i);
1027  eb = backwardPacked.GetError()->at(i);
1028  // check that there are indeed bins
1029  if (f+b != 0.0)
1030  asym = (f-b) / (f+b);
1031  else
1032  asym = 0.0;
1033  fData.AppendValue(asym);
1034  // calculate the error
1035  if (f+b != 0.0)
1036  error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1037  else
1038  error = 1.0;
1039  fData.AppendErrorValue(error);
1040  }
1041 
1042  CalcNoOfFitBins();
1043 
1044  // clean up
1045  fForward.clear();
1046  fForwardErr.clear();
1047  fBackward.clear();
1048  fBackwardErr.clear();
1049 
1050  return true;
1051 }
1052 
1053 //--------------------------------------------------------------------------
1054 // PrepareViewData (protected)
1055 //--------------------------------------------------------------------------
1072 Bool_t PRunAsymmetry::PrepareViewData(PRawRunData* runData, UInt_t histoNo[2])
1073 {
1074  // check if view_packing is wished
1075  Int_t packing = fPacking;
1076  if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
1077  packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1078  }
1079 
1080  // feed the parameter vector
1081  std::vector<Double_t> par;
1082  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1083  for (UInt_t i=0; i<paramList->size(); i++)
1084  par.push_back((*paramList)[i].fValue);
1085 
1086  // transform raw histo data. This is done the following way (for details see the manual):
1087  // first rebin the data, than calculate the asymmetry
1088 
1089  // first get start data, end data, and t0
1090  Int_t start[2] = {fGoodBins[0], fGoodBins[2]};
1091  Int_t end[2] = {fGoodBins[1], fGoodBins[3]};
1092  Int_t t0[2] = {static_cast<Int_t>(fT0s[0]), static_cast<Int_t>(fT0s[1])};
1093 
1094  // check if the data ranges and t0's between forward/backward are compatible
1095  Int_t fgb[2];
1096  if (start[0]-t0[0] != start[1]-t0[1]) { // wrong fgb aligning
1097  if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
1098  fgb[0] = start[0];
1099  fgb[1] = t0[1] + start[0]-t0[0];
1100  std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **WARNING** needed to shift backward fgb from ";
1101  std::cerr << start[1] << " to " << fgb[1] << std::endl;
1102  } else {
1103  fgb[0] = t0[0] + start[1]-t0[1];
1104  fgb[1] = start[1];
1105  std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **WARNING** needed to shift forward fgb from ";
1106  std::cerr << start[0] << " to " << fgb[0] << std::endl;
1107  }
1108  } else { // fgb aligning is correct
1109  fgb[0] = start[0];
1110  fgb[1] = start[1];
1111  }
1112 
1113  Int_t val = fgb[0]-packing*(fgb[0]/packing);
1114  do {
1115  if (fgb[1] - fgb[0] < 0)
1116  val += packing;
1117  } while (val + fgb[1] - fgb[0] < 0);
1118 
1119  start[0] = val;
1120  start[1] = val + fgb[1] - fgb[0];
1121 
1122  // make sure that there are equal number of rebinned bins in forward and backward
1123  UInt_t noOfBins0 = (runData->GetDataBin(histoNo[0])->size()-start[0])/packing;
1124  UInt_t noOfBins1 = (runData->GetDataBin(histoNo[1])->size()-start[1])/packing;
1125  if (noOfBins0 > noOfBins1)
1126  noOfBins0 = noOfBins1;
1127  end[0] = start[0] + noOfBins0 * packing;
1128  end[1] = start[1] + noOfBins0 * packing;
1129 
1130  // check if start, end, and t0 make any sense
1131  // 1st check if start and end are in proper order
1132  for (UInt_t i=0; i<2; i++) {
1133  if (end[i] < start[i]) { // need to swap them
1134  Int_t keep = end[i];
1135  end[i] = start[i];
1136  start[i] = keep;
1137  }
1138  // 2nd check if start is within proper bounds
1139  if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1140  std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
1141  std::cerr << std::endl;
1142  return false;
1143  }
1144  // 3rd check if end is within proper bounds
1145  if ((end[i] < 0) || (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1146  std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
1147  std::cerr << std::endl;
1148  return false;
1149  }
1150  // 4th check if t0 is within proper bounds
1151  if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1152  std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **ERROR** t0 data bin doesn't make any sense!";
1153  std::cerr << std::endl;
1154  return false;
1155  }
1156  }
1157 
1158  // everything looks fine, hence fill packed forward and backward histo
1159  PRunData forwardPacked;
1160  PRunData backwardPacked;
1161  Double_t value = 0.0;
1162  Double_t error = 0.0;
1163 
1164  // forward
1165  for (Int_t i=start[0]; i<end[0]; i++) {
1166  if (packing == 1) {
1167  forwardPacked.AppendValue(fForward[i]);
1168  forwardPacked.AppendErrorValue(fForwardErr[i]);
1169  } else { // packed data, i.e. packing > 1
1170  if (((i-start[0]) % packing == 0) && (i != start[0])) { // fill data
1171  // in order that after rebinning the fit does not need to be redone (important for plots)
1172  // the value is normalize to per bin
1173  value /= packing;
1174  forwardPacked.AppendValue(value);
1175  if (value == 0.0)
1176  forwardPacked.AppendErrorValue(1.0);
1177  else
1178  forwardPacked.AppendErrorValue(TMath::Sqrt(error)/packing);
1179  value = 0.0;
1180  error = 0.0;
1181  }
1182  value += fForward[i];
1183  error += fForwardErr[i]*fForwardErr[i];
1184  }
1185  }
1186 
1187  // backward
1188  for (Int_t i=start[1]; i<end[1]; i++) {
1189  if (packing == 1) {
1190  backwardPacked.AppendValue(fBackward[i]);
1191  backwardPacked.AppendErrorValue(fBackwardErr[i]);
1192  } else { // packed data, i.e. packing > 1
1193  if (((i-start[1]) % packing == 0) && (i != start[1])) { // fill data
1194  // in order that after rebinning the fit does not need to be redone (important for plots)
1195  // the value is normalize to per bin
1196  value /= packing;
1197  backwardPacked.AppendValue(value);
1198  if (value == 0.0)
1199  backwardPacked.AppendErrorValue(1.0);
1200  else
1201  backwardPacked.AppendErrorValue(TMath::Sqrt(error)/packing);
1202  value = 0.0;
1203  error = 0.0;
1204  }
1205  value += fBackward[i];
1206  error += fBackwardErr[i]*fBackwardErr[i];
1207  }
1208  }
1209 
1210  // check if packed forward and backward hist have the same size, otherwise take the minimum size
1211  UInt_t noOfBins = forwardPacked.GetValue()->size();
1212  if (forwardPacked.GetValue()->size() != backwardPacked.GetValue()->size()) {
1213  if (forwardPacked.GetValue()->size() > backwardPacked.GetValue()->size())
1214  noOfBins = backwardPacked.GetValue()->size();
1215  }
1216 
1217  // form asymmetry including error propagation
1218  Double_t asym;
1219  Double_t f, b, ef, eb, alpha = 1.0, beta = 1.0;
1220  // set data time start, and step
1221  // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1222  fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start[0])-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0[0])));
1223  fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(packing));
1224 
1225  // get the proper alpha and beta
1226  switch (fAlphaBetaTag) {
1227  case 1: // alpha == 1, beta == 1
1228  alpha = 1.0;
1229  beta = 1.0;
1230  break;
1231  case 2: // alpha != 1, beta == 1
1232  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1233  alpha = par[fRunInfo->GetAlphaParamNo()-1];
1234  } else { // alpha is function
1235  // get function number
1236  UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1237  // evaluate function
1238  alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1239  }
1240  beta = 1.0;
1241  break;
1242  case 3: // alpha == 1, beta != 1
1243  alpha = 1.0;
1244  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1245  beta = par[fRunInfo->GetBetaParamNo()-1];
1246  } else { // beta is a function
1247  // get function number
1248  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1249  // evaluate function
1250  beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1251  }
1252  break;
1253  case 4: // alpha != 1, beta != 1
1254  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1255  alpha = par[fRunInfo->GetAlphaParamNo()-1];
1256  } else { // alpha is function
1257  // get function number
1258  UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1259  // evaluate function
1260  alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1261  }
1262  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1263  beta = par[fRunInfo->GetBetaParamNo()-1];
1264  } else { // beta is a function
1265  // get function number
1266  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1267  // evaluate function
1268  beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1269  }
1270  break;
1271  default:
1272  break;
1273  }
1274 
1275  for (UInt_t i=0; i<forwardPacked.GetValue()->size(); i++) {
1276  // to make the formulae more readable
1277  f = forwardPacked.GetValue()->at(i);
1278  b = backwardPacked.GetValue()->at(i);
1279  ef = forwardPacked.GetError()->at(i);
1280  eb = backwardPacked.GetError()->at(i);
1281  // check that there are indeed bins
1282  if (f+b != 0.0)
1283  asym = (alpha*f-b) / (alpha*beta*f+b);
1284  else
1285  asym = 0.0;
1286  fData.AppendValue(asym);
1287  // calculate the error
1288  if (f+b != 0.0)
1289  error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1290  else
1291  error = 1.0;
1292  fData.AppendErrorValue(error);
1293  }
1294 
1295  CalcNoOfFitBins();
1296 
1297  // clean up
1298  fForward.clear();
1299  fForwardErr.clear();
1300  fBackward.clear();
1301  fBackwardErr.clear();
1302 
1303  // fill theory vector for kView
1304  // calculate functions
1305  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1307  }
1308 
1309  // calculate theory
1310  Double_t time;
1311  UInt_t size = runData->GetDataBin(histoNo[0])->size()/packing;
1312  Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1313 
1315  if (fTheoAsData) { // calculate theory only at the data points
1317  } else {
1318  // finer binning for the theory (8 times as many points = factor)
1319  size *= factor;
1320  fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1321  }
1322 
1323  for (UInt_t i=0; i<size; i++) {
1324  time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1325  value = fTheory->Func(time, par, fFuncValues);
1326  if (fabs(value) > 10.0) { // dirty hack needs to be fixed!!
1327  value = 0.0;
1328  }
1329  fData.AppendTheoryValue(value);
1330  }
1331 
1332  // clean up
1333  par.clear();
1334 
1335  return true;
1336 }
1337 
1338 //--------------------------------------------------------------------------
1339 // PrepareRRFViewData (protected)
1340 //--------------------------------------------------------------------------
1356 Bool_t PRunAsymmetry::PrepareRRFViewData(PRawRunData* runData, UInt_t histoNo[2])
1357 {
1358  // feed the parameter vector
1359  std::vector<Double_t> par;
1360  PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1361  for (UInt_t i=0; i<paramList->size(); i++)
1362  par.push_back((*paramList)[i].fValue);
1363 
1364  // ------------------------------------------------------------
1365  // 1. make all necessary checks
1366  // ------------------------------------------------------------
1367 
1368  // first get start data, end data, and t0
1369  Int_t start[2] = {fGoodBins[0], fGoodBins[2]};
1370  Int_t end[2] = {fGoodBins[1], fGoodBins[3]};
1371  Int_t t0[2] = {static_cast<Int_t>(fT0s[0]), static_cast<Int_t>(fT0s[1])};
1372  UInt_t packing = fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking;
1373 
1374  // check if the data ranges and t0's between forward/backward are compatible
1375  Int_t fgb[2];
1376  if (start[0]-t0[0] != start[1]-t0[1]) { // wrong fgb aligning
1377  if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
1378  fgb[0] = start[0];
1379  fgb[1] = t0[1] + start[0]-t0[0];
1380  std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **WARNING** needed to shift backward fgb from ";
1381  std::cerr << start[1] << " to " << fgb[1] << std::endl;
1382  } else {
1383  fgb[0] = t0[0] + start[1]-t0[1];
1384  fgb[1] = start[1];
1385  std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **WARNING** needed to shift forward fgb from ";
1386  std::cerr << start[1] << " to " << fgb[0] << std::endl;
1387  }
1388  } else { // fgb aligning is correct
1389  fgb[0] = start[0];
1390  fgb[1] = start[1];
1391  }
1392 
1393  Int_t val = fgb[0]-packing*(fgb[0]/packing);
1394  do {
1395  if (fgb[1] - fgb[0] < 0)
1396  val += packing;
1397  } while (val + fgb[1] - fgb[0] < 0);
1398 
1399  start[0] = val;
1400  start[1] = val + fgb[1] - fgb[0];
1401 
1402  // make sure that there are equal number of rebinned bins in forward and backward
1403  UInt_t noOfBins0 = runData->GetDataBin(histoNo[0])->size()-start[0];
1404  UInt_t noOfBins1 = runData->GetDataBin(histoNo[1])->size()-start[1];
1405  if (noOfBins0 > noOfBins1)
1406  noOfBins0 = noOfBins1;
1407  end[0] = start[0] + noOfBins0;
1408  end[1] = start[1] + noOfBins0;
1409 
1410  // check if start, end, and t0 make any sense
1411  // 1st check if start and end are in proper order
1412  for (UInt_t i=0; i<2; i++) {
1413  if (end[i] < start[i]) { // need to swap them
1414  Int_t keep = end[i];
1415  end[i] = start[i];
1416  start[i] = keep;
1417  }
1418  // 2nd check if start is within proper bounds
1419  if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1420  std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **ERROR** start data bin doesn't make any sense!";
1421  std::cerr << std::endl;
1422  return false;
1423  }
1424  // 3rd check if end is within proper bounds
1425  if ((end[i] < 0) || (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1426  std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **ERROR** end data bin doesn't make any sense!";
1427  std::cerr << std::endl;
1428  return false;
1429  }
1430  // 4th check if t0 is within proper bounds
1431  if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1432  std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **ERROR** t0 data bin doesn't make any sense!";
1433  std::cerr << std::endl;
1434  return false;
1435  }
1436  }
1437 
1438  // ------------------------------------------------------------
1439  // 2. build the asymmetry [A(t)] WITHOUT packing.
1440  // ------------------------------------------------------------
1441 
1442  PDoubleVector forward, forwardErr;
1443  PDoubleVector backward, backwardErr;
1444  Double_t error = 0.0;
1445  // forward
1446  for (Int_t i=start[0]; i<end[0]; i++) {
1447  forward.push_back(fForward[i]);
1448  forwardErr.push_back(fForwardErr[i]);
1449  }
1450  // backward
1451  for (Int_t i=start[1]; i<end[1]; i++) {
1452  backward.push_back(fBackward[i]);
1453  backwardErr.push_back(fBackwardErr[i]);
1454  }
1455 
1456  // check if packed forward and backward hist have the same size, otherwise take the minimum size
1457  UInt_t noOfBins = forward.size();
1458  if (forward.size() != backward.size()) {
1459  if (forward.size() > backward.size())
1460  noOfBins = backward.size();
1461  }
1462 
1463  // form asymmetry including error propagation
1464  PDoubleVector asymmetry, asymmetryErr;
1465  Double_t asym;
1466  Double_t f, b, ef, eb, alpha = 1.0, beta = 1.0;
1467 
1468  // get the proper alpha and beta
1469  switch (fAlphaBetaTag) {
1470  case 1: // alpha == 1, beta == 1
1471  alpha = 1.0;
1472  beta = 1.0;
1473  break;
1474  case 2: // alpha != 1, beta == 1
1475  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1476  alpha = par[fRunInfo->GetAlphaParamNo()-1];
1477  } else { // alpha is function
1478  // get function number
1479  UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1480  // evaluate function
1481  alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1482  }
1483  beta = 1.0;
1484  break;
1485  case 3: // alpha == 1, beta != 1
1486  alpha = 1.0;
1487  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1488  beta = par[fRunInfo->GetBetaParamNo()-1];
1489  } else { // beta is a function
1490  // get function number
1491  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1492  // evaluate function
1493  beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1494  }
1495  break;
1496  case 4: // alpha != 1, beta != 1
1497  if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1498  alpha = par[fRunInfo->GetAlphaParamNo()-1];
1499  } else { // alpha is function
1500  // get function number
1501  UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1502  // evaluate function
1503  alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1504  }
1505  if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1506  beta = par[fRunInfo->GetBetaParamNo()-1];
1507  } else { // beta is a function
1508  // get function number
1509  UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1510  // evaluate function
1511  beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1512  }
1513  break;
1514  default:
1515  break;
1516  }
1517 
1518  for (UInt_t i=0; i<noOfBins; i++) {
1519  // to make the formulae more readable
1520  f = forward[i];
1521  b = backward[i];
1522  ef = forwardErr[i];
1523  eb = backwardErr[i];
1524  // check that there are indeed bins
1525  if (f+b != 0.0)
1526  asym = (alpha*f-b) / (alpha*beta*f+b);
1527  else
1528  asym = 0.0;
1529  asymmetry.push_back(asym);
1530  // calculate the error
1531  if (f+b != 0.0)
1532  error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1533  else
1534  error = 1.0;
1535  asymmetryErr.push_back(error);
1536  }
1537 
1538  // clean up
1539  fForward.clear();
1540  fForwardErr.clear();
1541  fBackward.clear();
1542  fBackwardErr.clear();
1543 
1544 
1545  // ------------------------------------------------------------
1546  // 3. A_R(t) = A(t) * 2 cos(w_R t + phi_R)
1547  // ------------------------------------------------------------
1548 
1549  // check which units shall be used
1550  Double_t gammaRRF = 0.0, wRRF = 0.0, phaseRRF = 0.0;
1551  Double_t time;
1552 
1553  switch (fMsrInfo->GetMsrPlotList()->at(0).fRRFUnit) {
1554  case RRF_UNIT_kHz:
1555  gammaRRF = TMath::TwoPi()*1.0e-3;
1556  break;
1557  case RRF_UNIT_MHz:
1558  gammaRRF = TMath::TwoPi();
1559  break;
1560  case RRF_UNIT_Mcs:
1561  gammaRRF = 1.0;
1562  break;
1563  case RRF_UNIT_G:
1564  gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi();
1565  break;
1566  case RRF_UNIT_T:
1567  gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi()*1.0e4;
1568  break;
1569  default:
1570  gammaRRF = TMath::TwoPi();
1571  break;
1572  }
1573  wRRF = gammaRRF * fMsrInfo->GetMsrPlotList()->at(0).fRRFFreq;
1574  phaseRRF = fMsrInfo->GetMsrPlotList()->at(0).fRRFPhase / 180.0 * TMath::Pi();
1575 
1576  for (UInt_t i=0; i<asymmetry.size(); i++) {
1577  time = fTimeResolution*(static_cast<Double_t>(start[0])-t0[0]+static_cast<Double_t>(i));
1578  asymmetry[i] *= 2.0*TMath::Cos(wRRF*time+phaseRRF);
1579  }
1580 
1581  // ------------------------------------------------------------
1582  // 4. do the packing of A_R(t)
1583  // ------------------------------------------------------------
1584  Double_t value = 0.0;
1585  error = 0.0;
1586  for (UInt_t i=0; i<asymmetry.size(); i++) {
1587  if ((i % packing == 0) && (i != 0)) {
1588  value /= packing;
1589  fData.AppendValue(value);
1590  fData.AppendErrorValue(TMath::Sqrt(error)/packing);
1591 
1592  value = 0.0;
1593  error = 0.0;
1594  }
1595  value += asymmetry[i];
1596  error += asymmetryErr[i]*asymmetryErr[i];
1597  }
1598 
1599  // set data time start, and step
1600  // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1601  fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start[0])-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0[0])));
1602  fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(packing));
1603 
1604  // ------------------------------------------------------------
1605  // 5. calculate theory [T(t)] as close as possible to the time resolution [compatible with the RRF frequency]
1606  // 6. T_R(t) = T(t) * 2 cos(w_R t + phi_R)
1607  // ------------------------------------------------------------
1608  UInt_t rebinRRF = static_cast<UInt_t>((TMath::Pi()/2.0/wRRF)/fTimeResolution); // RRF time resolution / data time resolution
1610  fData.SetTheoryTimeStep(TMath::Pi()/2.0/wRRF/rebinRRF); // = theory time resolution as close as possible to the data time resolution compatible with wRRF
1611 
1612  // calculate functions
1613  for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1615  }
1616 
1617  Double_t theoryValue;
1618  for (UInt_t i=0; i<asymmetry.size(); i++) {
1619  time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1620  theoryValue = fTheory->Func(time, par, fFuncValues);
1621  theoryValue *= 2.0*TMath::Cos(wRRF * time + phaseRRF);
1622 
1623  if (fabs(theoryValue) > 10.0) { // dirty hack needs to be fixed!!
1624  theoryValue = 0.0;
1625  }
1626 
1627  fData.AppendTheoryValue(theoryValue);
1628  }
1629 
1630  // ------------------------------------------------------------
1631  // 7. do the packing of T_R(t)
1632  // ------------------------------------------------------------
1633 
1634  PDoubleVector theo;
1635  Double_t dval = 0.0;
1636  for (UInt_t i=0; i<fData.GetTheory()->size(); i++) {
1637  if ((i % rebinRRF == 0) && (i != 0)) {
1638  theo.push_back(dval/rebinRRF);
1639  dval = 0.0;
1640  }
1641  dval += fData.GetTheory()->at(i);
1642  }
1643  fData.ReplaceTheory(theo);
1644  theo.clear();
1645 
1646  // set the theory time start and step size
1647  fData.SetTheoryTimeStart(fData.GetTheoryTimeStart()+static_cast<Double_t>(rebinRRF-1)*fData.GetTheoryTimeStep()/2.0);
1649 
1650  // ------------------------------------------------------------
1651  // 8. calculate the Kaiser FIR filter coefficients
1652  // ------------------------------------------------------------
1653  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)
1654 
1655  // ------------------------------------------------------------
1656  // 9. filter T_R(t)
1657  // ------------------------------------------------------------
1658  FilterTheo();
1659 
1660  // clean up
1661  par.clear();
1662  forward.clear();
1663  forwardErr.clear();
1664  backward.clear();
1665  backwardErr.clear();
1666  asymmetry.clear();
1667  asymmetryErr.clear();
1668 
1669  return true;
1670 }
1671 
1672 //--------------------------------------------------------------------------
1673 // GetProperT0 (private)
1674 //--------------------------------------------------------------------------
1693 Bool_t PRunAsymmetry::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHistoNo, PUIntVector &backwardHistoNo)
1694 {
1695  // feed all T0's
1696  // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1697  fT0s.clear();
1698  // this strange fT0 size estimate is needed in case #forw histos != #back histos
1699  size_t size = 2*forwardHistoNo.size();
1700  if (backwardHistoNo.size() > forwardHistoNo.size())
1701  size = 2*backwardHistoNo.size();
1702  fT0s.resize(size);
1703  for (UInt_t i=0; i<fT0s.size(); i++) {
1704  fT0s[i] = -1.0;
1705  }
1706 
1707  // fill in the T0's from the msr-file (if present)
1708  for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1709  fT0s[i] = fRunInfo->GetT0Bin(i);
1710  }
1711 
1712  // fill in the missing T0's from the GLOBAL block section (if present)
1713  for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1714  if (fT0s[i] == -1) { // i.e. not given in the RUN block section
1715  fT0s[i] = globalBlock->GetT0Bin(i);
1716  }
1717  }
1718 
1719  // fill in the missing T0's from the data file, if not already present in the msr-file
1720  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1721  if (fT0s[2*i] == -1.0) // i.e. not present in the msr-file, try the data file
1722  if (runData->GetT0Bin(forwardHistoNo[i]) > 0.0) {
1723  fT0s[2*i] = runData->GetT0Bin(forwardHistoNo[i]);
1724  fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1725  }
1726  }
1727  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1728  if (fT0s[2*i+1] == -1.0) // i.e. not present in the msr-file, try the data file
1729  if (runData->GetT0Bin(backwardHistoNo[i]) > 0.0) {
1730  fT0s[2*i+1] = runData->GetT0Bin(backwardHistoNo[i]);
1731  fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1732  }
1733  }
1734 
1735  // 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
1736  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1737  if (fT0s[2*i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1738  fT0s[2*i] = runData->GetT0BinEstimated(forwardHistoNo[i]);
1739  fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1740 
1741  std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1742  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1743  std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(forwardHistoNo[i]);
1744  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1745  std::cerr << std::endl;
1746  }
1747  }
1748  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1749  if (fT0s[2*i+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1750  fT0s[2*i+1] = runData->GetT0BinEstimated(backwardHistoNo[i]);
1751  fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1752 
1753  std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1754  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1755  std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[i]);
1756  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1757  std::cerr << std::endl;
1758  }
1759  }
1760 
1761  // check if t0 is within proper bounds
1762  for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1763  if ((fT0s[2*i] < 0) || (fT0s[2*i] > static_cast<Int_t>(runData->GetDataBin(forwardHistoNo[i])->size()))) {
1764  std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **ERROR** t0 data bin (" << fT0s[2*i] << ") doesn't make any sense!";
1765  std::cerr << std::endl << ">> forwardHistoNo " << forwardHistoNo[i];
1766  std::cerr << std::endl;
1767  return false;
1768  }
1769  }
1770  for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1771  if ((fT0s[2*i+1] < 0) || (fT0s[2*i+1] > static_cast<Int_t>(runData->GetDataBin(backwardHistoNo[i])->size()))) {
1772  std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **ERROR** t0 data bin (" << fT0s[2*i+1] << ") doesn't make any sense!";
1773  std::cerr << std::endl << ">> backwardHistoNo " << backwardHistoNo[i];
1774  std::cerr << std::endl;
1775  return false;
1776  }
1777  }
1778 
1779  // check if addrun's are present, and if yes add the necessary t0's
1780  if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1781  PRawRunData *addRunData;
1782  fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1783  for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1784  // get run to be added to the main one
1785  addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
1786  if (addRunData == 0) { // couldn't get run
1787  std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1788  std::cerr << std::endl;
1789  return false;
1790  }
1791 
1792  // feed all T0's
1793  // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1794  fAddT0s[i-1].clear();
1795  fAddT0s[i-1].resize(2*forwardHistoNo.size());
1796  for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1797  fAddT0s[i-1][j] = -1.0;
1798  }
1799 
1800  // fill in the T0's from the msr-file (if present)
1801  for (Int_t j=0; j<fRunInfo->GetAddT0BinSize(i-1); j++) {
1802  fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1, j);
1803  }
1804 
1805  // fill in the T0's from the data file, if not already present in the msr-file
1806  for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1807  if (fAddT0s[i-1][2*j] == -1.0) // i.e. not present in the msr-file, try the data file
1808  if (addRunData->GetT0Bin(forwardHistoNo[j]) > 0.0) {
1809  fAddT0s[i-1][2*j] = addRunData->GetT0Bin(forwardHistoNo[j]);
1810  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1811  }
1812  }
1813  for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1814  if (fAddT0s[i-1][2*j+1] == -1.0) // i.e. not present in the msr-file, try the data file
1815  if (addRunData->GetT0Bin(backwardHistoNo[j]) > 0.0) {
1816  fAddT0s[i-1][2*j+1] = addRunData->GetT0Bin(backwardHistoNo[j]);
1817  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1818  }
1819  }
1820 
1821  // 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
1822  for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1823  if (fAddT0s[i-1][2*j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1824  fAddT0s[i-1][2*j] = addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1825  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1826 
1827  std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1828  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1829  std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1830  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1831  std::cerr << std::endl;
1832  }
1833  }
1834  for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1835  if (fAddT0s[i-1][2*j+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1836  fAddT0s[i-1][2*j+1] = addRunData->GetT0BinEstimated(backwardHistoNo[j]);
1837  fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1838 
1839  std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1840  std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1841  std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[j]);
1842  std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1843  std::cerr << std::endl;
1844  }
1845  }
1846  }
1847  }
1848 
1849  return true;
1850 }
1851 
1852 //--------------------------------------------------------------------------
1853 // GetProperDataRange (private)
1854 //--------------------------------------------------------------------------
1868 Bool_t PRunAsymmetry::GetProperDataRange(PRawRunData* runData, UInt_t histoNo[2])
1869 {
1870  // first get start/end data
1871  Int_t start[2] = {fRunInfo->GetDataRange(0), fRunInfo->GetDataRange(2)};
1872  Int_t end[2] = {fRunInfo->GetDataRange(1), fRunInfo->GetDataRange(3)};
1873  // check if data range has been provided in the RUN block. If not, try the GLOBAL block
1874  if (start[0] == -1) {
1875  start[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1876  }
1877  if (start[1] == -1) {
1878  start[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(2);
1879  }
1880  if (end[0] == -1) {
1881  end[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1882  }
1883  if (end[1] == -1) {
1884  end[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(3);
1885  }
1886 
1887  Double_t t0[2] = {fT0s[0], fT0s[1]};
1888  Int_t offset = static_cast<Int_t>((10.0e-3/fTimeResolution)); // needed in case first good bin is not given, default = 10ns
1889 
1890  // check if data range has been provided, and if not try to estimate them
1891  if (start[0] < 0) {
1892  start[0] = static_cast<Int_t>(t0[0])+offset;
1893  fRunInfo->SetDataRange(start[0], 0);
1894  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[0] << ".";
1895  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1896  std::cerr << std::endl;
1897  }
1898  if (start[1] < 0) {
1899  start[1] = static_cast<Int_t>(t0[1])+offset;
1900  fRunInfo->SetDataRange(start[1], 2);
1901  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[1] << ".";
1902  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1903  std::cerr << std::endl;
1904  }
1905  if (end[0] < 0) {
1906  end[0] = runData->GetDataBin(histoNo[0])->size();
1907  fRunInfo->SetDataRange(end[0], 1);
1908  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range end = " << end[0] << ".";
1909  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1910  std::cerr << std::endl;
1911  }
1912  if (end[1] < 0) {
1913  end[1] = runData->GetDataBin(histoNo[1])->size();
1914  fRunInfo->SetDataRange(end[1], 3);
1915  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range end = " << end[1] << ".";
1916  std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1917  std::cerr << std::endl;
1918  }
1919 
1920  // check if start, end, and t0 make any sense
1921  // 1st check if start and end are in proper order
1922  for (UInt_t i=0; i<2; i++) {
1923  if (end[i] < start[i]) { // need to swap them
1924  Int_t keep = end[i];
1925  end[i] = start[i];
1926  start[i] = keep;
1927  }
1928  // 2nd check if start is within proper bounds
1929  if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1930  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **ERROR** start data bin doesn't make any sense!";
1931  std::cerr << std::endl;
1932  return false;
1933  }
1934  // 3rd check if end is within proper bounds
1935  if (end[i] < 0) {
1936  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **ERROR** end data bin (" << end[i] << ") doesn't make any sense!";
1937  std::cerr << std::endl;
1938  return false;
1939  }
1940  if (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())) {
1941  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** end data bin (" << end[i] << ") > histo length (" << static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()) << ").";
1942  std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
1943  std::cerr << std::endl;
1944  end[i] = static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())-1;
1945  }
1946  // 4th check if t0 is within proper bounds
1947  if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1948  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **ERROR** t0 data bin doesn't make any sense!";
1949  std::cerr << std::endl;
1950  return false;
1951  }
1952  }
1953 
1954  // check that start-t0 is the same for forward as for backward, otherwise take max(start[i]-t0[i])
1955  if (fabs(static_cast<Double_t>(start[0])-t0[0]) > fabs(static_cast<Double_t>(start[1])-t0[1])){
1956  start[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(start[0]) - t0[0]);
1957  end[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(end[0]) - t0[0]);
1958  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift backward data range.";
1959  std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(2) << ", " << fRunInfo->GetDataRange(3);
1960  std::cerr << std::endl << ">> used : " << start[1] << ", " << end[1];
1961  std::cerr << std::endl;
1962  }
1963  if (fabs(static_cast<Double_t>(start[0])-t0[0]) < fabs(static_cast<Double_t>(start[1])-t0[1])){
1964  start[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(start[1]) - t0[1]);
1965  end[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(end[1]) - t0[1]);
1966  std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift forward data range.";
1967  std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(0) << ", " << fRunInfo->GetDataRange(1);
1968  std::cerr << std::endl << ">> used : " << start[0] << ", " << end[0];
1969  std::cerr << std::endl;
1970  }
1971 
1972  // keep good bins for potential latter use
1973  fGoodBins[0] = start[0];
1974  fGoodBins[1] = end[0];
1975  fGoodBins[2] = start[1];
1976  fGoodBins[3] = end[1];
1977 
1978  return true;
1979 }
1980 
1981 //--------------------------------------------------------------------------
1982 // GetProperFitRange (private)
1983 //--------------------------------------------------------------------------
1997 {
1998  // set fit start/end time; first check RUN Block
2001  // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
2002  if (fRunInfo->IsFitRangeInBin()) {
2003  fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
2004  fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
2005  // write these times back into the data structure. This way it is available when writting the log-file
2008  }
2009  if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
2010  fFitStartTime = globalBlock->GetFitRange(0);
2011  fFitEndTime = globalBlock->GetFitRange(1);
2012  // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
2013  if (globalBlock->IsFitRangeInBin()) {
2014  fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
2015  fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
2016  // write these times back into the data structure. This way it is available when writting the log-file
2017  globalBlock->SetFitRange(fFitStartTime, 0);
2018  globalBlock->SetFitRange(fFitEndTime, 1);
2019  }
2020  }
2022  fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
2023  fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
2024  std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
2025  std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
2026  }
2027 }
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
UInt_t fNoOfFitBins
number of bins to be be fitted
Definition: PRunAsymmetry.h:69
virtual Bool_t GetProperDataRange(PRawRunData *runData, UInt_t histoNo[2])
virtual const PDoubleVector * GetError()
Definition: PMusr.h:249
Int_t fStartTimeBin
bin at which the fit starts
Definition: PRunAsymmetry.h:80
Definition: PMusr.h:220
virtual Bool_t PrepareFitData()
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
virtual const Double_t GetTimeResolution()
Definition: PMusr.h:429
Int_t fEndTimeBin
bin at which the fit ends
Definition: PRunAsymmetry.h:81
virtual Int_t GetPacking()
Definition: PMusr.h:668
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
UInt_t fAlphaBetaTag
; ; ; .
Definition: PRunAsymmetry.h:68
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHisto, PUIntVector &backwardHistoNo)
virtual void SetTheoryTimeStart(Double_t dval)
Definition: PMusr.h:255
PMsrHandler * fMsrInfo
msr-file handler
Definition: PRunBase.h:71
virtual Bool_t IsFitRangeInBin()
Definition: PMusr.h:665
PDoubleVector fForwardErr
forward histo errors
Definition: PRunAsymmetry.h:74
virtual const PDoublePairVector * GetTemperature() const
Definition: PMusr.h:422
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
PDoubleVector fTemp
temperature(s) in (K)
Definition: PMusr.h:229
virtual void AppendErrorValue(Double_t dval)
Definition: PMusr.h:260
virtual Bool_t IsFitRangeInBin()
Definition: PMusr.h:588
virtual ~PRunAsymmetry()
virtual Double_t GetTheoryTimeStep()
Definition: PMusr.h:245
Int_t fGoodBins[4]
keep first/last good bins. 0=fgb, 1=lgb (forward); 2=fgb, 3=lgb (backward)
Definition: PRunAsymmetry.h:78
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 void CalcNoOfFitBins()
virtual const Double_t GetT0Bin(const UInt_t histoNo)
Definition: PMusr.h:431
Bool_t fTheoAsData
true=only calculate the theory points at the data points, false=calculate more points for the theory ...
Definition: PRunAsymmetry.h:71
virtual Int_t GetBackwardHistoNo(UInt_t idx=0)
Definition: PMusr.cpp:1444
Bool_t SubtractFixBkg()
PDoubleVector fBackwardErr
backward histo errors
Definition: PRunAsymmetry.h:76
virtual void SetT0Bin(Double_t dval, Int_t idx=-1)
Definition: PMusr.cpp:1712
#define RRF_UNIT_T
Definition: PMusr.h:158
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
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 UInt_t GetNoOfFitBins()
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 void SetAddT0Bin(Double_t dval, UInt_t addRunIdx, UInt_t histoNoIdx)
Definition: PMusr.cpp:1788
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
PDoubleVector fForward
forward histo data
Definition: PRunAsymmetry.h:73
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
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 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 Bool_t PrepareViewData(PRawRunData *runData, UInt_t histoNo[2])
virtual PMsrPlotList * GetMsrPlotList()
Definition: PMsrHandler.h:66
Bool_t SubtractEstimatedBkg()
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 SetFitRangeBin(const TString fitRange)
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
PDoubleVector fBackward
backward histo data
Definition: PRunAsymmetry.h:75
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
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
virtual void CalcTheory()
virtual Bool_t PrepareData()
#define ACCEL_PERIOD_PSI
Definition: PMusr.h:81
virtual Int_t GetAlphaParamNo()
Definition: PMusr.h:644
virtual Bool_t PrepareRRFViewData(PRawRunData *runData, UInt_t histoNo[2])
virtual void SetFitRange(Double_t dval, UInt_t idx)
Definition: PMusr.cpp:1068
Double_t fFitStartTime
fit start time
Definition: PRunBase.h:81
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 CalcChiSquare(const std::vector< Double_t > &par)
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
virtual void FilterTheo()
Definition: PRunBase.cpp:240
std::unique_ptr< PTheory > fTheory
theory needed to calculate chi-square
Definition: PRunBase.h:85
virtual Int_t GetPacking()
Definition: PMusr.h:591
virtual Double_t GetFitRange(UInt_t idx)
Definition: PMusr.cpp:1051
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
#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
Int_t fPacking
packing for this particular run. Either given in the RUN- or GLOBAL-block.
Definition: PRunAsymmetry.h:70
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