musrfit 1.10.0
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-2026 by Andreas Suter *
12 * andreas.suter@psi.ch *
13 * *
14 * This program is free software; you can redistribute it and/or modify *
15 * it under the terms of the GNU General Public License as published by *
16 * the Free Software Foundation; either version 2 of the License, or *
17 * (at your option) any later version. *
18 * *
19 * This program is distributed in the hope that it will be useful, *
20 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
21 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
22 * GNU General Public License for more details. *
23 * *
24 * You should have received a copy of the GNU General Public License *
25 * along with this program; if not, write to the *
26 * Free Software Foundation, Inc., *
27 * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
28 ***************************************************************************/
29
30#ifdef HAVE_CONFIG_H
31#include "config.h"
32#endif
33
34#ifdef HAVE_GOMP
35#include <omp.h>
36#endif
37
38#include <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//--------------------------------------------------------------------------
61{
62 fNoOfFitBins = 0;
63 fRRFPacking = -1;
64 fTheoAsData = false;
65
66 // the 2 following variables are need in case fit range is given in bins, and since
67 // the fit range can be changed in the command block, these variables need to be accessible
68 fGoodBins[0] = -1;
69 fGoodBins[1] = -1;
70}
71
72//--------------------------------------------------------------------------
73// Constructor
74//--------------------------------------------------------------------------
101PRunAsymmetryRRF::PRunAsymmetryRRF(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
102 PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
103{
104 // the 2 following variables are need in case fit range is given in bins, and since
105 // the fit range can be changed in the command block, these variables need to be accessible
106 fGoodBins[0] = -1;
107 fGoodBins[1] = -1;
108
109 fRRFPacking = fMsrInfo->GetMsrGlobal()->GetRRFPacking();
110 if (fRRFPacking == -1) { // this should NOT happen, somethin is severely wrong
111 std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **SEVERE ERROR**: Couldn't find any RRF packing information!";
112 std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
113 std::cerr << std::endl;
114 fValid = false;
115 return;
116 }
117
118 // check if alpha and/or beta is fixed --------------------
119
120 PMsrParamList *param = msrInfo->GetMsrParamList();
121
122 // check if alpha is given
123 if (fRunInfo->GetAlphaParamNo() == -1) { // no alpha given
124 std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **ERROR** no alpha parameter given! This is needed for an asymmetry fit!";
125 std::cerr << std::endl;
126 fValid = false;
127 return;
128 }
129 // check if alpha parameter is within proper bounds
130 if ((fRunInfo->GetAlphaParamNo() < 0) || (fRunInfo->GetAlphaParamNo() > static_cast<Int_t>(param->size()))) {
131 std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **ERROR** alpha parameter no = " << fRunInfo->GetAlphaParamNo();
132 std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
133 std::cerr << std::endl;
134 fValid = false;
135 return;
136 }
137 // check if alpha is fixed
138 Bool_t alphaFixedToOne = false;
139 if (((*param)[fRunInfo->GetAlphaParamNo()-1].fStep == 0.0) &&
140 ((*param)[fRunInfo->GetAlphaParamNo()-1].fValue == 1.0))
141 alphaFixedToOne = true;
142
143 // check if beta is given
144 Bool_t betaFixedToOne = false;
145 if (fRunInfo->GetBetaParamNo() == -1) { // no beta given hence assuming beta == 1
146 betaFixedToOne = true;
147 } else if ((fRunInfo->GetBetaParamNo() < 0) || (fRunInfo->GetBetaParamNo() > static_cast<Int_t>(param->size()))) { // check if beta parameter is within proper bounds
148 std::cerr << std::endl << ">> PRunAsymmetryRRF::PRunAsymmetryRRF(): **ERROR** beta parameter no = " << fRunInfo->GetBetaParamNo();
149 std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
150 std::cerr << std::endl;
151 fValid = false;
152 return;
153 } else { // check if beta is fixed
154 if (((*param)[fRunInfo->GetBetaParamNo()-1].fStep == 0.0) &&
155 ((*param)[fRunInfo->GetBetaParamNo()-1].fValue == 1.0))
156 betaFixedToOne = true;
157 }
158
159 // set fAlphaBetaTag
160 if (alphaFixedToOne && betaFixedToOne) // alpha == 1, beta == 1
161 fAlphaBetaTag = 1;
162 else if (!alphaFixedToOne && betaFixedToOne) // alpha != 1, beta == 1
163 fAlphaBetaTag = 2;
164 else if (alphaFixedToOne && !betaFixedToOne) // alpha == 1, beta != 1
165 fAlphaBetaTag = 3;
166 else
167 fAlphaBetaTag = 4;
168
169 // calculate fData
170 if (!PrepareData())
171 fValid = false;
172}
173
174//--------------------------------------------------------------------------
175// Destructor
176//--------------------------------------------------------------------------
184{
185 fForward.clear();
186 fForwardErr.clear();
187 fBackward.clear();
188 fBackwardErr.clear();
189}
190
191//--------------------------------------------------------------------------
192// CalcChiSquare (public)
193//--------------------------------------------------------------------------
214Double_t PRunAsymmetryRRF::CalcChiSquare(const std::vector<Double_t>& par)
215{
216 Double_t chisq = 0.0;
217 Double_t diff = 0.0;
218 Double_t asymFcnValue = 0.0;
219 Double_t a, b, f;
220
221 // calculate functions
222 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
223 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
224 }
225
226 // calculate chi square
227 Double_t time(1.0);
228 Int_t i;
229
230 // determine alpha/beta
231 switch (fAlphaBetaTag) {
232 case 1: // alpha == 1, beta == 1
233 a = 1.0;
234 b = 1.0;
235 break;
236 case 2: // alpha != 1, beta == 1
237 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
238 a = par[fRunInfo->GetAlphaParamNo()-1];
239 } else { // alpha is function
240 // get function number
241 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
242 // evaluate function
243 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
244 }
245 b = 1.0;
246 break;
247 case 3: // alpha == 1, beta != 1
248 a = 1.0;
249 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
250 b = par[fRunInfo->GetBetaParamNo()-1];
251 } else { // beta is a function
252 // get function number
253 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
254 // evaluate function
255 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
256 }
257 break;
258 case 4: // alpha != 1, beta != 1
259 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
260 a = par[fRunInfo->GetAlphaParamNo()-1];
261 } else { // alpha is function
262 // get function number
263 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
264 // evaluate function
265 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
266 }
267 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
268 b = par[fRunInfo->GetBetaParamNo()-1];
269 } else { // beta is a function
270 // get function number
271 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
272 // evaluate function
273 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
274 }
275 break;
276 default:
277 a = 1.0;
278 b = 1.0;
279 break;
280 }
281
282 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
283 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
284 // for a given set of parameters---which should be done outside of the parallelized loop.
285 // For all other functions it means a tiny and acceptable overhead.
286 asymFcnValue = fTheory->Func(time, par, fFuncValues);
287
288 #ifdef HAVE_GOMP
289 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
290 if (chunk < 10)
291 chunk = 10;
292 #pragma omp parallel for default(shared) private(i,time,diff,asymFcnValue,f) schedule(dynamic,chunk) reduction(+:chisq)
293 #endif
294 for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
295 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
296 f = fTheory->Func(time, par, fFuncValues);
297 asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
298 diff = fData.GetValue()->at(i) - asymFcnValue;
299 chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
300 }
301
302 return chisq;
303}
304
305//--------------------------------------------------------------------------
306// CalcChiSquareExpected (public)
307//--------------------------------------------------------------------------
321Double_t PRunAsymmetryRRF::CalcChiSquareExpected(const std::vector<Double_t>& par)
322{
323 return 0.0;
324}
325
326//--------------------------------------------------------------------------
327// CalcMaxLikelihood (public)
328//--------------------------------------------------------------------------
341Double_t PRunAsymmetryRRF::CalcMaxLikelihood(const std::vector<Double_t>& par)
342{
343 std::cout << std::endl << "PRunAsymmetryRRF::CalcMaxLikelihood(): not implemented yet ..." << std::endl;
344
345 return 1.0;
346}
347
348//--------------------------------------------------------------------------
349// GetNoOfFitBins (public)
350//--------------------------------------------------------------------------
361{
363
364 return fNoOfFitBins;
365}
366
367//--------------------------------------------------------------------------
368// SetFitRangeBin (public)
369//--------------------------------------------------------------------------
394void PRunAsymmetryRRF::SetFitRangeBin(const TString fitRange)
395{
396 TObjArray *tok = nullptr;
397 TObjString *ostr = nullptr;
398 TString str;
399 Ssiz_t idx = -1;
400 Int_t offset = 0;
401
402 tok = fitRange.Tokenize(" \t");
403
404 if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
405 // handle fgb+n0 entry
406 ostr = static_cast<TObjString*>(tok->At(1));
407 str = ostr->GetString();
408 // check if there is an offset present
409 idx = str.First("+");
410 if (idx != -1) { // offset present
411 str.Remove(0, idx+1);
412 if (str.IsFloat()) // if str is a valid number, convert is to an integer
413 offset = str.Atoi();
414 }
415 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
416
417 // handle lgb-n1 entry
418 ostr = static_cast<TObjString*>(tok->At(2));
419 str = ostr->GetString();
420 // check if there is an offset present
421 idx = str.First("-");
422 if (idx != -1) { // offset present
423 str.Remove(0, idx+1);
424 if (str.IsFloat()) // if str is a valid number, convert is to an integer
425 offset = str.Atoi();
426 }
427 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
428 } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
429 Int_t pos = 2*(fRunNo+1)-1;
430
431 if (pos + 1 >= tok->GetEntries()) {
432 std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
433 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
434 } else {
435 // handle fgb+n0 entry
436 ostr = static_cast<TObjString*>(tok->At(pos));
437 str = ostr->GetString();
438 // check if there is an offset present
439 idx = str.First("+");
440 if (idx != -1) { // offset present
441 str.Remove(0, idx+1);
442 if (str.IsFloat()) // if str is a valid number, convert is to an integer
443 offset = str.Atoi();
444 }
445 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
446
447 // handle lgb-n1 entry
448 ostr = static_cast<TObjString*>(tok->At(pos+1));
449 str = ostr->GetString();
450 // check if there is an offset present
451 idx = str.First("-");
452 if (idx != -1) { // offset present
453 str.Remove(0, idx+1);
454 if (str.IsFloat()) // if str is a valid number, convert is to an integer
455 offset = str.Atoi();
456 }
457 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
458 }
459 } else { // error
460 std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
461 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
462 }
463
464 // clean up
465 if (tok) {
466 delete tok;
467 }
468}
469
470//--------------------------------------------------------------------------
471// CalcNoOfFitBins (public)
472//--------------------------------------------------------------------------
486{
487 // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
488 fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
489 if (fStartTimeBin < 0)
490 fStartTimeBin = 0;
491 fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
492 if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
493 fEndTimeBin = fData.GetValue()->size();
494
497 else
498 fNoOfFitBins = 0;
499}
500
501//--------------------------------------------------------------------------
502// CalcTheory (protected)
503//--------------------------------------------------------------------------
525{
526 // feed the parameter vector
527 std::vector<Double_t> par;
528 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
529 for (UInt_t i=0; i<paramList->size(); i++)
530 par.push_back((*paramList)[i].fValue);
531
532 // calculate functions
533 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
534 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
535 }
536
537 // calculate asymmetry
538 Double_t asymFcnValue = 0.0;
539 Double_t a, b, f;
540 Double_t time;
541 for (UInt_t i=0; i<fData.GetValue()->size(); i++) {
542 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
543 switch (fAlphaBetaTag) {
544 case 1: // alpha == 1, beta == 1
545 asymFcnValue = fTheory->Func(time, par, fFuncValues);
546 break;
547 case 2: // alpha != 1, beta == 1
548 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
549 a = par[fRunInfo->GetAlphaParamNo()-1];
550 } else { // alpha is function
551 // get function number
552 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
553 // evaluate function
554 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
555 }
556 f = fTheory->Func(time, par, fFuncValues);
557 asymFcnValue = (f*(a+1.0)-(a-1.0))/((a+1.0)-f*(a-1.0));
558 break;
559 case 3: // alpha == 1, beta != 1
560 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
561 b = par[fRunInfo->GetBetaParamNo()-1];
562 } else { // beta is a function
563 // get function number
564 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
565 // evaluate function
566 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
567 }
568 f = fTheory->Func(time, par, fFuncValues);
569 asymFcnValue = f*(b+1.0)/(2.0-f*(b-1.0));
570 break;
571 case 4: // alpha != 1, beta != 1
572 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
573 a = par[fRunInfo->GetAlphaParamNo()-1];
574 } else { // alpha is function
575 // get function number
576 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
577 // evaluate function
578 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
579 }
580 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
581 b = par[fRunInfo->GetBetaParamNo()-1];
582 } else { // beta is a function
583 // get function number
584 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
585 // evaluate function
586 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
587 }
588 f = fTheory->Func(time, par, fFuncValues);
589 asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
590 break;
591 default:
592 asymFcnValue = 0.0;
593 break;
594 }
595 fData.AppendTheoryValue(asymFcnValue);
596 }
597
598 // clean up
599 par.clear();
600}
601
602//--------------------------------------------------------------------------
603// PrepareData (protected)
604//--------------------------------------------------------------------------
635{
636 if (!fValid)
637 return false;
638
639 // keep the Global block info
640 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
641
642 // get forward/backward histo from PRunDataHandler object ------------------------
643 // get the correct run
644 PRawRunData *runData = fRawData->GetRunData(*(fRunInfo->GetRunName()));
645 if (!runData) { // run not found
646 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
647 std::cerr << std::endl;
648 return false;
649 }
650
651 // keep the field from the meta-data from the data-file
652 fMetaData.fField = runData->GetField();
653
654 // keep the energy from the meta-data from the data-file
655 fMetaData.fEnergy = runData->GetEnergy();
656
657 // keep the temperature(s) from the meta-data from the data-file
658 for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
659 fMetaData.fTemp.push_back(runData->GetTemperature(i));
660
661 // collect histogram numbers
662 PUIntVector forwardHistoNo;
663 PUIntVector backwardHistoNo;
664 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
665 forwardHistoNo.push_back(fRunInfo->GetForwardHistoNo(i));
666
667 if (!runData->IsPresent(forwardHistoNo[i])) {
668 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **PANIC ERROR**:";
669 std::cerr << std::endl << ">> forwardHistoNo found = " << forwardHistoNo[i] << ", which is NOT present in the data file!?!?";
670 std::cerr << std::endl << ">> Will quit :-(";
671 std::cerr << std::endl;
672 // clean up
673 forwardHistoNo.clear();
674 backwardHistoNo.clear();
675 return false;
676 }
677 }
678 for (UInt_t i=0; i<fRunInfo->GetBackwardHistoNoSize(); i++) {
679 backwardHistoNo.push_back(fRunInfo->GetBackwardHistoNo(i));
680
681 if (!runData->IsPresent(backwardHistoNo[i])) {
682 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **PANIC ERROR**:";
683 std::cerr << std::endl << ">> backwardHistoNo found = " << backwardHistoNo[i] << ", which is NOT present in the data file!?!?";
684 std::cerr << std::endl << ">> Will quit :-(";
685 std::cerr << std::endl;
686 // clean up
687 forwardHistoNo.clear();
688 backwardHistoNo.clear();
689 return false;
690 }
691 }
692
693 // keep the time resolution in (us)
694 fTimeResolution = runData->GetTimeResolution()/1.0e3;
695 std::cout.precision(10);
696 std::cout << std::endl << ">> PRunAsymmetryRRF::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
697
698 // get all the proper t0's and addt0's for the current RUN block
699 if (!GetProperT0(runData, globalBlock, forwardHistoNo, backwardHistoNo)) {
700 return false;
701 }
702
703 // keep the histo of each group at this point (addruns handled below)
704 std::vector<PDoubleVector> forward, backward;
705 forward.resize(forwardHistoNo.size()); // resize to number of groups
706 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
707 forward[i].resize(runData->GetDataBin(forwardHistoNo[i])->size());
708 forward[i] = *runData->GetDataBin(forwardHistoNo[i]);
709 }
710 // check if a dead time correction has to be done
711 // this will be done automatically in the function itself, which also
712 // checks in the global and run section
713 DeadTimeCorrection(forward, forwardHistoNo);
714
715 backward.resize(backwardHistoNo.size()); // resize to number of groups
716 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
717 backward[i].resize(runData->GetDataBin(backwardHistoNo[i])->size());
718 backward[i] = *runData->GetDataBin(backwardHistoNo[i]);
719 }
720 // check if a dead time correction has to be done
721 // this will be done automatically in the function itself, which also
722 // checks in the global and run section
723 DeadTimeCorrection(forward, forwardHistoNo);
724
725 // check if addrun's are present, and if yes add data
726 // check if there are runs to be added to the current one
727 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
728 PRawRunData *addRunData;
729 std::vector<PDoubleVector> addForward, addBackward;
730 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
731 // get run to be added to the main one
732 addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
733 if (addRunData == nullptr) { // couldn't get run
734 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
735 std::cerr << std::endl;
736 return false;
737 }
738
739 // dead time correction handling
740 addForward.clear();
741 addForward.resize(forwardHistoNo.size());
742 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
743 addForward[j].resize(addRunData->GetDataBin(forwardHistoNo[j])->size());
744 addForward[j] = *addRunData->GetDataBin(forwardHistoNo[j]);
745 }
746 DeadTimeCorrection(addForward, forwardHistoNo);
747 addBackward.clear();
748 addBackward.resize(backwardHistoNo.size());
749 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
750 addBackward[j].resize(addRunData->GetDataBin(backwardHistoNo[j])->size());
751 addBackward[j] = *addRunData->GetDataBin(backwardHistoNo[j]);
752 }
753 DeadTimeCorrection(addBackward, backwardHistoNo);
754
755 // add forward run
756 UInt_t addRunSize;
757 for (UInt_t k=0; k<forwardHistoNo.size(); k++) { // fill each group
758 addRunSize = addRunData->GetDataBin(forwardHistoNo[k])->size();
759 for (UInt_t j=0; j<addRunData->GetDataBin(forwardHistoNo[k])->size(); j++) { // loop over the bin indices
760 // make sure that the index stays in the proper range
761 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)) {
762 forward[k][j] += addForward[k][j+(Int_t)fAddT0s[i-1][2*k]-(Int_t)fT0s[2*k]];
763 }
764 }
765 }
766
767 // add backward run
768 for (UInt_t k=0; k<backwardHistoNo.size(); k++) { // fill each group
769 addRunSize = addRunData->GetDataBin(backwardHistoNo[k])->size();
770 for (UInt_t j=0; j<addRunData->GetDataBin(backwardHistoNo[k])->size(); j++) { // loop over the bin indices
771 // make sure that the index stays in the proper range
772 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)) {
773 backward[k][j] += addBackward[k][j+(Int_t)fAddT0s[i-1][2*k+1]-(Int_t)fT0s[2*k+1]];
774 }
775 }
776 }
777 }
778 }
779
780 // set forward/backward histo data of the first group
781 fForward.resize(forward[0].size());
782 fBackward.resize(backward[0].size());
783 for (UInt_t i=0; i<fForward.size(); i++) {
784 fForward[i] = forward[0][i];
785 fBackward[i] = backward[0][i];
786 }
787
788 // group histograms, add all the remaining forward histograms of the group
789 for (UInt_t i=1; i<forwardHistoNo.size(); i++) { // loop over the groupings
790 for (UInt_t j=0; j<runData->GetDataBin(forwardHistoNo[i])->size(); j++) { // loop over the bin indices
791 // make sure that the index stays within proper range
792 if (((Int_t)j+fT0s[2*i]-fT0s[0] >= 0) && (j+fT0s[2*i]-fT0s[0] < runData->GetDataBin(forwardHistoNo[i])->size())) {
793 fForward[j] += forward[i][j+(Int_t)fT0s[2*i]-(Int_t)fT0s[0]];
794 }
795 }
796 }
797
798 // group histograms, add all the remaining backward histograms of the group
799 for (UInt_t i=1; i<backwardHistoNo.size(); i++) { // loop over the groupings
800 for (UInt_t j=0; j<runData->GetDataBin(backwardHistoNo[i])->size(); j++) { // loop over the bin indices
801 // make sure that the index stays within proper range
802 if ((j+fT0s[2*i+1]-fT0s[1] >= 0) && (j+fT0s[2*i+1]-fT0s[1] < runData->GetDataBin(backwardHistoNo[i])->size())) {
803 fBackward[j] += backward[i][j+(Int_t)fT0s[2*i+1]-(Int_t)fT0s[1]];
804 }
805 }
806 }
807
808 // subtract background from histogramms ------------------------------------------
809 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
810 if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
812 return false;
813 } else { // no background given to do the job, try to estimate it
814 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
815 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
816 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.1), 2);
817 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.6), 3);
818 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **WARNING** Neither fix background nor background bins are given!";
819 std::cerr << std::endl << ">> Will try the following:";
820 std::cerr << std::endl << ">> forward: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
821 std::cerr << std::endl << ">> backward: bkg start = " << fRunInfo->GetBkgRange(2) << ", bkg end = " << fRunInfo->GetBkgRange(3);
822 std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
823 std::cerr << std::endl;
825 return false;
826 }
827 } else { // fixed background given
828 if (!SubtractFixBkg())
829 return false;
830 }
831
832 UInt_t histoNo[2] = {forwardHistoNo[0], backwardHistoNo[0]};
833
834 // get the data range (fgb/lgb) for the current RUN block
835 if (!GetProperDataRange(runData, histoNo)) {
836 return false;
837 }
838
839 // get the fit range for the current RUN block
840 GetProperFitRange(globalBlock);
841
842 // everything looks fine, hence fill data set
843 Bool_t status;
844 switch(fHandleTag) {
845 case kFit:
847 break;
848 case kView:
849 status = PrepareViewData(runData, histoNo);
850 break;
851 default:
852 status = false;
853 break;
854 }
855
856 // clean up
857 forwardHistoNo.clear();
858 backwardHistoNo.clear();
859
860 return status;
861}
862
863//--------------------------------------------------------------------------
864// SubtractFixBkg (private)
865//--------------------------------------------------------------------------
883{
884 Double_t dval;
885 for (UInt_t i=0; i<fForward.size(); i++) {
886 // keep the error, and make sure that the bin is NOT empty
887 if (fForward[i] != 0.0)
888 dval = TMath::Sqrt(fForward[i]);
889 else
890 dval = 1.0;
891 fForwardErr.push_back(dval);
892 fForward[i] -= fRunInfo->GetBkgFix(0);
893
894 // keep the error, and make sure that the bin is NOT empty
895 if (fBackward[i] != 0.0)
896 dval = TMath::Sqrt(fBackward[i]);
897 else
898 dval = 1.0;
899 fBackwardErr.push_back(dval);
900 fBackward[i] -= fRunInfo->GetBkgFix(1);
901 }
902
903 return true;
904}
905
906//--------------------------------------------------------------------------
907// SubtractEstimatedBkg (private)
908//--------------------------------------------------------------------------
927{
928 Double_t beamPeriod = 0.0;
929
930 // check if data are from PSI, RAL, or TRIUMF
931 if (fRunInfo->GetInstitute()->Contains("psi"))
932 beamPeriod = ACCEL_PERIOD_PSI;
933 else if (fRunInfo->GetInstitute()->Contains("ral"))
934 beamPeriod = ACCEL_PERIOD_RAL;
935 else if (fRunInfo->GetInstitute()->Contains("triumf"))
936 beamPeriod = ACCEL_PERIOD_TRIUMF;
937 else
938 beamPeriod = 0.0;
939
940 // check if start and end are in proper order
941 UInt_t start[2] = {static_cast<UInt_t>(fRunInfo->GetBkgRange(0)), static_cast<UInt_t>(fRunInfo->GetBkgRange(2))};
942 UInt_t end[2] = {static_cast<UInt_t>(fRunInfo->GetBkgRange(1)), static_cast<UInt_t>(fRunInfo->GetBkgRange(3))};
943 for (UInt_t i=0; i<2; i++) {
944 if (end[i] < start[i]) {
945 std::cout << std::endl << "PRunAsymmetryRRF::SubtractEstimatedBkg(): end = " << end[i] << " > start = " << start[i] << "! Will swap them!";
946 UInt_t keep = end[i];
947 end[i] = start[i];
948 start[i] = keep;
949 }
950 }
951
952 // calculate proper background range
953 for (UInt_t i=0; i<2; i++) {
954 if (beamPeriod != 0.0) {
955 Double_t timeBkg = static_cast<Double_t>(end[i]-start[i])*fTimeResolution; // length of the background intervall in time
956 UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
957 // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
958 end[i] = start[i] + static_cast<UInt_t>((fullCycles*beamPeriod)/fTimeResolution);
959 std::cout << "PRunAsymmetryRRF::SubtractEstimatedBkg(): Background " << start[i] << ", " << end[i] << std::endl;
960 if (end[i] == start[i])
961 end[i] = fRunInfo->GetBkgRange(2*i+1);
962 }
963 }
964
965 // check if start is within histogram bounds
966 if ((start[0] >= fForward.size()) || (start[1] >= fBackward.size())) {
967 std::cerr << std::endl << ">> PRunAsymmetryRRF::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
968 std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
969 std::cerr << std::endl << ">> background start (f/b) = (" << start[0] << "/" << start[1] << ").";
970 return false;
971 }
972
973 // check if end is within histogram bounds
974 if ((end[0] >= fForward.size()) || (end[1] >= fBackward.size())) {
975 std::cerr << std::endl << ">> PRunAsymmetryRRF::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
976 std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
977 std::cerr << std::endl << ">> background end (f/b) = (" << end[0] << "/" << end[1] << ").";
978 return false;
979 }
980
981 // calculate background
982 Double_t bkg[2] = {0.0, 0.0};
983 Double_t errBkg[2] = {0.0, 0.0};
984
985 // forward
986 for (UInt_t i=start[0]; i<=end[0]; i++)
987 bkg[0] += fForward[i];
988 errBkg[0] = TMath::Sqrt(bkg[0])/(end[0] - start[0] + 1);
989 bkg[0] /= static_cast<Double_t>(end[0] - start[0] + 1);
990 std::cout << std::endl << ">> estimated forward histo background: " << bkg[0];
991
992 // backward
993 for (UInt_t i=start[1]; i<=end[1]; i++)
994 bkg[1] += fBackward[i];
995 errBkg[1] = TMath::Sqrt(bkg[1])/(end[1] - start[1] + 1);
996 bkg[1] /= static_cast<Double_t>(end[1] - start[1] + 1);
997 std::cout << std::endl << ">> estimated backward histo background: " << bkg[1] << std::endl;
998
999 // correct error for forward, backward
1000 Double_t errVal = 0.0;
1001 for (UInt_t i=0; i<fForward.size(); i++) {
1002 if (fForward[i] > 0.0)
1003 errVal = TMath::Sqrt(fForward[i]+errBkg[0]*errBkg[0]);
1004 else
1005 errVal = 1.0;
1006 fForwardErr.push_back(errVal);
1007 if (fBackward[i] > 0.0)
1008 errVal = TMath::Sqrt(fBackward[i]+errBkg[1]*errBkg[1]);
1009 else
1010 errVal = 1.0;
1011 fBackwardErr.push_back(errVal);
1012 }
1013
1014 // subtract background from data
1015 for (UInt_t i=0; i<fForward.size(); i++) {
1016 fForward[i] -= bkg[0];
1017 fBackward[i] -= bkg[1];
1018 }
1019
1020 fRunInfo->SetBkgEstimated(bkg[0], 0);
1021 fRunInfo->SetBkgEstimated(bkg[1], 1);
1022
1023 return true;
1024}
1025
1026//--------------------------------------------------------------------------
1027// PrepareFitData (protected)
1028//--------------------------------------------------------------------------
1051{
1052 // transform raw histo data. At this point, the raw data are already background corrected.
1053
1054 // 1st: form the asymmetry of the original data
1055
1056 // forward and backward detectors might have different fgb-t0 offset. Take the maximum of both.
1057 Int_t fgbOffset = fGoodBins[0]-static_cast<Int_t>(fT0s[0]);
1058 if (fgbOffset < fGoodBins[2]-static_cast<Int_t>(fT0s[1]))
1059 fgbOffset = fGoodBins[2]-static_cast<Int_t>(fT0s[1]);
1060 // last good bin (lgb) is the minimum of forward/backward lgb
1061 Int_t lgb_offset = fGoodBins[1]-static_cast<Int_t>(fT0s[0])+fgbOffset;
1062 if (lgb_offset < fGoodBins[3]-static_cast<Int_t>(fT0s[1])+fgbOffset)
1063 lgb_offset = fGoodBins[3]-static_cast<Int_t>(fT0s[1])+fgbOffset;
1064
1065 Int_t fgb = static_cast<Int_t>(fT0s[0])+fgbOffset;
1066 if (fgb < 0)
1067 fgb=0;
1068 Int_t lgb = fgb + lgb_offset;
1069 if (lgb > fForward.size())
1070 lgb = fForward.size();
1071 if (lgb > fBackward.size())
1072 lgb = fBackward.size();
1073 Int_t dt0 = static_cast<Int_t>(fT0s[0])-static_cast<Int_t>(fT0s[1]);
1074
1075 PDoubleVector asym;
1076 PDoubleVector asymErr;
1077 Double_t asymVal, asymValErr;
1078 Double_t ff, bb, eff, ebb;
1079 Int_t idx=0;
1080 for (Int_t i=fgb; i<lgb; i++) {
1081 ff = fForward[i];
1082 idx = i-dt0;
1083 if (idx < 0)
1084 idx = 0;
1085 if (idx >= fBackward.size())
1086 idx = fBackward.size()-1;
1087 bb = fBackward[idx];
1088 if (ff+bb != 0.0)
1089 asymVal = (ff-bb)/(ff+bb);
1090 else
1091 asymVal = 0.0;
1092 asym.push_back(asymVal);
1093 eff = fForwardErr[i];
1094 ebb = fBackwardErr[i-dt0];
1095 if ((asymVal != 0.0) && (ff+bb) > 0.0)
1096 asymValErr = 2.0/pow((ff+bb),2.0)*sqrt(bb*bb*eff*eff+ff*ff*ebb*ebb);
1097 else
1098 asymValErr = 1.0;
1099 asymErr.push_back(asymValErr);
1100 }
1101
1102 // 2nd: a_rrf = a * 2*cos(w_rrf*t + phi_rrf)
1103 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
1104 Double_t wRRF = globalBlock->GetRRFFreq("Mc");
1105 Double_t phaseRRF = globalBlock->GetRRFPhase()*TMath::TwoPi()/180.0;
1106
1107 Double_t startTime = fTimeResolution * static_cast<Double_t>(fgbOffset);
1108 Double_t time=0.0;
1109 for (UInt_t i=0; i<asym.size(); i++) {
1110 time = startTime + i*fTimeResolution;
1111 asym[i] *= 2.0*cos(wRRF*time+phaseRRF);
1112 }
1113
1114 // 3rd: rrf packing
1115 PDoubleVector asymRRF;
1116 asymVal = 0.0;
1117 for (UInt_t i=0; i<asym.size(); i++) {
1118 if ((i+1) % fRRFPacking == 0) {
1119 asymRRF.push_back(asymVal/fRRFPacking);
1120 asymVal = 0.0;
1121 }
1122 asymVal += asym[i];
1123 }
1124
1125 // 4th: rrf packing error
1126 PDoubleVector asymRRFErr;
1127 asymValErr = 0.0;
1128 for (UInt_t i=0; i<asymErr.size(); i++) {
1129 if ((i+1) % fRRFPacking == 0) {
1130 asymRRFErr.push_back(sqrt(2.0*asymValErr)/fRRFPacking); // factor of two is needed due to the rescaling
1131 asymValErr = 0.0;
1132 }
1133 asymValErr += asymErr[i]*asymErr[i];
1134 }
1135
1136 fData.SetDataTimeStart(startTime+fTimeResolution*(static_cast<Double_t>(fRRFPacking-1)/2.0));
1137 fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(fRRFPacking));
1138
1139 for (UInt_t i=0; i<asymRRF.size(); i++) {
1140 fData.AppendValue(asymRRF[i]);
1141 fData.AppendErrorValue(asymRRFErr[i]);
1142 }
1143
1145
1146 // clean up
1147 fForward.clear();
1148 fForwardErr.clear();
1149 fBackward.clear();
1150 fBackwardErr.clear();
1151 asym.clear();
1152 asymErr.clear();
1153 asymRRF.clear();
1154 asymRRFErr.clear();
1155
1156 return true;
1157}
1158
1159//--------------------------------------------------------------------------
1160// PrepareViewData (protected)
1161//--------------------------------------------------------------------------
1187Bool_t PRunAsymmetryRRF::PrepareViewData(PRawRunData* runData, UInt_t histoNo[2])
1188{
1189 // feed the parameter vector
1190 std::vector<Double_t> par;
1191 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1192 for (UInt_t i=0; i<paramList->size(); i++)
1193 par.push_back((*paramList)[i].fValue);
1194
1195 // first get start data, end data, and t0
1196 Int_t start[2] = {fGoodBins[0], fGoodBins[2]};
1197 Int_t end[2] = {fGoodBins[1], fGoodBins[3]};
1198 Int_t t0[2] = {static_cast<Int_t>(fT0s[0]), static_cast<Int_t>(fT0s[1])};
1199
1200 // check if the data ranges and t0's between forward/backward are compatible
1201 Int_t fgb[2];
1202 if (start[0]-t0[0] != start[1]-t0[1]) { // wrong fgb aligning
1203 if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
1204 fgb[0] = start[0];
1205 fgb[1] = t0[1] + start[0]-t0[0];
1206 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **WARNING** needed to shift backward fgb from ";
1207 std::cerr << start[1] << " to " << fgb[1] << std::endl;
1208 } else {
1209 fgb[0] = t0[0] + start[1]-t0[1];
1210 fgb[1] = start[1];
1211 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **WARNING** needed to shift forward fgb from ";
1212 std::cerr << start[0] << " to " << fgb[0] << std::endl;
1213 }
1214 } else { // fgb aligning is correct
1215 fgb[0] = start[0];
1216 fgb[1] = start[1];
1217 }
1218 start[0] = fgb[0];
1219 start[1] = fgb[1];
1220
1221 // make sure that there are equal number of bins in forward and backward
1222 UInt_t noOfBins0 = runData->GetDataBin(histoNo[0])->size()-start[0];
1223 UInt_t noOfBins1 = runData->GetDataBin(histoNo[1])->size()-start[1];
1224 if (noOfBins0 > noOfBins1)
1225 noOfBins0 = noOfBins1;
1226 end[0] = start[0] + noOfBins0;
1227 end[1] = start[1] + noOfBins0;
1228
1229 // check if start, end, and t0 make any sense
1230 for (UInt_t i=0; i<2; i++) {
1231 // 1st check if start is within proper bounds
1232 if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1233 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
1234 std::cerr << std::endl;
1235 return false;
1236 }
1237 // 2nd check if end is within proper bounds
1238 if ((end[i] < 0) || (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1239 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
1240 std::cerr << std::endl;
1241 return false;
1242 }
1243 // 3rd check if t0 is within proper bounds
1244 if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1245 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareViewData(): **ERROR** t0 data bin doesn't make any sense!";
1246 std::cerr << std::endl;
1247 return false;
1248 }
1249 }
1250
1251 // check if forward and backward histo have the same size, otherwise take the minimum size
1252 UInt_t noOfBins = fForward.size();
1253 if (noOfBins > fBackward.size()) {
1254 noOfBins = fBackward.size();
1255 }
1256
1257 // form asymmetry including error propagation
1258 Double_t asym, error;
1259 Double_t f, b, ef, eb, alpha = 1.0, beta = 1.0;
1260
1261 // get the proper alpha and beta
1262 switch (fAlphaBetaTag) {
1263 case 1: // alpha == 1, beta == 1
1264 alpha = 1.0;
1265 beta = 1.0;
1266 break;
1267 case 2: // alpha != 1, beta == 1
1268 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1269 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1270 } else { // alpha is function
1271 // get function number
1272 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1273 // evaluate function
1274 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1275 }
1276 beta = 1.0;
1277 break;
1278 case 3: // alpha == 1, beta != 1
1279 alpha = 1.0;
1280 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1281 beta = par[fRunInfo->GetBetaParamNo()-1];
1282 } else { // beta is a function
1283 // get function number
1284 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1285 // evaluate function
1286 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1287 }
1288 break;
1289 case 4: // alpha != 1, beta != 1
1290 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1291 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1292 } else { // alpha is function
1293 // get function number
1294 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1295 // evaluate function
1296 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1297 }
1298 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1299 beta = par[fRunInfo->GetBetaParamNo()-1];
1300 } else { // beta is a function
1301 // get function number
1302 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1303 // evaluate function
1304 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1305 }
1306 break;
1307 default:
1308 break;
1309 }
1310
1311 PDoubleVector asymVec, asymErr;
1312 Int_t dtBin = start[1]-start[0];
1313 for (Int_t i=start[0]; i<end[0]; i++) {
1314 // to make the formulae more readable
1315 f = fForward[i];
1316 b = fBackward[i+dtBin];
1317 ef = fForwardErr[i];
1318 eb = fBackwardErr[i+dtBin];
1319 // check that there are indeed bins
1320 if (f+b != 0.0)
1321 asym = (alpha*f-b) / (alpha*beta*f+b);
1322 else
1323 asym = 0.0;
1324 asymVec.push_back(asym);
1325 // calculate the error
1326 if (f+b != 0.0)
1327 error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1328 else
1329 error = 1.0;
1330 asymErr.push_back(error);
1331 }
1332
1333 // RRF transform
1334 // a_rrf = a * 2*cos(w_rrf*t + phi_rrf)
1335 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
1336 Double_t wRRF = globalBlock->GetRRFFreq("Mc");
1337 Double_t phaseRRF = globalBlock->GetRRFPhase()*TMath::TwoPi()/180.0;
1338 Double_t startTime=fTimeResolution*(static_cast<Double_t>(start[0])-t0[0]);
1339 Double_t time = 0.0;
1340 for (UInt_t i=0; i<asymVec.size(); i++) {
1341 time = startTime + i * fTimeResolution;
1342 asymVec[i] *= 2.0*cos(wRRF*time+phaseRRF); // factor of 2 needed to keep the asymmetry
1343 }
1344
1345 Double_t dval = 0.0;
1346 PDoubleVector asymRRF;
1347 for (UInt_t i=0; i<asymVec.size(); i++) {
1348 if ((i+1) % fRRFPacking == 0) {
1349 asymRRF.push_back(dval/fRRFPacking);
1350 dval = 0.0;
1351 }
1352 dval += asymVec[i];
1353 }
1354
1355 // RRF packing error
1356 PDoubleVector asymRRFErr;
1357 dval = 0.0;
1358 for (UInt_t i=0; i<asymErr.size(); i++) {
1359 if ((i+1) % fRRFPacking == 0) {
1360 asymRRFErr.push_back(sqrt(2.0*dval)/fRRFPacking); // factor of two is needed due to the rescaling
1361 dval = 0.0;
1362 }
1363 dval += asymErr[i]*asymErr[i];
1364 }
1365
1366 fData.SetDataTimeStart(startTime+fTimeResolution*(static_cast<Double_t>(fRRFPacking-1)/2.0));
1367 fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(fRRFPacking));
1368
1369 for (UInt_t i=0; i<asymRRF.size(); i++) {
1370 fData.AppendValue(asymRRF[i]);
1371 fData.AppendErrorValue(asymRRFErr[i]);
1372 }
1373
1375
1376 // clean up
1377 fForward.clear();
1378 fForwardErr.clear();
1379 fBackward.clear();
1380 fBackwardErr.clear();
1381 asymVec.clear();
1382 asymErr.clear();
1383 asymRRF.clear();
1384 asymRRFErr.clear();
1385
1386 // fill theory vector for kView
1387 // calculate functions
1388 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1389 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
1390 }
1391
1392 // calculate theory
1393 UInt_t size = runData->GetDataBin(histoNo[0])->size();
1394 Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1395 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1396 if (fTheoAsData) { // calculate theory only at the data points
1397 fData.SetTheoryTimeStep(fData.GetDataTimeStep());
1398 } else {
1399 // finer binning for the theory (8 times as many points = factor)
1400 size *= factor;
1401 fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1402 }
1403
1404 for (UInt_t i=0; i<size; i++) {
1405 time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1406 dval = fTheory->Func(time, par, fFuncValues);
1407 if (fabs(dval) > 10.0) { // dirty hack needs to be fixed!!
1408 dval = 0.0;
1409 }
1410 fData.AppendTheoryValue(dval);
1411 }
1412
1413 // clean up
1414 par.clear();
1415
1416 return true;
1417}
1418
1419//--------------------------------------------------------------------------
1420// GetProperT0 (private)
1421//--------------------------------------------------------------------------
1440Bool_t PRunAsymmetryRRF::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHistoNo, PUIntVector &backwardHistoNo)
1441{
1442 // feed all T0's
1443 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1444 fT0s.clear();
1445 // this strange fT0 size estimate is needed in case #forw histos != #back histos
1446 size_t size = 2*forwardHistoNo.size();
1447 if (backwardHistoNo.size() > forwardHistoNo.size())
1448 size = 2*backwardHistoNo.size();
1449 fT0s.resize(size);
1450 for (UInt_t i=0; i<fT0s.size(); i++) {
1451 fT0s[i] = -1.0;
1452 }
1453
1454 // fill in the T0's from the msr-file (if present)
1455 for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1456 fT0s[i] = fRunInfo->GetT0Bin(i);
1457 }
1458
1459 // fill in the missing T0's from the GLOBAL block section (if present)
1460 for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1461 if (fT0s[i] == -1.0) { // i.e. not given in the RUN block section
1462 fT0s[i] = globalBlock->GetT0Bin(i);
1463 }
1464 }
1465
1466 // fill in the missing T0's from the data file, if not already present in the msr-file
1467 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1468 if (fT0s[2*i] == -1.0) // i.e. not present in the msr-file, try the data file
1469 if (runData->GetT0Bin(forwardHistoNo[i]) > 0.0) {
1470 fT0s[2*i] = runData->GetT0Bin(forwardHistoNo[i]);
1471 fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1472 }
1473 }
1474 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1475 if (fT0s[2*i+1] == -1.0) // i.e. not present in the msr-file, try the data file
1476 if (runData->GetT0Bin(backwardHistoNo[i]) > 0.0) {
1477 fT0s[2*i+1] = runData->GetT0Bin(backwardHistoNo[i]);
1478 fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1479 }
1480 }
1481
1482 // 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
1483 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1484 if (fT0s[2*i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1485 fT0s[2*i] = runData->GetT0BinEstimated(forwardHistoNo[i]);
1486 fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1487
1488 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1489 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1490 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(forwardHistoNo[i]);
1491 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1492 std::cerr << std::endl;
1493 }
1494 }
1495 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1496 if (fT0s[2*i+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1497 fT0s[2*i+1] = runData->GetT0BinEstimated(backwardHistoNo[i]);
1498 fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1499
1500 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1501 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1502 std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[i]);
1503 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1504 std::cerr << std::endl;
1505 }
1506 }
1507
1508 // check if t0 is within proper bounds
1509 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1510 if ((fT0s[2*i] < 0) || (fT0s[2*i] > static_cast<Int_t>(runData->GetDataBin(forwardHistoNo[i])->size()))) {
1511 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **ERROR** t0 data bin (" << fT0s[2*i] << ") doesn't make any sense!";
1512 std::cerr << std::endl << ">> forwardHistoNo " << forwardHistoNo[i];
1513 std::cerr << std::endl;
1514 return false;
1515 }
1516 }
1517 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1518 if ((fT0s[2*i+1] < 0) || (fT0s[2*i+1] > static_cast<Int_t>(runData->GetDataBin(backwardHistoNo[i])->size()))) {
1519 std::cerr << std::endl << ">> PRunAsymmetryRRF::PrepareData(): **ERROR** t0 data bin (" << fT0s[2*i+1] << ") doesn't make any sense!";
1520 std::cerr << std::endl << ">> backwardHistoNo " << backwardHistoNo[i];
1521 std::cerr << std::endl;
1522 return false;
1523 }
1524 }
1525
1526 // check if addrun's are present, and if yes add the necessary t0's
1527 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1528 PRawRunData *addRunData;
1529 fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1530 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1531 // get run to be added to the main one
1532 addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
1533 if (addRunData == nullptr) { // couldn't get run
1534 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1535 std::cerr << std::endl;
1536 return false;
1537 }
1538
1539 // feed all T0's
1540 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1541 fAddT0s[i-1].clear();
1542 fAddT0s[i-1].resize(2*forwardHistoNo.size());
1543 for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1544 fAddT0s[i-1][j] = -1.0;
1545 }
1546
1547 // fill in the T0's from the msr-file (if present)
1548 for (Int_t j=0; j<fRunInfo->GetAddT0BinSize(i-1); j++) {
1549 fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1, j);
1550 }
1551
1552 // fill in the T0's from the data file, if not already present in the msr-file
1553 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1554 if (fAddT0s[i-1][2*j] == -1.0) // i.e. not present in the msr-file, try the data file
1555 if (addRunData->GetT0Bin(forwardHistoNo[j]) > 0.0) {
1556 fAddT0s[i-1][2*j] = addRunData->GetT0Bin(forwardHistoNo[j]);
1557 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1558 }
1559 }
1560 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1561 if (fAddT0s[i-1][2*j+1] == -1.0) // i.e. not present in the msr-file, try the data file
1562 if (addRunData->GetT0Bin(backwardHistoNo[j]) > 0.0) {
1563 fAddT0s[i-1][2*j+1] = addRunData->GetT0Bin(backwardHistoNo[j]);
1564 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1565 }
1566 }
1567
1568 // 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
1569 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1570 if (fAddT0s[i-1][2*j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1571 fAddT0s[i-1][2*j] = addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1572 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1573
1574 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1575 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1576 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1577 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1578 std::cerr << std::endl;
1579 }
1580 }
1581 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1582 if (fAddT0s[i-1][2*j+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1583 fAddT0s[i-1][2*j+1] = addRunData->GetT0BinEstimated(backwardHistoNo[j]);
1584 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1585
1586 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1587 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1588 std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[j]);
1589 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1590 std::cerr << std::endl;
1591 }
1592 }
1593 }
1594 }
1595
1596 return true;
1597}
1598
1599//--------------------------------------------------------------------------
1600// GetProperDataRange (private)
1601//--------------------------------------------------------------------------
1615Bool_t PRunAsymmetryRRF::GetProperDataRange(PRawRunData* runData, UInt_t histoNo[2])
1616{
1617 // first get start/end data
1618 Int_t start[2] = {fRunInfo->GetDataRange(0), fRunInfo->GetDataRange(2)};
1619 Int_t end[2] = {fRunInfo->GetDataRange(1), fRunInfo->GetDataRange(3)};
1620 // check if data range has been provided in the RUN block. If not, try the GLOBAL block
1621 if (start[0] == -1) {
1622 start[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1623 }
1624 if (start[1] == -1) {
1625 start[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(2);
1626 }
1627 if (end[0] == -1) {
1628 end[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1629 }
1630 if (end[1] == -1) {
1631 end[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(3);
1632 }
1633
1634 Double_t t0[2] = {fT0s[0], fT0s[1]};
1635 Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution); // needed in case first good bin is not given, default = 10ns
1636
1637 // check if data range has been provided, and if not try to estimate them
1638 if (start[0] < 0) {
1639 start[0] = static_cast<Int_t>(t0[0])+offset;
1640 fRunInfo->SetDataRange(start[0], 0);
1641 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[0] << ".";
1642 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1643 std::cerr << std::endl;
1644 }
1645 if (start[1] < 0) {
1646 start[1] = static_cast<Int_t>(t0[1])+offset;
1647 fRunInfo->SetDataRange(start[1], 2);
1648 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[1] << ".";
1649 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1650 std::cerr << std::endl;
1651 }
1652 if (end[0] < 0) {
1653 end[0] = runData->GetDataBin(histoNo[0])->size();
1654 fRunInfo->SetDataRange(end[0], 1);
1655 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range end = " << end[0] << ".";
1656 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1657 std::cerr << std::endl;
1658 }
1659 if (end[1] < 0) {
1660 end[1] = runData->GetDataBin(histoNo[1])->size();
1661 fRunInfo->SetDataRange(end[1], 3);
1662 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range end = " << end[1] << ".";
1663 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1664 std::cerr << std::endl;
1665 }
1666
1667 // check if start, end, and t0 make any sense
1668 // 1st check if start and end are in proper order
1669 for (UInt_t i=0; i<2; i++) {
1670 if (end[i] < start[i]) { // need to swap them
1671 Int_t keep = end[i];
1672 end[i] = start[i];
1673 start[i] = keep;
1674 }
1675 // 2nd check if start is within proper bounds
1676 if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1677 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **ERROR** start data bin doesn't make any sense!";
1678 std::cerr << std::endl;
1679 return false;
1680 }
1681 // 3rd check if end is within proper bounds
1682 if (end[i] < 0) {
1683 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **ERROR** end data bin (" << end[i] << ") doesn't make any sense!";
1684 std::cerr << std::endl;
1685 return false;
1686 }
1687 if (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())) {
1688 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **WARNING** end data bin (" << end[i] << ") > histo length (" << (Int_t)runData->GetDataBin(histoNo[i])->size() << ").";
1689 std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
1690 std::cerr << std::endl;
1691 end[i] = static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())-1;
1692 }
1693 // 4th check if t0 is within proper bounds
1694 if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1695 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange(): **ERROR** t0 data bin doesn't make any sense!";
1696 std::cerr << std::endl;
1697 return false;
1698 }
1699 }
1700
1701 // check that start-t0 is the same for forward as for backward, otherwise take max(start[i]-t0[i])
1702 if (fabs(static_cast<Double_t>(start[0])-t0[0]) > fabs(static_cast<Double_t>(start[1])-t0[1])){
1703 start[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(start[0]) - t0[0]);
1704 end[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(end[0]) - t0[0]);
1705 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange **WARNING** needed to shift backward data range.";
1706 std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(2) << ", " << fRunInfo->GetDataRange(3);
1707 std::cerr << std::endl << ">> used : " << start[1] << ", " << end[1];
1708 std::cerr << std::endl;
1709 }
1710 if (fabs(static_cast<Double_t>(start[0])-t0[0]) < fabs(static_cast<Double_t>(start[1])-t0[1])){
1711 start[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(start[1]) - t0[1]);
1712 end[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(end[1]) - t0[1]);
1713 std::cerr << std::endl << ">> PRunAsymmetryRRF::GetProperDataRange **WARNING** needed to shift forward data range.";
1714 std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(0) << ", " << fRunInfo->GetDataRange(1);
1715 std::cerr << std::endl << ">> used : " << start[0] << ", " << end[0];
1716 std::cerr << std::endl;
1717 }
1718
1719 // keep good bins for potential latter use
1720 fGoodBins[0] = start[0];
1721 fGoodBins[1] = end[0];
1722 fGoodBins[2] = start[1];
1723 fGoodBins[3] = end[1];
1724
1725 // make sure that fGoodBins are in proper range for fForward and fBackward
1726 if (fGoodBins[0] < 0)
1727 fGoodBins[0]=0;
1728 if (fGoodBins[1] > fForward.size()) {
1729 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift forward lgb,";
1730 std::cerr << std::endl << ">> from " << fGoodBins[1] << " to " << fForward.size()-1 << std::endl;
1731 fGoodBins[1]=fForward.size()-1;
1732 }
1733 if (fGoodBins[2] < 0)
1734 fGoodBins[2]=0;
1735 if (fGoodBins[3] > fBackward.size()) {
1736 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift backward lgb,";
1737 std::cerr << std::endl << ">> from " << fGoodBins[1] << " to " << fForward.size()-1 << std::endl;
1738 fGoodBins[3]=fBackward.size()-1;
1739 }
1740
1741 return true;
1742}
1743
1744//--------------------------------------------------------------------------
1745// GetProperFitRange (private)
1746//--------------------------------------------------------------------------
1760{
1761 // set fit start/end time; first check RUN Block
1762 fFitStartTime = fRunInfo->GetFitRange(0);
1763 fFitEndTime = fRunInfo->GetFitRange(1);
1764 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1765 if (fRunInfo->IsFitRangeInBin()) {
1766 fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1767 fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1768 // write these times back into the data structure. This way it is available when writting the log-file
1769 fRunInfo->SetFitRange(fFitStartTime, 0);
1770 fRunInfo->SetFitRange(fFitEndTime, 1);
1771 }
1772 if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
1773 fFitStartTime = globalBlock->GetFitRange(0);
1774 fFitEndTime = globalBlock->GetFitRange(1);
1775 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1776 if (globalBlock->IsFitRangeInBin()) {
1777 fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1778 fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1779 // write these times back into the data structure. This way it is available when writting the log-file
1780 globalBlock->SetFitRange(fFitStartTime, 0);
1781 globalBlock->SetFitRange(fFitEndTime, 1);
1782 }
1783 }
1785 fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
1786 fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
1787 std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
1788 std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
1789 }
1790}
std::vector< UInt_t > PUIntVector
Definition PMusr.h:361
#define ACCEL_PERIOD_TRIUMF
TRIUMF accelerator cycle: 43.37 ns.
Definition PMusr.h:151
EPMusrHandleTag
Definition PMusr.h:413
@ kFit
Fitting mode - perform least-squares fit to data.
Definition PMusr.h:415
@ kView
Viewing mode - display data and theory without fitting.
Definition PMusr.h:416
#define PMUSR_UNDEFINED
Definition PMusr.h:172
#define MSR_PARAM_FUN_OFFSET
Offset added to function indices for parameter parsing.
Definition PMusr.h:260
std::vector< PMsrParamStructure > PMsrParamList
Definition PMusr.h:1022
#define ACCEL_PERIOD_PSI
PSI (Paul Scherrer Institute) accelerator cycle: 19.75 ns.
Definition PMusr.h:149
std::vector< Double_t > PDoubleVector
Definition PMusr.h:385
#define ACCEL_PERIOD_RAL
RAL (Rutherford Appleton Lab) - pulsed beam.
Definition PMusr.h:153
if(xmlFile.is_open())
return status
virtual void SetFitRange(Double_t dval, UInt_t idx)
Definition PMusr.cpp:1137
virtual Double_t GetFitRange(UInt_t idx)
Definition PMusr.cpp:1120
virtual UInt_t GetT0BinSize()
Definition PMusr.h:1050
virtual Bool_t IsFitRangeInBin()
Definition PMusr.h:1055
virtual Int_t GetFitRangeOffset(UInt_t idx)
Definition PMusr.cpp:1157
virtual Double_t GetT0Bin(UInt_t idx=0)
Definition PMusr.cpp:1004
virtual Double_t GetRRFPhase()
Definition PMusr.h:1046
virtual Double_t GetRRFFreq(const char *unit)
Definition PMusr.cpp:831
MSR file parser and manager for the musrfit framework.
virtual PMsrParamList * GetMsrParamList()
Returns pointer to fit parameter list.
virtual const PDoubleVector * GetDataBin(const UInt_t histoNo)
Definition PMusr.h:880
virtual const Double_t GetT0Bin(const UInt_t histoNo)
Definition PMusr.h:870
virtual const Double_t GetTimeResolution()
Definition PMusr.h:868
virtual const Bool_t IsPresent(UInt_t histoNo)
Definition PMusr.h:869
virtual const UInt_t GetNoOfTemperatures()
Definition PMusr.h:860
virtual const Double_t GetField()
Definition PMusr.h:859
virtual const Double_t GetEnergy()
Definition PMusr.h:864
virtual const Double_t GetT0BinEstimated(const UInt_t histoNo)
Definition PMusr.h:871
virtual const PDoublePairVector * GetTemperature() const
Definition PMusr.h:861
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
Determines the proper fit range from global block.
PRunAsymmetryRRF()
Default constructor.
virtual UInt_t GetNoOfFitBins()
Returns the number of bins used in the fit.
Int_t fGoodBins[4]
Good bin boundaries: [0]=forward first, [1]=forward last, [2]=backward first, [3]=backward last.
PDoubleVector fForwardErr
Forward detector histogram errors.
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
Calculates expected chi-square (for statistical analysis).
virtual Bool_t GetProperDataRange(PRawRunData *runData, UInt_t histoNo[2])
Retrieves proper data range for histograms.
Bool_t fTheoAsData
If true, theory calculated only at data points; if false, extra points for nicer Fourier transforms.
virtual Bool_t PrepareViewData(PRawRunData *runData, UInt_t histoNo[2])
Prepares RRF data for viewing/plotting.
Bool_t SubtractEstimatedBkg()
Estimates and subtracts background from histograms.
virtual ~PRunAsymmetryRRF()
Destructor.
virtual Bool_t PrepareData()
Prepares all data for RRF fitting or viewing.
Int_t fEndTimeBin
Last bin index for fitting (after RRF transformation)
Int_t fRRFPacking
RRF packing factor from GLOBAL block (required for RRF analysis)
Bool_t SubtractFixBkg()
Subtracts fixed background from histograms.
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHisto, PUIntVector &backwardHistoNo)
Retrieves proper t0 values for all histograms.
virtual void SetFitRangeBin(const TString fitRange)
Sets the fit range in bins (can be changed dynamically via COMMAND block).
PDoubleVector fForward
Forward detector histogram data.
UInt_t fNoOfFitBins
Number of bins included in the fit after RRF packing.
PDoubleVector fBackward
Backward detector histogram data.
UInt_t fAlphaBetaTag
Tag indicating α/β configuration: 1=both unity, 2=α free/β unity, 3=α unity/β free,...
virtual Bool_t PrepareFitData()
Prepares RRF data specifically for fitting.
virtual void CalcTheory()
Calculates theoretical RRF asymmetry function.
virtual Double_t CalcChiSquare(const std::vector< Double_t > &par)
Calculates chi-square for the RRF asymmetry fit.
virtual void CalcNoOfFitBins()
Calculates the number of bins to be fitted.
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
Calculates maximum likelihood estimator.
Int_t fStartTimeBin
First bin index for fitting (after RRF transformation)
PDoubleVector fBackwardErr
Backward detector histogram errors.
Double_t fTimeResolution
Time resolution of raw histogram data in microseconds (μs), e.g., 0.01953125 μs for PSI GPS.
Definition PRunBase.h:276
Bool_t fValid
Flag indicating if run object initialized successfully; false if any error occurred.
Definition PRunBase.h:266
Double_t fFitEndTime
Fit range end time in microseconds (μs) relative to t0.
Definition PRunBase.h:282
PDoubleVector fFuncValues
Cached values of user-defined functions from FUNCTIONS block, evaluated at current parameters.
Definition PRunBase.h:284
PMsrHandler * fMsrInfo
Pointer to MSR file handler (owned externally, not deleted here)
Definition PRunBase.h:271
virtual void DeadTimeCorrection(std::vector< PDoubleVector > &histos, PUIntVector &histoNo)
carry out dead time correction
Definition PRunBase.cpp:171
PMetaData fMetaData
Experimental metadata extracted from data file header (magnetic field, temperature,...
Definition PRunBase.h:277
std::unique_ptr< PTheory > fTheory
Theory function evaluator (smart pointer, automatically deleted)
Definition PRunBase.h:285
std::vector< PDoubleVector > fAddT0s
Time-zero bin values for additional runs to be added to main run.
Definition PRunBase.h:279
EPMusrHandleTag fHandleTag
Operation mode: kFit (fitting), kView (display only), kEmpty (uninitialized)
Definition PRunBase.h:268
PRunData fData
Processed data container: background-corrected, packed, with theory values.
Definition PRunBase.h:275
PRunDataHandler * fRawData
Pointer to raw data handler (owned externally, not deleted here)
Definition PRunBase.h:273
PDoubleVector fT0s
Time-zero bin values for all histograms in this run (forward, backward, etc.)
Definition PRunBase.h:278
PRunBase()
Default constructor.
Definition PRunBase.cpp:54
Int_t fRunNo
Run number (0-based index in MSR file RUN blocks)
Definition PRunBase.h:270
PMsrRunBlock * fRunInfo
Pointer to this run's RUN block settings within fMsrInfo.
Definition PRunBase.h:272
Double_t fFitStartTime
Fit range start time in microseconds (μs) relative to t0.
Definition PRunBase.h:281
Raw data file reader and format converter for μSR data.