musrfit 1.10.0
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-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
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//--------------------------------------------------------------------------
60{
61 fNoOfFitBins = 0;
62 fPacking = -1;
63 fTheoAsData = false;
64
65 // the 2 following variables are need in case fit range is given in bins, and since
66 // the fit range can be changed in the command block, these variables need to be accessible
67 fGoodBins[0] = -1;
68 fGoodBins[1] = -1;
69
70 fStartTimeBin = -1;
71 fEndTimeBin = -1;
72}
73
74//--------------------------------------------------------------------------
75// Constructor
76//--------------------------------------------------------------------------
99PRunAsymmetry::PRunAsymmetry(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
100 PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
101{
102 // the following variables are need in case fit range is given in bins, and since
103 // the fit range can be changed in the command block, these variables need to be accessible
104 fGoodBins[0] = -1;
105 fGoodBins[1] = -1;
106 fGoodBins[2] = -1;
107 fGoodBins[3] = -1;
108
109 fStartTimeBin = -1;
110 fEndTimeBin = -1;
111
112 fPacking = fRunInfo->GetPacking();
113 if (fPacking == -1) { // i.e. packing is NOT given in the RUN-block, it must be given in the GLOBAL-block
114 fPacking = fMsrInfo->GetMsrGlobal()->GetPacking();
115 }
116 if (fPacking == -1) { // this should NOT happen, somethin is severely wrong
117 std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **SEVERE ERROR**: Couldn't find any packing information!";
118 std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
119 std::cerr << std::endl;
120 fValid = false;
121 return;
122 }
123
124 // check if alpha and/or beta is fixed --------------------
125
126 PMsrParamList *param = msrInfo->GetMsrParamList();
127
128 // check if alpha is given
129 if (fRunInfo->GetAlphaParamNo() == -1) { // no alpha given
130 std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** no alpha parameter given! This is needed for an asymmetry fit!";
131 std::cerr << std::endl;
132 fValid = false;
133 return;
134 }
135 // check if alpha parameter is within proper bounds
136 if (fRunInfo->GetAlphaParamNo() >= MSR_PARAM_FUN_OFFSET) { // alpha seems to be a function
137 if ((fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET < 0) ||
138 (fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET > msrInfo->GetNoOfFuncs())) {
139 std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** alpha parameter is a function with no = " << fRunInfo->GetAlphaParamNo();
140 std::cerr << std::endl << ">> This is out of bound, since there are only " << msrInfo->GetNoOfFuncs() << " functions.";
141 std::cerr << std::endl;
142 fValid = false;
143 return;
144 }
145 } else if ((fRunInfo->GetAlphaParamNo() < 0) || (fRunInfo->GetAlphaParamNo() > static_cast<Int_t>(param->size()))) {
146 std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** alpha parameter no = " << fRunInfo->GetAlphaParamNo();
147 std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
148 std::cerr << std::endl;
149 fValid = false;
150 return;
151 }
152 // check if alpha is fixed
153 Bool_t alphaFixedToOne = false;
154 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is parameter
155 if (((*param)[fRunInfo->GetAlphaParamNo()-1].fStep == 0.0) &&
156 ((*param)[fRunInfo->GetAlphaParamNo()-1].fValue == 1.0))
157 alphaFixedToOne = true;
158 }
159 // check if beta is given
160 Bool_t betaFixedToOne = false;
161 if (fRunInfo->GetBetaParamNo() == -1) { // no beta given hence assuming beta == 1
162 betaFixedToOne = true;
163 } else if ((fRunInfo->GetBetaParamNo() < 0) || (fRunInfo->GetBetaParamNo() > static_cast<Int_t>(param->size()))) { // check if beta parameter is within proper bounds
164 std::cerr << std::endl << ">> PRunAsymmetry::PRunAsymmetry(): **ERROR** beta parameter no = " << fRunInfo->GetBetaParamNo();
165 std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
166 std::cerr << std::endl;
167 fValid = false;
168 return;
169 } else { // check if beta is fixed
170 if (((*param)[fRunInfo->GetBetaParamNo()-1].fStep == 0.0) &&
171 ((*param)[fRunInfo->GetBetaParamNo()-1].fValue == 1.0))
172 betaFixedToOne = true;
173 }
174
175 // set fAlphaBetaTag
176 if (alphaFixedToOne && betaFixedToOne) // alpha == 1, beta == 1
177 fAlphaBetaTag = 1;
178 else if (!alphaFixedToOne && betaFixedToOne) // alpha != 1, beta == 1
179 fAlphaBetaTag = 2;
180 else if (alphaFixedToOne && !betaFixedToOne) // alpha == 1, beta != 1
181 fAlphaBetaTag = 3;
182 else
183 fAlphaBetaTag = 4;
184
185 // calculate fData
186 if (!PrepareData())
187 fValid = false;
188}
189
190//--------------------------------------------------------------------------
191// Destructor
192//--------------------------------------------------------------------------
200{
201 fForward.clear();
202 fForwardErr.clear();
203 fBackward.clear();
204 fBackwardErr.clear();
205}
206
207//--------------------------------------------------------------------------
208// CalcChiSquare (public)
209//--------------------------------------------------------------------------
227Double_t PRunAsymmetry::CalcChiSquare(const std::vector<Double_t>& par)
228{
229 Double_t chisq = 0.0;
230 Double_t diff = 0.0;
231 Double_t asymFcnValue = 0.0;
232 Double_t a, b, f;
233
234 // calculate functions
235 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
236 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
237 }
238
239 // calculate chi square
240 Double_t time(1.0);
241 Int_t i;
242
243 // determine alpha/beta
244 switch (fAlphaBetaTag) {
245 case 1: // alpha == 1, beta == 1
246 a = 1.0;
247 b = 1.0;
248 break;
249 case 2: // alpha != 1, beta == 1
250 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
251 a = par[fRunInfo->GetAlphaParamNo()-1];
252 } else { // alpha is function
253 // get function number
254 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
255 // evaluate function
256 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
257 }
258 b = 1.0;
259 break;
260 case 3: // alpha == 1, beta != 1
261 a = 1.0;
262 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
263 b = par[fRunInfo->GetBetaParamNo()-1];
264 } else { // beta is a function
265 // get function number
266 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
267 // evaluate function
268 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
269 }
270 break;
271 case 4: // alpha != 1, beta != 1
272 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
273 a = par[fRunInfo->GetAlphaParamNo()-1];
274 } else { // alpha is function
275 // get function number
276 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
277 // evaluate function
278 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
279 }
280 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
281 b = par[fRunInfo->GetBetaParamNo()-1];
282 } else { // beta is a function
283 // get function number
284 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
285 // evaluate function
286 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
287 }
288 break;
289 default:
290 a = 1.0;
291 b = 1.0;
292 break;
293 }
294
295 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
296 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
297 // for a given set of parameters---which should be done outside of the parallelized loop.
298 // For all other functions it means a tiny and acceptable overhead.
299 asymFcnValue = fTheory->Func(time, par, fFuncValues);
300
301 #ifdef HAVE_GOMP
302 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
303 if (chunk < 10)
304 chunk = 10;
305 #pragma omp parallel for default(shared) private(i,time,diff,asymFcnValue,f) schedule(dynamic,chunk) reduction(+:chisq)
306 #endif
307 for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
308 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
309 f = fTheory->Func(time, par, fFuncValues);
310 asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
311 diff = fData.GetValue()->at(i) - asymFcnValue;
312 chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
313 }
314
315 return chisq;
316}
317
318//--------------------------------------------------------------------------
319// CalcChiSquareExpected (public)
320//--------------------------------------------------------------------------
334Double_t PRunAsymmetry::CalcChiSquareExpected(const std::vector<Double_t>& par)
335{
336 return 0.0;
337}
338
339//--------------------------------------------------------------------------
340// CalcMaxLikelihood (public)
341//--------------------------------------------------------------------------
354Double_t PRunAsymmetry::CalcMaxLikelihood(const std::vector<Double_t>& par)
355{
356 std::cout << std::endl << "PRunAsymmetry::CalcMaxLikelihood(): not implemented yet ..." << std::endl;
357
358 return 1.0;
359}
360
361//--------------------------------------------------------------------------
362// GetNoOfFitBins (public)
363//--------------------------------------------------------------------------
374{
376
377 return fNoOfFitBins;
378}
379
380//--------------------------------------------------------------------------
381// SetFitRangeBin (public)
382//--------------------------------------------------------------------------
405void PRunAsymmetry::SetFitRangeBin(const TString fitRange)
406{
407 TObjArray *tok = nullptr;
408 TObjString *ostr = nullptr;
409 TString str;
410 Ssiz_t idx = -1;
411 Int_t offset = 0;
412
413 tok = fitRange.Tokenize(" \t");
414
415 if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
416 // handle fgb+n0 entry
417 ostr = dynamic_cast<TObjString*>(tok->At(1));
418 str = ostr->GetString();
419 // check if there is an offset present
420 idx = str.First("+");
421 if (idx != -1) { // offset present
422 str.Remove(0, idx+1);
423 if (str.IsFloat()) // if str is a valid number, convert is to an integer
424 offset = str.Atoi();
425 }
426 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
427
428 // handle lgb-n1 entry
429 ostr = dynamic_cast<TObjString*>(tok->At(2));
430 str = ostr->GetString();
431 // check if there is an offset present
432 idx = str.First("-");
433 if (idx != -1) { // offset present
434 str.Remove(0, idx+1);
435 if (str.IsFloat()) // if str is a valid number, convert is to an integer
436 offset = str.Atoi();
437 }
438 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
439 } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
440 Int_t pos = 2*(fRunNo+1)-1;
441
442 if (pos + 1 >= tok->GetEntries()) {
443 std::cerr << std::endl << ">> PRunAsymmetry::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
444 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
445 } else {
446 // handle fgb+n0 entry
447 ostr = static_cast<TObjString*>(tok->At(pos));
448 str = ostr->GetString();
449 // check if there is an offset present
450 idx = str.First("+");
451 if (idx != -1) { // offset present
452 str.Remove(0, idx+1);
453 if (str.IsFloat()) // if str is a valid number, convert is to an integer
454 offset = str.Atoi();
455 }
456 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
457
458 // handle lgb-n1 entry
459 ostr = static_cast<TObjString*>(tok->At(pos+1));
460 str = ostr->GetString();
461 // check if there is an offset present
462 idx = str.First("-");
463 if (idx != -1) { // offset present
464 str.Remove(0, idx+1);
465 if (str.IsFloat()) // if str is a valid number, convert is to an integer
466 offset = str.Atoi();
467 }
468 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
469 }
470 } else { // error
471 std::cerr << std::endl << ">> PRunAsymmetry::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
472 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
473 }
474
475 // clean up
476 if (tok) {
477 delete tok;
478 }
479}
480
481//--------------------------------------------------------------------------
482// CalcNoOfFitBins (public)
483//--------------------------------------------------------------------------
497{
498 // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
499 fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
500 if (fStartTimeBin < 0)
501 fStartTimeBin = 0;
502 fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
503 if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
504 fEndTimeBin = fData.GetValue()->size();
505
507 fNoOfFitBins = static_cast<UInt_t>(fEndTimeBin - fStartTimeBin);
508 else
509 fNoOfFitBins = 0;
510}
511
512//--------------------------------------------------------------------------
513// CalcTheory (protected)
514//--------------------------------------------------------------------------
532{
533 // feed the parameter vector
534 std::vector<Double_t> par;
535 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
536 for (UInt_t i=0; i<paramList->size(); i++)
537 par.push_back((*paramList)[i].fValue);
538
539 // calculate functions
540 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
541 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
542 }
543
544 // calculate asymmetry
545 Double_t asymFcnValue = 0.0;
546 Double_t a, b, f;
547 Double_t time;
548 for (UInt_t i=0; i<fData.GetValue()->size(); i++) {
549 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
550 switch (fAlphaBetaTag) {
551 case 1: // alpha == 1, beta == 1
552 asymFcnValue = fTheory->Func(time, par, fFuncValues);
553 break;
554 case 2: // alpha != 1, beta == 1
555 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
556 a = par[fRunInfo->GetAlphaParamNo()-1];
557 } else { // alpha is function
558 // get function number
559 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
560 // evaluate function
561 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
562 }
563 f = fTheory->Func(time, par, fFuncValues);
564 asymFcnValue = (f*(a+1.0)-(a-1.0))/((a+1.0)-f*(a-1.0));
565 break;
566 case 3: // alpha == 1, beta != 1
567 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
568 b = par[fRunInfo->GetBetaParamNo()-1];
569 } else { // beta is a function
570 // get function number
571 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
572 // evaluate function
573 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
574 }
575 f = fTheory->Func(time, par, fFuncValues);
576 asymFcnValue = f*(b+1.0)/(2.0-f*(b-1.0));
577 break;
578 case 4: // alpha != 1, beta != 1
579 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
580 a = par[fRunInfo->GetAlphaParamNo()-1];
581 } else { // alpha is function
582 // get function number
583 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
584 // evaluate function
585 a = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
586 }
587 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
588 b = par[fRunInfo->GetBetaParamNo()-1];
589 } else { // beta is a function
590 // get function number
591 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
592 // evaluate function
593 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
594 }
595 f = fTheory->Func(time, par, fFuncValues);
596 asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0));
597 break;
598 default:
599 asymFcnValue = 0.0;
600 break;
601 }
602 fData.AppendTheoryValue(asymFcnValue);
603 }
604
605 // clean up
606 par.clear();
607}
608
609//--------------------------------------------------------------------------
610// PrepareData (protected)
611//--------------------------------------------------------------------------
640{
641 if (!fValid)
642 return false;
643
644 // keep the Global block info
645 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
646
647 // get forward/backward histo from PRunDataHandler object ------------------------
648 // get the correct run
649 PRawRunData *runData = fRawData->GetRunData(*(fRunInfo->GetRunName()));
650 if (!runData) { // run not found
651 std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
652 std::cerr << std::endl;
653 return false;
654 }
655
656 // keep the field from the meta-data from the data-file
657 fMetaData.fField = runData->GetField();
658
659 // keep the energy from the meta-data from the data-file
660 fMetaData.fEnergy = runData->GetEnergy();
661
662 // keep the temperature(s) from the meta-data from the data-file
663 for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
664 fMetaData.fTemp.push_back(runData->GetTemperature(i));
665
666 // collect histogram numbers
667 PUIntVector forwardHistoNo;
668 PUIntVector backwardHistoNo;
669 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
670 forwardHistoNo.push_back(fRunInfo->GetForwardHistoNo(i));
671
672 if (!runData->IsPresent(forwardHistoNo[i])) {
673 std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **PANIC ERROR**:";
674 std::cerr << std::endl << ">> forwardHistoNo found = " << forwardHistoNo[i] << ", which is NOT present in the data file!?!?";
675 std::cerr << std::endl << ">> Will quit :-(";
676 std::cerr << std::endl;
677 // clean up
678 forwardHistoNo.clear();
679 backwardHistoNo.clear();
680 return false;
681 }
682 }
683 for (UInt_t i=0; i<fRunInfo->GetBackwardHistoNoSize(); i++) {
684 backwardHistoNo.push_back(fRunInfo->GetBackwardHistoNo(i));
685
686 if (!runData->IsPresent(backwardHistoNo[i])) {
687 std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **PANIC ERROR**:";
688 std::cerr << std::endl << ">> backwardHistoNo found = " << backwardHistoNo[i] << ", which is NOT present in the data file!?!?";
689 std::cerr << std::endl << ">> Will quit :-(";
690 std::cerr << std::endl;
691 // clean up
692 forwardHistoNo.clear();
693 backwardHistoNo.clear();
694 return false;
695 }
696 }
697
698 // keep the time resolution in (us)
699 fTimeResolution = runData->GetTimeResolution()/1.0e3;
700 std::cout.precision(10);
701 std::cout << std::endl << ">> PRunAsymmetry::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
702
703 // get all the proper t0's and addt0's for the current RUN block
704 if (!GetProperT0(runData, globalBlock, forwardHistoNo, backwardHistoNo)) {
705 return false;
706 }
707
708 // keep the histo of each group at this point (addruns handled below)
709 std::vector<PDoubleVector> forward, backward;
710 forward.resize(forwardHistoNo.size()); // resize to number of groups
711 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
712 forward[i].resize(runData->GetDataBin(forwardHistoNo[i])->size());
713 forward[i] = *runData->GetDataBin(forwardHistoNo[i]);
714 }
715 // check if a dead time correction has to be done
716 // this will be done automatically in the function itself, which also
717 // checks in the global and run section
718 DeadTimeCorrection(forward, forwardHistoNo);
719
720 backward.resize(backwardHistoNo.size()); // resize to number of groups
721 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
722 backward[i].resize(runData->GetDataBin(backwardHistoNo[i])->size());
723 backward[i] = *runData->GetDataBin(backwardHistoNo[i]);
724 }
725 // check if a dead time correction has to be done
726 // this will be done automatically in the function itself, which also
727 // checks in the global and run section
728 DeadTimeCorrection(backward, backwardHistoNo);
729
730 // check if addrun's are present, and if yes add data
731 // check if there are runs to be added to the current one
732 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
733 PRawRunData *addRunData;
734 std::vector<PDoubleVector> addForward, addBackward;
735 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
736 // get run to be added to the main one
737 addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
738 if (addRunData == nullptr) { // couldn't get run
739 std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
740 std::cerr << std::endl;
741 return false;
742 }
743
744 // dead time correction handling
745 addForward.clear();
746 addForward.resize(forwardHistoNo.size());
747 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
748 addForward[j].resize(addRunData->GetDataBin(forwardHistoNo[j])->size());
749 addForward[j] = *addRunData->GetDataBin(forwardHistoNo[j]);
750 }
751 DeadTimeCorrection(addForward, forwardHistoNo);
752 addBackward.clear();
753 addBackward.resize(backwardHistoNo.size());
754 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
755 addBackward[j].resize(addRunData->GetDataBin(backwardHistoNo[j])->size());
756 addBackward[j] = *addRunData->GetDataBin(backwardHistoNo[j]);
757 }
758 DeadTimeCorrection(addBackward, backwardHistoNo);
759
760 // add forward run
761 UInt_t addRunSize;
762 for (UInt_t k=0; k<forwardHistoNo.size(); k++) { // fill each group
763 addRunSize = addRunData->GetDataBin(forwardHistoNo[k])->size();
764 for (UInt_t j=0; j<addRunData->GetDataBin(forwardHistoNo[k])->size(); j++) { // loop over the bin indices
765 // make sure that the index stays in the proper range
766 if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]) >= 0) &&
767 (j+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]) < addRunSize)) {
768 forward[k][j] += addForward[k][j+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k])];
769 }
770 }
771 }
772
773 // add backward run
774 for (UInt_t k=0; k<backwardHistoNo.size(); k++) { // fill each group
775 addRunSize = addRunData->GetDataBin(backwardHistoNo[k])->size();
776 for (UInt_t j=0; j<addRunData->GetDataBin(backwardHistoNo[k])->size(); j++) { // loop over the bin indices
777 // make sure that the index stays in the proper range
778 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) &&
779 (j+static_cast<Int_t>(fAddT0s[i-1][2*k+1])-static_cast<Int_t>(fT0s[2*k+1]) < addRunSize)) {
780 backward[k][j] += addBackward[k][j+static_cast<Int_t>(fAddT0s[i-1][2*k+1])-static_cast<Int_t>(fT0s[2*k+1])];
781 }
782 }
783 }
784 }
785 }
786
787 // set forward/backward histo data of the first group
788 fForward.resize(forward[0].size());
789 for (UInt_t i=0; i<fForward.size(); i++) {
790 fForward[i] = forward[0][i];
791 }
792 fBackward.resize(backward[0].size());
793 for (UInt_t i=0; i<fBackward.size(); i++) {
794 fBackward[i] = backward[0][i];
795 }
796
797 // group histograms, add all the remaining forward histograms of the group
798 for (UInt_t i=1; i<forwardHistoNo.size(); i++) { // loop over the groupings
799 for (UInt_t j=0; j<runData->GetDataBin(forwardHistoNo[i])->size(); j++) { // loop over the bin indices
800 // make sure that the index stays within proper range
801 if ((j+fT0s[2*i]-fT0s[0] >= 0) && (j+fT0s[2*i]-fT0s[0] < runData->GetDataBin(forwardHistoNo[i])->size())) {
802 fForward[j] += forward[i][j+static_cast<Int_t>(fT0s[2*i])-static_cast<Int_t>(fT0s[0])];
803 }
804 }
805 }
806
807 // group histograms, add all the remaining backward histograms of the group
808 for (UInt_t i=1; i<backwardHistoNo.size(); i++) { // loop over the groupings
809 for (UInt_t j=0; j<runData->GetDataBin(backwardHistoNo[i])->size(); j++) { // loop over the bin indices
810 // make sure that the index stays within proper range
811 if ((j+fT0s[2*i+1]-fT0s[1] >= 0) && (j+fT0s[2*i+1]-fT0s[1] < runData->GetDataBin(backwardHistoNo[i])->size())) {
812 fBackward[j] += backward[i][j+static_cast<Int_t>(fT0s[2*i+1])-static_cast<Int_t>(fT0s[1])];
813 }
814 }
815 }
816
817 // subtract background from histogramms ------------------------------------------
818 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
819 if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
821 return false;
822 } else { // no background given to do the job, try to estimate it
823 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
824 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
825 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.1), 2);
826 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.6), 3);
827 std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **WARNING** Neither fix background nor background bins are given!";
828 std::cerr << std::endl << ">> Will try the following:";
829 std::cerr << std::endl << ">> forward: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
830 std::cerr << std::endl << ">> backward: bkg start = " << fRunInfo->GetBkgRange(2) << ", bkg end = " << fRunInfo->GetBkgRange(3);
831 std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
832 std::cerr << std::endl;
834 return false;
835 }
836 } else { // fixed background given
837 if (!SubtractFixBkg())
838 return false;
839 }
840
841 UInt_t histoNo[2] = {forwardHistoNo[0], backwardHistoNo[0]};
842
843 // get the data range (fgb/lgb) for the current RUN block
844 if (!GetProperDataRange(runData, histoNo)) {
845 return false;
846 }
847
848 // get the fit range for the current RUN block
849 GetProperFitRange(globalBlock);
850
851 // everything looks fine, hence fill data set
852 Bool_t status;
853 switch(fHandleTag) {
854 case kFit:
856 break;
857 case kView:
858 if (fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking == 0)
859 status = PrepareViewData(runData, histoNo);
860 else
861 status = PrepareRRFViewData(runData, histoNo);
862 break;
863 default:
864 status = false;
865 break;
866 }
867
868 // clean up
869 forwardHistoNo.clear();
870 backwardHistoNo.clear();
871
872 return status;
873}
874
875//--------------------------------------------------------------------------
876// SubtractFixBkg (private)
877//--------------------------------------------------------------------------
895{
896 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) {
897 std::cerr << "PRunAsymmetry::SubtractFixBkg(): **ERROR** no fixed bkg for forward set. Will do nothing here." << std::endl;
898 return false;
899 }
900 if (fRunInfo->GetBkgFix(1) == PMUSR_UNDEFINED) {
901 std::cerr << "PRunAsymmetry::SubtractFixBkg(): **ERROR** no fixed bkg for backward set. Will do nothing here." << std::endl;
902 std::cerr << " you need an entry like:" << std::endl;
903 std::cerr << " backgr.fix 2 3" << std::endl;
904 std::cerr << " i.e. two entries for forward and backward." << std::endl;
905 return false;
906 }
907
908
909 Double_t dval;
910
911 for (UInt_t i=0; i<fForward.size(); i++) {
912 // keep the error, and make sure that the bin is NOT empty
913 if (fForward[i] != 0.0)
914 dval = TMath::Sqrt(fForward[i]);
915 else
916 dval = 1.0;
917 fForwardErr.push_back(dval);
918 fForward[i] -= fRunInfo->GetBkgFix(0);
919
920 // keep the error, and make sure that the bin is NOT empty
921 if (fBackward[i] != 0.0)
922 dval = TMath::Sqrt(fBackward[i]);
923 else
924 dval = 1.0;
925 fBackwardErr.push_back(dval);
926 fBackward[i] -= fRunInfo->GetBkgFix(1);
927 }
928
929 return true;
930}
931
932//--------------------------------------------------------------------------
933// SubtractEstimatedBkg (private)
934//--------------------------------------------------------------------------
953{
954 Double_t beamPeriod = 0.0;
955
956 // check if data are from PSI, RAL, or TRIUMF
957 if (fRunInfo->GetInstitute()->Contains("psi"))
958 beamPeriod = ACCEL_PERIOD_PSI;
959 else if (fRunInfo->GetInstitute()->Contains("ral"))
960 beamPeriod = ACCEL_PERIOD_RAL;
961 else if (fRunInfo->GetInstitute()->Contains("triumf"))
962 beamPeriod = ACCEL_PERIOD_TRIUMF;
963 else
964 beamPeriod = 0.0;
965
966 // check if start and end are in proper order
967 Int_t start[2] = {fRunInfo->GetBkgRange(0), fRunInfo->GetBkgRange(2)};
968 Int_t end[2] = {fRunInfo->GetBkgRange(1), fRunInfo->GetBkgRange(3)};
969 for (UInt_t i=0; i<2; i++) {
970 if (end[i] < start[i]) {
971 std::cout << std::endl << ">> PRunAsymmetry::SubtractEstimatedBkg(): end = " << end[i] << " > start = " << start[i] << "! Will swap them!";
972 UInt_t keep = end[i];
973 end[i] = start[i];
974 start[i] = keep;
975 }
976 }
977
978 // calculate proper background range
979 for (UInt_t i=0; i<2; i++) {
980 if (beamPeriod != 0.0) {
981 Double_t timeBkg = static_cast<Double_t>(end[i]-start[i])*(fTimeResolution*fPacking); // length of the background intervall in time
982 UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
983 // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
984 end[i] = start[i] + static_cast<UInt_t>((fullCycles*beamPeriod)/(fTimeResolution*fPacking));
985 std::cout << ">> PRunAsymmetry::SubtractEstimatedBkg(): Background " << start[i] << ", " << end[i] << std::endl;
986 if (end[i] == start[i])
987 end[i] = fRunInfo->GetBkgRange(2*i+1);
988 }
989 }
990
991 // check if start is within histogram bounds
992 if ((start[0] < 0) || (start[0] >= fForward.size()) ||
993 (start[1] < 0) || (start[1] >= fBackward.size())) {
994 std::cerr << std::endl << ">> PRunAsymmetry::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
995 std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
996 std::cerr << std::endl << ">> background start (f/b) = (" << start[0] << "/" << start[1] << ").";
997 return false;
998 }
999
1000 // check if end is within histogram bounds
1001 if ((end[0] < 0) || (end[0] >= fForward.size()) ||
1002 (end[1] < 0) || (end[1] >= fBackward.size())) {
1003 std::cerr << std::endl << ">> PRunAsymmetry::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
1004 std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForward.size() << "/" << fBackward.size() << ").";
1005 std::cerr << std::endl << ">> background end (f/b) = (" << end[0] << "/" << end[1] << ").";
1006 return false;
1007 }
1008
1009 // calculate background
1010 Double_t bkg[2] = {0.0, 0.0};
1011 Double_t errBkg[2] = {0.0, 0.0};
1012
1013 // forward
1014 for (UInt_t i=start[0]; i<=end[0]; i++)
1015 bkg[0] += fForward[i];
1016 errBkg[0] = TMath::Sqrt(bkg[0])/(end[0] - start[0] + 1);
1017 bkg[0] /= static_cast<Double_t>(end[0] - start[0] + 1);
1018 std::cout << std::endl << ">> estimated forward histo background: " << bkg[0];
1019
1020 // backward
1021 for (UInt_t i=start[1]; i<=end[1]; i++)
1022 bkg[1] += fBackward[i];
1023 errBkg[1] = TMath::Sqrt(bkg[1])/(end[1] - start[1] + 1);
1024 bkg[1] /= static_cast<Double_t>(end[1] - start[1] + 1);
1025 std::cout << std::endl << ">> estimated backward histo background: " << bkg[1] << std::endl;
1026
1027 // correct error for forward, backward
1028 Double_t errVal = 0.0;
1029 for (UInt_t i=0; i<fForward.size(); i++) {
1030 if (fForward[i] > 0.0)
1031 errVal = TMath::Sqrt(fForward[i]+errBkg[0]*errBkg[0]);
1032 else
1033 errVal = 1.0;
1034 fForwardErr.push_back(errVal);
1035 if (fBackward[i] > 0.0)
1036 errVal = TMath::Sqrt(fBackward[i]+errBkg[1]*errBkg[1]);
1037 else
1038 errVal = 1.0;
1039 fBackwardErr.push_back(errVal);
1040 }
1041
1042 // subtract background from data
1043 for (UInt_t i=0; i<fForward.size(); i++) {
1044 fForward[i] -= bkg[0];
1045 fBackward[i] -= bkg[1];
1046 }
1047
1048 fRunInfo->SetBkgEstimated(bkg[0], 0);
1049 fRunInfo->SetBkgEstimated(bkg[1], 1);
1050
1051 return true;
1052}
1053
1054//--------------------------------------------------------------------------
1055// PrepareFitData (protected)
1056//--------------------------------------------------------------------------
1078{
1079 // transform raw histo data. This is done the following way (for details see the manual):
1080 // first rebin the data, than calculate the asymmetry
1081
1082 // everything looks fine, hence fill packed forward and backward histo
1083 PRunData forwardPacked;
1084 PRunData backwardPacked;
1085 Double_t value = 0.0;
1086 Double_t error = 0.0;
1087 // forward
1088 for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
1089 if (fPacking == 1) {
1090 forwardPacked.AppendValue(fForward[i]);
1091 forwardPacked.AppendErrorValue(fForwardErr[i]);
1092 } else { // packed data, i.e. fPacking > 1
1093 if (((i-fGoodBins[0]) % fPacking == 0) && (i != fGoodBins[0])) { // fill data
1094 // in order that after rebinning the fit does not need to be redone (important for plots)
1095 // the value is normalize to per bin
1096 value /= fPacking;
1097 forwardPacked.AppendValue(value);
1098 if (value == 0.0)
1099 forwardPacked.AppendErrorValue(1.0);
1100 else
1101 forwardPacked.AppendErrorValue(TMath::Sqrt(error)/fPacking);
1102 value = 0.0;
1103 error = 0.0;
1104 }
1105 value += fForward[i];
1106 error += fForwardErr[i]*fForwardErr[i];
1107 }
1108 }
1109 // backward
1110 for (Int_t i=fGoodBins[2]; i<fGoodBins[3]; i++) {
1111 if (fPacking == 1) {
1112 backwardPacked.AppendValue(fBackward[i]);
1113 backwardPacked.AppendErrorValue(fBackwardErr[i]);
1114 } else { // packed data, i.e. fPacking > 1
1115 if (((i-fGoodBins[2]) % fPacking == 0) && (i != fGoodBins[2])) { // fill data
1116 // in order that after rebinning the fit does not need to be redone (important for plots)
1117 // the value is normalize to per bin
1118 value /= fPacking;
1119 backwardPacked.AppendValue(value);
1120 if (value == 0.0)
1121 backwardPacked.AppendErrorValue(1.0);
1122 else
1123 backwardPacked.AppendErrorValue(TMath::Sqrt(error)/fPacking);
1124 value = 0.0;
1125 error = 0.0;
1126 }
1127 value += fBackward[i];
1128 error += fBackwardErr[i]*fBackwardErr[i];
1129 }
1130 }
1131
1132 // check if packed forward and backward hist have the same size, otherwise take the minimum size
1133 UInt_t noOfBins = forwardPacked.GetValue()->size();
1134 if (forwardPacked.GetValue()->size() != backwardPacked.GetValue()->size()) {
1135 if (forwardPacked.GetValue()->size() > backwardPacked.GetValue()->size())
1136 noOfBins = backwardPacked.GetValue()->size();
1137 }
1138
1139 // form asymmetry including error propagation
1140 Double_t asym;
1141 Double_t f, b, ef, eb;
1142 // fill data time start, and step
1143 // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1144 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(fGoodBins[0])-0.5) + static_cast<Double_t>(fPacking)/2.0 - static_cast<Double_t>(fT0s[0])));
1145 fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(fPacking));
1146 for (UInt_t i=0; i<noOfBins; i++) {
1147 // to make the formulae more readable
1148 f = forwardPacked.GetValue()->at(i);
1149 b = backwardPacked.GetValue()->at(i);
1150 ef = forwardPacked.GetError()->at(i);
1151 eb = backwardPacked.GetError()->at(i);
1152 // check that there are indeed bins
1153 if (f+b != 0.0)
1154 asym = (f-b) / (f+b);
1155 else
1156 asym = 0.0;
1157 fData.AppendValue(asym);
1158 // calculate the error
1159 if (f+b != 0.0)
1160 error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1161 else
1162 error = 1.0;
1163 fData.AppendErrorValue(error);
1164 }
1165
1167
1168 // clean up
1169 fForward.clear();
1170 fForwardErr.clear();
1171 fBackward.clear();
1172 fBackwardErr.clear();
1173
1174 return true;
1175}
1176
1177//--------------------------------------------------------------------------
1178// PrepareViewData (protected)
1179//--------------------------------------------------------------------------
1201Bool_t PRunAsymmetry::PrepareViewData(PRawRunData* runData, UInt_t histoNo[2])
1202{
1203 // check if view_packing is wished
1204 Int_t packing = fPacking;
1205 if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
1206 packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1207 }
1208
1209 // feed the parameter vector
1210 std::vector<Double_t> par;
1211 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1212 for (UInt_t i=0; i<paramList->size(); i++)
1213 par.push_back((*paramList)[i].fValue);
1214
1215 // transform raw histo data. This is done the following way (for details see the manual):
1216 // first rebin the data, than calculate the asymmetry
1217
1218 // first get start data, end data, and t0
1219 Int_t start[2] = {fGoodBins[0], fGoodBins[2]};
1220 Int_t end[2] = {fGoodBins[1], fGoodBins[3]};
1221 Int_t t0[2] = {static_cast<Int_t>(fT0s[0]), static_cast<Int_t>(fT0s[1])};
1222
1223 // check if the data ranges and t0's between forward/backward are compatible
1224 Int_t fgb[2];
1225 if (start[0]-t0[0] != start[1]-t0[1]) { // wrong fgb aligning
1226 if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
1227 fgb[0] = start[0];
1228 fgb[1] = t0[1] + start[0]-t0[0];
1229 std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **WARNING** needed to shift backward fgb from ";
1230 std::cerr << start[1] << " to " << fgb[1] << std::endl;
1231 } else {
1232 fgb[0] = t0[0] + start[1]-t0[1];
1233 fgb[1] = start[1];
1234 std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **WARNING** needed to shift forward fgb from ";
1235 std::cerr << start[0] << " to " << fgb[0] << std::endl;
1236 }
1237 } else { // fgb aligning is correct
1238 fgb[0] = start[0];
1239 fgb[1] = start[1];
1240 }
1241
1242 Int_t val = fgb[0]-packing*(fgb[0]/packing);
1243 do {
1244 if (fgb[1] - fgb[0] < 0)
1245 val += packing;
1246 } while (val + fgb[1] - fgb[0] < 0);
1247
1248 start[0] = val;
1249 start[1] = val + fgb[1] - fgb[0];
1250
1251 // make sure that there are equal number of rebinned bins in forward and backward
1252 UInt_t noOfBins0 = (runData->GetDataBin(histoNo[0])->size()-start[0])/packing;
1253 UInt_t noOfBins1 = (runData->GetDataBin(histoNo[1])->size()-start[1])/packing;
1254 if (noOfBins0 > noOfBins1)
1255 noOfBins0 = noOfBins1;
1256 end[0] = start[0] + noOfBins0 * packing;
1257 end[1] = start[1] + noOfBins0 * packing;
1258
1259 // check if start, end, and t0 make any sense
1260 // 1st check if start and end are in proper order
1261 for (UInt_t i=0; i<2; i++) {
1262 if (end[i] < start[i]) { // need to swap them
1263 Int_t keep = end[i];
1264 end[i] = start[i];
1265 start[i] = keep;
1266 }
1267 // 2nd check if start is within proper bounds
1268 if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1269 std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
1270 std::cerr << std::endl;
1271 return false;
1272 }
1273 // 3rd check if end is within proper bounds
1274 if ((end[i] < 0) || (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1275 std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
1276 std::cerr << std::endl;
1277 return false;
1278 }
1279 // 4th check if t0 is within proper bounds
1280 if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1281 std::cerr << std::endl << ">> PRunAsymmetry::PrepareViewData(): **ERROR** t0 data bin doesn't make any sense!";
1282 std::cerr << std::endl;
1283 return false;
1284 }
1285 }
1286
1287 // everything looks fine, hence fill packed forward and backward histo
1288 PRunData forwardPacked;
1289 PRunData backwardPacked;
1290 Double_t value = 0.0;
1291 Double_t error = 0.0;
1292
1293 // forward
1294 for (Int_t i=start[0]; i<end[0]; i++) {
1295 if (packing == 1) {
1296 forwardPacked.AppendValue(fForward[i]);
1297 forwardPacked.AppendErrorValue(fForwardErr[i]);
1298 } else { // packed data, i.e. packing > 1
1299 if (((i-start[0]) % packing == 0) && (i != start[0])) { // fill data
1300 // in order that after rebinning the fit does not need to be redone (important for plots)
1301 // the value is normalize to per bin
1302 value /= packing;
1303 forwardPacked.AppendValue(value);
1304 if (value == 0.0)
1305 forwardPacked.AppendErrorValue(1.0);
1306 else
1307 forwardPacked.AppendErrorValue(TMath::Sqrt(error)/packing);
1308 value = 0.0;
1309 error = 0.0;
1310 }
1311 value += fForward[i];
1312 error += fForwardErr[i]*fForwardErr[i];
1313 }
1314 }
1315
1316 // backward
1317 for (Int_t i=start[1]; i<end[1]; i++) {
1318 if (packing == 1) {
1319 backwardPacked.AppendValue(fBackward[i]);
1320 backwardPacked.AppendErrorValue(fBackwardErr[i]);
1321 } else { // packed data, i.e. packing > 1
1322 if (((i-start[1]) % packing == 0) && (i != start[1])) { // fill data
1323 // in order that after rebinning the fit does not need to be redone (important for plots)
1324 // the value is normalize to per bin
1325 value /= packing;
1326 backwardPacked.AppendValue(value);
1327 if (value == 0.0)
1328 backwardPacked.AppendErrorValue(1.0);
1329 else
1330 backwardPacked.AppendErrorValue(TMath::Sqrt(error)/packing);
1331 value = 0.0;
1332 error = 0.0;
1333 }
1334 value += fBackward[i];
1335 error += fBackwardErr[i]*fBackwardErr[i];
1336 }
1337 }
1338
1339 // check if packed forward and backward hist have the same size, otherwise take the minimum size
1340 UInt_t noOfBins = forwardPacked.GetValue()->size();
1341 if (forwardPacked.GetValue()->size() != backwardPacked.GetValue()->size()) {
1342 if (forwardPacked.GetValue()->size() > backwardPacked.GetValue()->size())
1343 noOfBins = backwardPacked.GetValue()->size();
1344 }
1345
1346 // form asymmetry including error propagation
1347 Double_t asym;
1348 Double_t f, b, ef, eb, alpha = 1.0, beta = 1.0;
1349 // set data time start, and step
1350 // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1351 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start[0])-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0[0])));
1352 fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(packing));
1353
1354 // get the proper alpha and beta
1355 switch (fAlphaBetaTag) {
1356 case 1: // alpha == 1, beta == 1
1357 alpha = 1.0;
1358 beta = 1.0;
1359 break;
1360 case 2: // alpha != 1, beta == 1
1361 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1362 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1363 } else { // alpha is function
1364 // get function number
1365 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1366 // evaluate function
1367 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1368 }
1369 beta = 1.0;
1370 break;
1371 case 3: // alpha == 1, beta != 1
1372 alpha = 1.0;
1373 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1374 beta = par[fRunInfo->GetBetaParamNo()-1];
1375 } else { // beta is a function
1376 // get function number
1377 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1378 // evaluate function
1379 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1380 }
1381 break;
1382 case 4: // alpha != 1, beta != 1
1383 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1384 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1385 } else { // alpha is function
1386 // get function number
1387 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1388 // evaluate function
1389 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1390 }
1391 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1392 beta = par[fRunInfo->GetBetaParamNo()-1];
1393 } else { // beta is a function
1394 // get function number
1395 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1396 // evaluate function
1397 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1398 }
1399 break;
1400 default:
1401 break;
1402 }
1403
1404 for (UInt_t i=0; i<forwardPacked.GetValue()->size(); i++) {
1405 // to make the formulae more readable
1406 f = forwardPacked.GetValue()->at(i);
1407 b = backwardPacked.GetValue()->at(i);
1408 ef = forwardPacked.GetError()->at(i);
1409 eb = backwardPacked.GetError()->at(i);
1410 // check that there are indeed bins
1411 if (f+b != 0.0)
1412 asym = (alpha*f-b) / (alpha*beta*f+b);
1413 else
1414 asym = 0.0;
1415 fData.AppendValue(asym);
1416 // calculate the error
1417 if (f+b != 0.0)
1418 error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1419 else
1420 error = 1.0;
1421 fData.AppendErrorValue(error);
1422 }
1423
1425
1426 // clean up
1427 fForward.clear();
1428 fForwardErr.clear();
1429 fBackward.clear();
1430 fBackwardErr.clear();
1431
1432 // fill theory vector for kView
1433 // calculate functions
1434 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1435 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
1436 }
1437
1438 // calculate theory
1439 Double_t time;
1440 UInt_t size = runData->GetDataBin(histoNo[0])->size()/packing;
1441 Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1442
1443 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1444 if (fTheoAsData) { // calculate theory only at the data points
1445 fData.SetTheoryTimeStep(fData.GetDataTimeStep());
1446 } else {
1447 // finer binning for the theory (8 times as many points = factor)
1448 size *= factor;
1449 fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1450 }
1451
1452 for (UInt_t i=0; i<size; i++) {
1453 time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1454 value = fTheory->Func(time, par, fFuncValues);
1455 if (fabs(value) > 10.0) { // dirty hack needs to be fixed!!
1456 value = 0.0;
1457 }
1458 fData.AppendTheoryValue(value);
1459 }
1460
1461 // clean up
1462 par.clear();
1463
1464 return true;
1465}
1466
1467//--------------------------------------------------------------------------
1468// PrepareRRFViewData (protected)
1469//--------------------------------------------------------------------------
1493Bool_t PRunAsymmetry::PrepareRRFViewData(PRawRunData* runData, UInt_t histoNo[2])
1494{
1495 // feed the parameter vector
1496 std::vector<Double_t> par;
1497 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1498 for (UInt_t i=0; i<paramList->size(); i++)
1499 par.push_back((*paramList)[i].fValue);
1500
1501 // ------------------------------------------------------------
1502 // 1. make all necessary checks
1503 // ------------------------------------------------------------
1504
1505 // first get start data, end data, and t0
1506 Int_t start[2] = {fGoodBins[0], fGoodBins[2]};
1507 Int_t end[2] = {fGoodBins[1], fGoodBins[3]};
1508 Int_t t0[2] = {static_cast<Int_t>(fT0s[0]), static_cast<Int_t>(fT0s[1])};
1509 UInt_t packing = fMsrInfo->GetMsrPlotList()->at(0).fRRFPacking;
1510
1511 // check if the data ranges and t0's between forward/backward are compatible
1512 Int_t fgb[2];
1513 if (start[0]-t0[0] != start[1]-t0[1]) { // wrong fgb aligning
1514 if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
1515 fgb[0] = start[0];
1516 fgb[1] = t0[1] + start[0]-t0[0];
1517 std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **WARNING** needed to shift backward fgb from ";
1518 std::cerr << start[1] << " to " << fgb[1] << std::endl;
1519 } else {
1520 fgb[0] = t0[0] + start[1]-t0[1];
1521 fgb[1] = start[1];
1522 std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **WARNING** needed to shift forward fgb from ";
1523 std::cerr << start[1] << " to " << fgb[0] << std::endl;
1524 }
1525 } else { // fgb aligning is correct
1526 fgb[0] = start[0];
1527 fgb[1] = start[1];
1528 }
1529
1530 Int_t val = fgb[0]-packing*(fgb[0]/packing);
1531 do {
1532 if (fgb[1] - fgb[0] < 0)
1533 val += packing;
1534 } while (val + fgb[1] - fgb[0] < 0);
1535
1536 start[0] = val;
1537 start[1] = val + fgb[1] - fgb[0];
1538
1539 // make sure that there are equal number of rebinned bins in forward and backward
1540 UInt_t noOfBins0 = runData->GetDataBin(histoNo[0])->size()-start[0];
1541 UInt_t noOfBins1 = runData->GetDataBin(histoNo[1])->size()-start[1];
1542 if (noOfBins0 > noOfBins1)
1543 noOfBins0 = noOfBins1;
1544 end[0] = start[0] + noOfBins0;
1545 end[1] = start[1] + noOfBins0;
1546
1547 // check if start, end, and t0 make any sense
1548 // 1st check if start and end are in proper order
1549 for (UInt_t i=0; i<2; i++) {
1550 if (end[i] < start[i]) { // need to swap them
1551 Int_t keep = end[i];
1552 end[i] = start[i];
1553 start[i] = keep;
1554 }
1555 // 2nd check if start is within proper bounds
1556 if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1557 std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **ERROR** start data bin doesn't make any sense!";
1558 std::cerr << std::endl;
1559 return false;
1560 }
1561 // 3rd check if end is within proper bounds
1562 if ((end[i] < 0) || (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1563 std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **ERROR** end data bin doesn't make any sense!";
1564 std::cerr << std::endl;
1565 return false;
1566 }
1567 // 4th check if t0 is within proper bounds
1568 if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1569 std::cerr << std::endl << ">> PRunAsymmetry::PrepareRRFViewData(): **ERROR** t0 data bin doesn't make any sense!";
1570 std::cerr << std::endl;
1571 return false;
1572 }
1573 }
1574
1575 // ------------------------------------------------------------
1576 // 2. build the asymmetry [A(t)] WITHOUT packing.
1577 // ------------------------------------------------------------
1578
1579 PDoubleVector forward, forwardErr;
1580 PDoubleVector backward, backwardErr;
1581 Double_t error = 0.0;
1582 // forward
1583 for (Int_t i=start[0]; i<end[0]; i++) {
1584 forward.push_back(fForward[i]);
1585 forwardErr.push_back(fForwardErr[i]);
1586 }
1587 // backward
1588 for (Int_t i=start[1]; i<end[1]; i++) {
1589 backward.push_back(fBackward[i]);
1590 backwardErr.push_back(fBackwardErr[i]);
1591 }
1592
1593 // check if packed forward and backward hist have the same size, otherwise take the minimum size
1594 UInt_t noOfBins = forward.size();
1595 if (forward.size() != backward.size()) {
1596 if (forward.size() > backward.size())
1597 noOfBins = backward.size();
1598 }
1599
1600 // form asymmetry including error propagation
1601 PDoubleVector asymmetry, asymmetryErr;
1602 Double_t asym;
1603 Double_t f, b, ef, eb, alpha = 1.0, beta = 1.0;
1604
1605 // get the proper alpha and beta
1606 switch (fAlphaBetaTag) {
1607 case 1: // alpha == 1, beta == 1
1608 alpha = 1.0;
1609 beta = 1.0;
1610 break;
1611 case 2: // alpha != 1, beta == 1
1612 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1613 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1614 } else { // alpha is function
1615 // get function number
1616 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1617 // evaluate function
1618 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1619 }
1620 beta = 1.0;
1621 break;
1622 case 3: // alpha == 1, beta != 1
1623 alpha = 1.0;
1624 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1625 beta = par[fRunInfo->GetBetaParamNo()-1];
1626 } else { // beta is a function
1627 // get function number
1628 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1629 // evaluate function
1630 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1631 }
1632 break;
1633 case 4: // alpha != 1, beta != 1
1634 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1635 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1636 } else { // alpha is function
1637 // get function number
1638 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1639 // evaluate function
1640 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1641 }
1642 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1643 beta = par[fRunInfo->GetBetaParamNo()-1];
1644 } else { // beta is a function
1645 // get function number
1646 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1647 // evaluate function
1648 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1649 }
1650 break;
1651 default:
1652 break;
1653 }
1654
1655 for (UInt_t i=0; i<noOfBins; i++) {
1656 // to make the formulae more readable
1657 f = forward[i];
1658 b = backward[i];
1659 ef = forwardErr[i];
1660 eb = backwardErr[i];
1661 // check that there are indeed bins
1662 if (f+b != 0.0)
1663 asym = (alpha*f-b) / (alpha*beta*f+b);
1664 else
1665 asym = 0.0;
1666 asymmetry.push_back(asym);
1667 // calculate the error
1668 if (f+b != 0.0)
1669 error = 2.0/((f+b)*(f+b))*TMath::Sqrt(b*b*ef*ef+eb*eb*f*f);
1670 else
1671 error = 1.0;
1672 asymmetryErr.push_back(error);
1673 }
1674
1675 // clean up
1676 fForward.clear();
1677 fForwardErr.clear();
1678 fBackward.clear();
1679 fBackwardErr.clear();
1680
1681
1682 // ------------------------------------------------------------
1683 // 3. A_R(t) = A(t) * 2 cos(w_R t + phi_R)
1684 // ------------------------------------------------------------
1685
1686 // check which units shall be used
1687 Double_t gammaRRF = 0.0, wRRF = 0.0, phaseRRF = 0.0;
1688 Double_t time;
1689
1690 switch (fMsrInfo->GetMsrPlotList()->at(0).fRRFUnit) {
1691 case RRF_UNIT_kHz:
1692 gammaRRF = TMath::TwoPi()*1.0e-3;
1693 break;
1694 case RRF_UNIT_MHz:
1695 gammaRRF = TMath::TwoPi();
1696 break;
1697 case RRF_UNIT_Mcs:
1698 gammaRRF = 1.0;
1699 break;
1700 case RRF_UNIT_G:
1701 gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi();
1702 break;
1703 case RRF_UNIT_T:
1704 gammaRRF = GAMMA_BAR_MUON*TMath::TwoPi()*1.0e4;
1705 break;
1706 default:
1707 gammaRRF = TMath::TwoPi();
1708 break;
1709 }
1710 wRRF = gammaRRF * fMsrInfo->GetMsrPlotList()->at(0).fRRFFreq;
1711 phaseRRF = fMsrInfo->GetMsrPlotList()->at(0).fRRFPhase / 180.0 * TMath::Pi();
1712
1713 for (UInt_t i=0; i<asymmetry.size(); i++) {
1714 time = fTimeResolution*(static_cast<Double_t>(start[0])-t0[0]+static_cast<Double_t>(i));
1715 asymmetry[i] *= 2.0*TMath::Cos(wRRF*time+phaseRRF);
1716 }
1717
1718 // ------------------------------------------------------------
1719 // 4. do the packing of A_R(t)
1720 // ------------------------------------------------------------
1721 Double_t value = 0.0;
1722 error = 0.0;
1723 for (UInt_t i=0; i<asymmetry.size(); i++) {
1724 if ((i % packing == 0) && (i != 0)) {
1725 value /= packing;
1726 fData.AppendValue(value);
1727 fData.AppendErrorValue(TMath::Sqrt(error)/packing);
1728
1729 value = 0.0;
1730 error = 0.0;
1731 }
1732 value += asymmetry[i];
1733 error += asymmetryErr[i]*asymmetryErr[i];
1734 }
1735
1736 // set data time start, and step
1737 // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
1738 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start[0])-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0[0])));
1739 fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(packing));
1740
1741 // ------------------------------------------------------------
1742 // 5. calculate theory [T(t)] as close as possible to the time resolution [compatible with the RRF frequency]
1743 // 6. T_R(t) = T(t) * 2 cos(w_R t + phi_R)
1744 // ------------------------------------------------------------
1745 UInt_t rebinRRF = static_cast<UInt_t>((TMath::Pi()/2.0/wRRF)/fTimeResolution); // RRF time resolution / data time resolution
1746 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1747 fData.SetTheoryTimeStep(TMath::Pi()/2.0/wRRF/rebinRRF); // = theory time resolution as close as possible to the data time resolution compatible with wRRF
1748
1749 // calculate functions
1750 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1751 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
1752 }
1753
1754 Double_t theoryValue;
1755 for (UInt_t i=0; i<asymmetry.size(); i++) {
1756 time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1757 theoryValue = fTheory->Func(time, par, fFuncValues);
1758 theoryValue *= 2.0*TMath::Cos(wRRF * time + phaseRRF);
1759
1760 if (fabs(theoryValue) > 10.0) { // dirty hack needs to be fixed!!
1761 theoryValue = 0.0;
1762 }
1763
1764 fData.AppendTheoryValue(theoryValue);
1765 }
1766
1767 // ------------------------------------------------------------
1768 // 7. do the packing of T_R(t)
1769 // ------------------------------------------------------------
1770
1771 PDoubleVector theo;
1772 Double_t dval = 0.0;
1773 for (UInt_t i=0; i<fData.GetTheory()->size(); i++) {
1774 if ((i % rebinRRF == 0) && (i != 0)) {
1775 theo.push_back(dval/rebinRRF);
1776 dval = 0.0;
1777 }
1778 dval += fData.GetTheory()->at(i);
1779 }
1780 fData.ReplaceTheory(theo);
1781 theo.clear();
1782
1783 // set the theory time start and step size
1784 fData.SetTheoryTimeStart(fData.GetTheoryTimeStart()+static_cast<Double_t>(rebinRRF-1)*fData.GetTheoryTimeStep()/2.0);
1785 fData.SetTheoryTimeStep(rebinRRF*fData.GetTheoryTimeStep());
1786
1787 // ------------------------------------------------------------
1788 // 8. calculate the Kaiser FIR filter coefficients
1789 // ------------------------------------------------------------
1790 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)
1791
1792 // ------------------------------------------------------------
1793 // 9. filter T_R(t)
1794 // ------------------------------------------------------------
1795 FilterTheo();
1796
1797 // clean up
1798 par.clear();
1799 forward.clear();
1800 forwardErr.clear();
1801 backward.clear();
1802 backwardErr.clear();
1803 asymmetry.clear();
1804 asymmetryErr.clear();
1805
1806 return true;
1807}
1808
1809//--------------------------------------------------------------------------
1810// GetProperT0 (private)
1811//--------------------------------------------------------------------------
1840Bool_t PRunAsymmetry::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHistoNo, PUIntVector &backwardHistoNo)
1841{
1842 // feed all T0's
1843 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1844 fT0s.clear();
1845 // this strange fT0 size estimate is needed in case #forw histos != #back histos
1846 size_t size = 2*forwardHistoNo.size();
1847 if (backwardHistoNo.size() > forwardHistoNo.size())
1848 size = 2*backwardHistoNo.size();
1849 fT0s.resize(size);
1850 for (UInt_t i=0; i<fT0s.size(); i++) {
1851 fT0s[i] = -1.0;
1852 }
1853
1854 // fill in the T0's from the msr-file (if present)
1855 for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1856 fT0s[i] = fRunInfo->GetT0Bin(i);
1857 }
1858
1859 // fill in the missing T0's from the GLOBAL block section (if present)
1860 for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1861 if (fT0s[i] == -1) { // i.e. not given in the RUN block section
1862 fT0s[i] = globalBlock->GetT0Bin(i);
1863 }
1864 }
1865
1866 // fill in the missing T0's from the data file, if not already present in the msr-file
1867 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1868 if (fT0s[2*i] == -1.0) // i.e. not present in the msr-file, try the data file
1869 if (runData->GetT0Bin(forwardHistoNo[i]) > 0.0) {
1870 fT0s[2*i] = runData->GetT0Bin(forwardHistoNo[i]);
1871 fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1872 }
1873 }
1874 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1875 if (fT0s[2*i+1] == -1.0) // i.e. not present in the msr-file, try the data file
1876 if (runData->GetT0Bin(backwardHistoNo[i]) > 0.0) {
1877 fT0s[2*i+1] = runData->GetT0Bin(backwardHistoNo[i]);
1878 fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1879 }
1880 }
1881
1882 // 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
1883 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1884 if (fT0s[2*i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1885 fT0s[2*i] = runData->GetT0BinEstimated(forwardHistoNo[i]);
1886 fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1887
1888 std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1889 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1890 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(forwardHistoNo[i]);
1891 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1892 std::cerr << std::endl;
1893 }
1894 }
1895 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1896 if (fT0s[2*i+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1897 fT0s[2*i+1] = runData->GetT0BinEstimated(backwardHistoNo[i]);
1898 fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1899
1900 std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1901 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1902 std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[i]);
1903 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1904 std::cerr << std::endl;
1905 }
1906 }
1907
1908 // check if t0 is within proper bounds
1909 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1910 if ((fT0s[2*i] < 0) || (fT0s[2*i] > static_cast<Int_t>(runData->GetDataBin(forwardHistoNo[i])->size()))) {
1911 std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **ERROR** t0 data bin (" << fT0s[2*i] << ") doesn't make any sense!";
1912 std::cerr << std::endl << ">> forwardHistoNo " << forwardHistoNo[i];
1913 std::cerr << std::endl;
1914 return false;
1915 }
1916 }
1917 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1918 if ((fT0s[2*i+1] < 0) || (fT0s[2*i+1] > static_cast<Int_t>(runData->GetDataBin(backwardHistoNo[i])->size()))) {
1919 std::cerr << std::endl << ">> PRunAsymmetry::PrepareData(): **ERROR** t0 data bin (" << fT0s[2*i+1] << ") doesn't make any sense!";
1920 std::cerr << std::endl << ">> backwardHistoNo " << backwardHistoNo[i];
1921 std::cerr << std::endl;
1922 return false;
1923 }
1924 }
1925
1926 // check if addrun's are present, and if yes add the necessary t0's
1927 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1928 PRawRunData *addRunData;
1929 fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1930 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1931 // get run to be added to the main one
1932 addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
1933 if (addRunData == 0) { // couldn't get run
1934 std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1935 std::cerr << std::endl;
1936 return false;
1937 }
1938
1939 // feed all T0's
1940 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1941 fAddT0s[i-1].clear();
1942 fAddT0s[i-1].resize(2*forwardHistoNo.size());
1943 for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1944 fAddT0s[i-1][j] = -1.0;
1945 }
1946
1947 // fill in the T0's from the msr-file (if present)
1948 for (Int_t j=0; j<fRunInfo->GetAddT0BinSize(i-1); j++) {
1949 fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1, j);
1950 }
1951
1952 // fill in the T0's from the data file, if not already present in the msr-file
1953 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1954 if (fAddT0s[i-1][2*j] == -1.0) // i.e. not present in the msr-file, try the data file
1955 if (addRunData->GetT0Bin(forwardHistoNo[j]) > 0.0) {
1956 fAddT0s[i-1][2*j] = addRunData->GetT0Bin(forwardHistoNo[j]);
1957 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1958 }
1959 }
1960 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1961 if (fAddT0s[i-1][2*j+1] == -1.0) // i.e. not present in the msr-file, try the data file
1962 if (addRunData->GetT0Bin(backwardHistoNo[j]) > 0.0) {
1963 fAddT0s[i-1][2*j+1] = addRunData->GetT0Bin(backwardHistoNo[j]);
1964 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1965 }
1966 }
1967
1968 // 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
1969 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1970 if (fAddT0s[i-1][2*j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1971 fAddT0s[i-1][2*j] = addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1972 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1973
1974 std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1975 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1976 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1977 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1978 std::cerr << std::endl;
1979 }
1980 }
1981 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1982 if (fAddT0s[i-1][2*j+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1983 fAddT0s[i-1][2*j+1] = addRunData->GetT0BinEstimated(backwardHistoNo[j]);
1984 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1985
1986 std::cerr << std::endl << ">> PRunAsymmetry::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1987 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1988 std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[j]);
1989 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1990 std::cerr << std::endl;
1991 }
1992 }
1993 }
1994 }
1995
1996 return true;
1997}
1998
1999//--------------------------------------------------------------------------
2000// GetProperDataRange (private)
2001//--------------------------------------------------------------------------
2028Bool_t PRunAsymmetry::GetProperDataRange(PRawRunData* runData, UInt_t histoNo[2])
2029{
2030 // first get start/end data
2031 Int_t start[2] = {fRunInfo->GetDataRange(0), fRunInfo->GetDataRange(2)};
2032 Int_t end[2] = {fRunInfo->GetDataRange(1), fRunInfo->GetDataRange(3)};
2033 // check if data range has been provided in the RUN block. If not, try the GLOBAL block
2034 if (start[0] == -1) {
2035 start[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
2036 }
2037 if (start[1] == -1) {
2038 start[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(2);
2039 }
2040 if (end[0] == -1) {
2041 end[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
2042 }
2043 if (end[1] == -1) {
2044 end[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(3);
2045 }
2046
2047 Double_t t0[2] = {fT0s[0], fT0s[1]};
2048 Int_t offset = static_cast<Int_t>((10.0e-3/fTimeResolution)); // needed in case first good bin is not given, default = 10ns
2049
2050 // check if data range has been provided, and if not try to estimate them
2051 if (start[0] < 0) {
2052 start[0] = static_cast<Int_t>(t0[0])+offset;
2053 fRunInfo->SetDataRange(start[0], 0);
2054 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[0] << ".";
2055 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
2056 std::cerr << std::endl;
2057 }
2058 if (start[1] < 0) {
2059 start[1] = static_cast<Int_t>(t0[1])+offset;
2060 fRunInfo->SetDataRange(start[1], 2);
2061 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[1] << ".";
2062 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
2063 std::cerr << std::endl;
2064 }
2065 if (end[0] < 0) {
2066 end[0] = runData->GetDataBin(histoNo[0])->size();
2067 fRunInfo->SetDataRange(end[0], 1);
2068 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range end = " << end[0] << ".";
2069 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
2070 std::cerr << std::endl;
2071 }
2072 if (end[1] < 0) {
2073 end[1] = runData->GetDataBin(histoNo[1])->size();
2074 fRunInfo->SetDataRange(end[1], 3);
2075 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range end = " << end[1] << ".";
2076 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
2077 std::cerr << std::endl;
2078 }
2079
2080 // check if start, end, and t0 make any sense
2081 // 1st check if start and end are in proper order
2082 for (UInt_t i=0; i<2; i++) {
2083 if (end[i] < start[i]) { // need to swap them
2084 Int_t keep = end[i];
2085 end[i] = start[i];
2086 start[i] = keep;
2087 }
2088 // 2nd check if start is within proper bounds
2089 if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
2090 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **ERROR** start data bin doesn't make any sense!";
2091 std::cerr << std::endl;
2092 return false;
2093 }
2094 // 3rd check if end is within proper bounds
2095 if (end[i] < 0) {
2096 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **ERROR** end data bin (" << end[i] << ") doesn't make any sense!";
2097 std::cerr << std::endl;
2098 return false;
2099 }
2100 if (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())) {
2101 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **WARNING** end data bin (" << end[i] << ") > histo length (" << static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()) << ").";
2102 std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
2103 std::cerr << std::endl;
2104 end[i] = static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())-1;
2105 }
2106 // 4th check if t0 is within proper bounds
2107 if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
2108 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange(): **ERROR** t0 data bin doesn't make any sense!";
2109 std::cerr << std::endl;
2110 return false;
2111 }
2112 }
2113
2114 // check that start-t0 is the same for forward as for backward, otherwise take max(start[i]-t0[i])
2115 if (fabs(static_cast<Double_t>(start[0])-t0[0]) > fabs(static_cast<Double_t>(start[1])-t0[1])){
2116 start[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(start[0]) - t0[0]);
2117 end[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(end[0]) - t0[0]);
2118 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift backward data range.";
2119 std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(2) << ", " << fRunInfo->GetDataRange(3);
2120 std::cerr << std::endl << ">> used : " << start[1] << ", " << end[1];
2121 std::cerr << std::endl;
2122 }
2123 if (fabs(static_cast<Double_t>(start[0])-t0[0]) < fabs(static_cast<Double_t>(start[1])-t0[1])){
2124 start[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(start[1]) - t0[1]);
2125 end[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(end[1]) - t0[1]);
2126 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift forward data range.";
2127 std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(0) << ", " << fRunInfo->GetDataRange(1);
2128 std::cerr << std::endl << ">> used : " << start[0] << ", " << end[0];
2129 std::cerr << std::endl;
2130 }
2131
2132 // keep good bins for potential latter use
2133 fGoodBins[0] = start[0];
2134 fGoodBins[1] = end[0];
2135 fGoodBins[2] = start[1];
2136 fGoodBins[3] = end[1];
2137
2138 // make sure that fGoodBins are in proper range for fForward and fBackward
2139 if (fGoodBins[0] < 0)
2140 fGoodBins[0]=0;
2141 if (fGoodBins[1] > fForward.size()) {
2142 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift forward lgb,";
2143 std::cerr << std::endl << ">> from " << fGoodBins[1] << " to " << fForward.size()-1 << std::endl;
2144 fGoodBins[1]=fForward.size()-1;
2145 }
2146 if (fGoodBins[2] < 0)
2147 fGoodBins[2]=0;
2148 if (fGoodBins[3] > fBackward.size()) {
2149 std::cerr << std::endl << ">> PRunAsymmetry::GetProperDataRange **WARNING** needed to shift backward lgb,";
2150 std::cerr << std::endl << ">> from " << fGoodBins[1] << " to " << fForward.size()-1 << std::endl;
2151 fGoodBins[3]=fBackward.size()-1;
2152 }
2153
2154 return true;
2155}
2156
2157//--------------------------------------------------------------------------
2158// GetProperFitRange (private)
2159//--------------------------------------------------------------------------
2183{
2184 // set fit start/end time; first check RUN Block
2185 fFitStartTime = fRunInfo->GetFitRange(0);
2186 fFitEndTime = fRunInfo->GetFitRange(1);
2187 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
2188 if (fRunInfo->IsFitRangeInBin()) {
2189 fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
2190 fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
2191 // write these times back into the data structure. This way it is available when writting the log-file
2192 fRunInfo->SetFitRange(fFitStartTime, 0);
2193 fRunInfo->SetFitRange(fFitEndTime, 1);
2194 }
2195 if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
2196 fFitStartTime = globalBlock->GetFitRange(0);
2197 fFitEndTime = globalBlock->GetFitRange(1);
2198 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
2199 if (globalBlock->IsFitRangeInBin()) {
2200 fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
2201 fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
2202 // write these times back into the data structure. This way it is available when writting the log-file
2203 globalBlock->SetFitRange(fFitStartTime, 0);
2204 globalBlock->SetFitRange(fFitEndTime, 1);
2205 }
2206 }
2208 fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
2209 fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
2210 std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
2211 std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
2212 }
2213}
std::vector< UInt_t > PUIntVector
Definition PMusr.h:361
#define ACCEL_PERIOD_TRIUMF
TRIUMF accelerator cycle: 43.37 ns.
Definition PMusr.h:151
EPMusrHandleTag
Definition PMusr.h:413
@ kFit
Fitting mode - perform least-squares fit to data.
Definition PMusr.h:415
@ kView
Viewing mode - display data and theory without fitting.
Definition PMusr.h:416
#define RRF_UNIT_MHz
Frequency in MHz (megahertz)
Definition PMusr.h:335
#define PMUSR_UNDEFINED
Definition PMusr.h:172
#define GAMMA_BAR_MUON
Definition PMusr.h:138
#define 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 RRF_UNIT_Mcs
Angular frequency in Mc/s (Mega-cycles per second)
Definition PMusr.h:337
#define RRF_UNIT_G
Equivalent magnetic field in Gauss (G)
Definition PMusr.h:339
#define RRF_UNIT_kHz
Frequency in kHz (kilohertz)
Definition PMusr.h:333
#define ACCEL_PERIOD_PSI
PSI (Paul Scherrer Institute) accelerator cycle: 19.75 ns.
Definition PMusr.h:149
#define RRF_UNIT_T
Equivalent magnetic field in Tesla (T)
Definition PMusr.h:341
std::vector< Double_t > PDoubleVector
Definition PMusr.h:385
#define ACCEL_PERIOD_RAL
RAL (Rutherford Appleton Lab) - pulsed beam.
Definition PMusr.h:153
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
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
UInt_t fNoOfFitBins
Number of bins included in the fit.
PRunAsymmetry()
Default constructor.
Int_t fPacking
Bin packing factor from RUN or GLOBAL block.
PDoubleVector fForwardErr
Forward detector histogram errors.
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
Determines the proper fit range from global block.
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHisto, PUIntVector &backwardHistoNo)
Retrieves proper t0 values for all histograms.
Bool_t SubtractFixBkg()
Subtracts fixed background from histograms.
PDoubleVector fBackwardErr
Backward detector histogram errors.
virtual void CalcNoOfFitBins()
Calculates the number of bins to be fitted.
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
Calculates expected chi-square (for statistical analysis).
virtual void CalcTheory()
Calculates theoretical asymmetry function.
virtual Bool_t PrepareViewData(PRawRunData *runData, UInt_t histoNo[2])
Prepares data for viewing/plotting.
UInt_t fAlphaBetaTag
Tag indicating α/β configuration: 1=both unity, 2=α free/β unity, 3=α unity/β free,...
virtual Bool_t PrepareFitData()
Prepares data specifically for fitting.
Int_t fGoodBins[4]
Good bin boundaries: [0]=forward first, [1]=forward last, [2]=backward first, [3]=backward last.
Int_t fEndTimeBin
Last bin index for fitting.
virtual Double_t CalcChiSquare(const std::vector< Double_t > &par)
Calculates chi-square for the current parameter set.
virtual ~PRunAsymmetry()
Destructor.
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
Calculates maximum likelihood estimator.
virtual Bool_t PrepareRRFViewData(PRawRunData *runData, UInt_t histoNo[2])
Prepares rotating reference frame (RRF) data for viewing.
virtual UInt_t GetNoOfFitBins()
Returns the number of bins used in the fit.
virtual Bool_t GetProperDataRange(PRawRunData *runData, UInt_t histoNo[2])
Retrieves proper data range for histograms.
virtual Bool_t PrepareData()
Prepares all data for fitting or viewing.
PDoubleVector fForward
Forward detector histogram data.
Bool_t SubtractEstimatedBkg()
Estimates and subtracts background from histograms.
PDoubleVector fBackward
Backward detector histogram data.
Int_t fStartTimeBin
First bin index for fitting.
Bool_t fTheoAsData
If true, theory calculated only at data points; if false, extra points for nicer Fourier transforms.
virtual void SetFitRangeBin(const TString fitRange)
Sets the fit range in bins.
Double_t fTimeResolution
Time resolution of raw histogram data in microseconds (μs), e.g., 0.01953125 μs for PSI GPS.
Definition PRunBase.h:276
Bool_t fValid
Flag indicating if run object initialized successfully; false if any error occurred.
Definition PRunBase.h:266
Double_t fFitEndTime
Fit range end time in microseconds (μs) relative to t0.
Definition PRunBase.h:282
PDoubleVector fFuncValues
Cached values of user-defined functions from FUNCTIONS block, evaluated at current parameters.
Definition PRunBase.h:284
PMsrHandler * fMsrInfo
Pointer to MSR file handler (owned externally, not deleted here)
Definition PRunBase.h:271
virtual void DeadTimeCorrection(std::vector< PDoubleVector > &histos, PUIntVector &histoNo)
carry out dead time correction
Definition PRunBase.cpp:171
PMetaData fMetaData
Experimental metadata extracted from data file header (magnetic field, temperature,...
Definition PRunBase.h:277
std::unique_ptr< PTheory > fTheory
Theory function evaluator (smart pointer, automatically deleted)
Definition PRunBase.h:285
std::vector< PDoubleVector > fAddT0s
Time-zero bin values for additional runs to be added to main run.
Definition PRunBase.h:279
EPMusrHandleTag fHandleTag
Operation mode: kFit (fitting), kView (display only), kEmpty (uninitialized)
Definition PRunBase.h:268
virtual void CalculateKaiserFilterCoeff(Double_t wc, Double_t A, Double_t dw)
Calculates Kaiser window FIR filter coefficients for RRF smoothing.
Definition PRunBase.cpp:339
PRunData fData
Processed data container: background-corrected, packed, with theory values.
Definition PRunBase.h:275
PRunDataHandler * fRawData
Pointer to raw data handler (owned externally, not deleted here)
Definition PRunBase.h:273
PDoubleVector fT0s
Time-zero bin values for all histograms in this run (forward, backward, etc.)
Definition PRunBase.h:278
PRunBase()
Default constructor.
Definition PRunBase.cpp:54
Int_t fRunNo
Run number (0-based index in MSR file RUN blocks)
Definition PRunBase.h:270
PMsrRunBlock * fRunInfo
Pointer to this run's RUN block settings within fMsrInfo.
Definition PRunBase.h:272
Double_t fFitStartTime
Fit range start time in microseconds (μs) relative to t0.
Definition PRunBase.h:281
virtual void FilterTheo()
Applies Kaiser FIR filter to theory values for RRF fits.
Definition PRunBase.cpp:405
Raw data file reader and format converter for μSR data.
virtual void AppendErrorValue(Double_t dval)
Definition PMusr.h:497
virtual void AppendValue(Double_t dval)
Definition PMusr.h:494
virtual const PDoubleVector * GetValue()
Returns pointer to data value vector (asymmetry, counts, or y-data)
Definition PMusr.h:468
virtual const PDoubleVector * GetError()
Returns pointer to data error vector (statistical uncertainties)
Definition PMusr.h:470