musrfit 1.10.0
PRunAsymmetryBNMR.cpp
Go to the documentation of this file.
1/***************************************************************************
2
3 PRunAsymmetryBNMR.cpp
4
5 Author: Zaher Salman
6 Based on PRunAsymmetry.cpp by Andreas Suter
7 e-mail: zaher.salman@psi.ch
8
9***************************************************************************/
10
11/***************************************************************************
12 * Copyright (C) 2018-2026 by Zaher Salman *
13 * zaher.salman@psi.ch *
14 * *
15 * This program is free software; you can redistribute it and/or modify *
16 * it under the terms of the GNU General Public License as published by *
17 * the Free Software Foundation; either version 2 of the License, or *
18 * (at your option) any later version. *
19 * *
20 * This program is distributed in the hope that it will be useful, *
21 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
22 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
23 * GNU General Public License for more details. *
24 * *
25 * You should have received a copy of the GNU General Public License *
26 * along with this program; if not, write to the *
27 * Free Software Foundation, Inc., *
28 * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
29 ***************************************************************************/
30
31#ifdef HAVE_CONFIG_H
32#include "config.h"
33#endif
34
35#ifdef HAVE_GOMP
36#include <omp.h>
37#endif
38
39#include <stdio.h>
40
41#include <iostream>
42
43#include <TString.h>
44#include <TObjArray.h>
45#include <TObjString.h>
46
47#include "PMusr.h"
48#include "PRunAsymmetryBNMR.h"
49
50//--------------------------------------------------------------------------
51// Constructor
52//--------------------------------------------------------------------------
61{
62 fNoOfFitBins = 0;
63 fPacking = -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 fStartTimeBin = -1;
72 fEndTimeBin = -1;
73}
74
75//--------------------------------------------------------------------------
76// Constructor
77//--------------------------------------------------------------------------
102PRunAsymmetryBNMR::PRunAsymmetryBNMR(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
103 PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
104{
105 // the 2 following variables are need in case fit range is given in bins, and since
106 // the fit range can be changed in the command block, these variables need to be accessible
107 fGoodBins[0] = -1;
108 fGoodBins[1] = -1;
109
110 fStartTimeBin = -1;
111 fEndTimeBin = -1;
112
113 fPacking = fRunInfo->GetPacking();
114 if (fPacking == -1) { // i.e. packing is NOT given in the RUN-block, it must be given in the GLOBAL-block
115 fPacking = fMsrInfo->GetMsrGlobal()->GetPacking();
116 }
117 if (fPacking == -1) { // this should NOT happen, somethin is severely wrong
118 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PRunAsymmetryBNMR(): **SEVERE ERROR**: Couldn't find any packing information!";
119 std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
120 std::cerr << std::endl;
121 fValid = false;
122 return;
123 }
124
125 // check if alpha and/or beta is fixed --------------------
126
127 PMsrParamList *param = msrInfo->GetMsrParamList();
128
129 // check if alpha is given
130 Bool_t alphaFixedToOne = false;
131 Bool_t alphaSetDefault = false;
132 if (fRunInfo->GetAlphaParamNo() == -1) { // no alpha given
133 // std::cerr << std::endl << ">> PRunAsymmetryBNMR::PRunAsymmetryBNMR(): **ERROR** no alpha parameter given! This is needed for an asymmetry fit!";
134 // std::cerr << std::endl;
135 // fValid = false;
136 // return;
137 alphaSetDefault = true;
138 } else if ((fRunInfo->GetAlphaParamNo() < 0) || (fRunInfo->GetAlphaParamNo() > (Int_t)param->size())) { // check if alpha parameter is within proper bounds
139 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PRunAsymmetryBNMR(): **ERROR** alpha parameter no = " << fRunInfo->GetAlphaParamNo();
140 std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
141 std::cerr << std::endl;
142 fValid = false;
143 return;
144 } else { // check if alpha is fixed
145 if (((*param)[fRunInfo->GetAlphaParamNo()-1].fStep == 0.0) &&
146 ((*param)[fRunInfo->GetAlphaParamNo()-1].fValue == 1.0))
147 alphaFixedToOne = true;
148 }
149
150 // check if beta is given
151 Bool_t betaFixedToOne = false;
152 if (fRunInfo->GetBetaParamNo() == -1) { // no beta given hence assuming beta == 1
153 betaFixedToOne = true;
154 } else if ((fRunInfo->GetBetaParamNo() < 0) || (fRunInfo->GetBetaParamNo() > (Int_t)param->size())) { // check if beta parameter is within proper bounds
155 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PRunAsymmetryBNMR(): **ERROR** beta parameter no = " << fRunInfo->GetBetaParamNo();
156 std::cerr << std::endl << ">> This is out of bound, since there are only " << param->size() << " parameters.";
157 std::cerr << std::endl;
158 fValid = false;
159 return;
160 } else { // check if beta is fixed
161 if (((*param)[fRunInfo->GetBetaParamNo()-1].fStep == 0.0) &&
162 ((*param)[fRunInfo->GetBetaParamNo()-1].fValue == 1.0))
163 betaFixedToOne = true;
164 }
165
166 // set fAlphaBetaTag
167 if (alphaFixedToOne && betaFixedToOne) // alpha == 1, beta == 1
168 fAlphaBetaTag = 1;
169 else if (!alphaFixedToOne && betaFixedToOne & !alphaSetDefault) // alpha != 1, beta == 1
170 fAlphaBetaTag = 2;
171 else if (alphaFixedToOne && !betaFixedToOne) // alpha == 1, beta != 1
172 fAlphaBetaTag = 3;
173 else if (!alphaFixedToOne && betaFixedToOne & alphaSetDefault) // alpha ??, beta == 1
174 fAlphaBetaTag = 5;
175 else if (!alphaFixedToOne && !betaFixedToOne & alphaSetDefault) // alpha ??, beta != 1
176 fAlphaBetaTag = 6;
177 else
178 fAlphaBetaTag = 4;
179
180 // calculate fData
181 if (!PrepareData())
182 fValid = false;
183}
184
185//--------------------------------------------------------------------------
186// Destructor
187//--------------------------------------------------------------------------
195{
196 fForwardp.clear();
197 fForwardpErr.clear();
198 fBackwardp.clear();
199 fBackwardpErr.clear();
200
201 fForwardm.clear();
202 fForwardmErr.clear();
203 fBackwardm.clear();
204 fBackwardmErr.clear();
205}
206
207//--------------------------------------------------------------------------
208// CalcChiSquare (public)
209//--------------------------------------------------------------------------
226Double_t PRunAsymmetryBNMR::CalcChiSquare(const std::vector<Double_t>& par)
227{
228 Double_t chisq = 0.0;
229 Double_t diff = 0.0;
230 Double_t asymFcnValue = 0.0;
231 Double_t a, b, f;
232
233 // calculate functions
234 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
235 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
236 }
237
238 // calculate chi square
239 Double_t time(1.0),alphaest;
240 Int_t i;
241
242 // determine alpha/beta
243 alphaest = fRunInfo->GetEstimatedAlpha();
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 case 5: // alpha ?? , beta == 1
290 a = alphaest;
291 b = 1.0;
292 break;
293 case 6: // alpha ??, beta != 1
294 a = alphaest;
295 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
296 b = par[fRunInfo->GetBetaParamNo()-1];
297 } else { // beta is a function
298 // get function number
299 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
300 // evaluate function
301 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
302 }
303 break;
304 default:
305 a = 1.0;
306 b = 1.0;
307 break;
308 }
309
310 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
311 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
312 // for a given set of parameters---which should be done outside of the parallelized loop.
313 // For all other functions it means a tiny and acceptable overhead.
314 asymFcnValue = fTheory->Func(time, par, fFuncValues);
315
316 #ifdef HAVE_GOMP
317 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
318 if (chunk < 10)
319 chunk = 10;
320 #pragma omp parallel for default(shared) private(i,time,diff,asymFcnValue,f) schedule(dynamic,chunk) reduction(+:chisq)
321 #endif
322 for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
323 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
324 f = fTheory->Func(time, par, fFuncValues)/2.0;
325 asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0))-(-f*(a*b+1.0)-(a-1.0))/((a+1.0)+f*(a*b-1.0));
326 diff = fData.GetValue()->at(i) - asymFcnValue;
327 chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
328 }
329
330 return chisq;
331}
332
333//--------------------------------------------------------------------------
334// CalcChiSquareExpected (public)
335//--------------------------------------------------------------------------
345Double_t PRunAsymmetryBNMR::CalcChiSquareExpected(const std::vector<Double_t>& par)
346{
347 return 0.0;
348}
349
350//--------------------------------------------------------------------------
351// CalcMaxLikelihood (public)
352//--------------------------------------------------------------------------
362Double_t PRunAsymmetryBNMR::CalcMaxLikelihood(const std::vector<Double_t>& par)
363{
364 std::cout << std::endl << "PRunAsymmetryBNMR::CalcMaxLikelihood(): not implemented yet ..." << std::endl;
365
366 return 1.0;
367}
368
369//--------------------------------------------------------------------------
370// GetNoOfFitBins (public)
371//--------------------------------------------------------------------------
380{
382
383 return fNoOfFitBins;
384}
385
386//--------------------------------------------------------------------------
387// SetFitRangeBin (public)
388//--------------------------------------------------------------------------
407void PRunAsymmetryBNMR::SetFitRangeBin(const TString fitRange)
408{
409 TObjArray *tok = nullptr;
410 TObjString *ostr = nullptr;
411 TString str;
412 Ssiz_t idx = -1;
413 Int_t offset = 0;
414
415 tok = fitRange.Tokenize(" \t");
416
417 if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
418 // handle fgb+n0 entry
419 ostr = dynamic_cast<TObjString*>(tok->At(1));
420 str = ostr->GetString();
421 // check if there is an offset present
422 idx = str.First("+");
423 if (idx != -1) { // offset present
424 str.Remove(0, idx+1);
425 if (str.IsFloat()) // if str is a valid number, convert is to an integer
426 offset = str.Atoi();
427 }
428 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
429
430 // handle lgb-n1 entry
431 ostr = dynamic_cast<TObjString*>(tok->At(2));
432 str = ostr->GetString();
433 // check if there is an offset present
434 idx = str.First("-");
435 if (idx != -1) { // offset present
436 str.Remove(0, idx+1);
437 if (str.IsFloat()) // if str is a valid number, convert is to an integer
438 offset = str.Atoi();
439 }
440 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
441 } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
442 Int_t pos = 2*(fRunNo+1)-1;
443
444 if (pos + 1 >= tok->GetEntries()) {
445 std::cerr << std::endl << ">> PRunAsymmetryBNMR::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
446 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
447 } else {
448 // handle fgb+n0 entry
449 ostr = static_cast<TObjString*>(tok->At(pos));
450 str = ostr->GetString();
451 // check if there is an offset present
452 idx = str.First("+");
453 if (idx != -1) { // offset present
454 str.Remove(0, idx+1);
455 if (str.IsFloat()) // if str is a valid number, convert is to an integer
456 offset = str.Atoi();
457 }
458 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
459
460 // handle lgb-n1 entry
461 ostr = static_cast<TObjString*>(tok->At(pos+1));
462 str = ostr->GetString();
463 // check if there is an offset present
464 idx = str.First("-");
465 if (idx != -1) { // offset present
466 str.Remove(0, idx+1);
467 if (str.IsFloat()) // if str is a valid number, convert is to an integer
468 offset = str.Atoi();
469 }
470 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
471 }
472 } else { // error
473 std::cerr << std::endl << ">> PRunAsymmetryBNMR::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
474 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
475 }
476
477 // clean up
478 if (tok) {
479 delete tok;
480 }
481}
482
483//--------------------------------------------------------------------------
484// CalcNoOfFitBins (public)
485//--------------------------------------------------------------------------
494{
495 // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
496 fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
497 if (fStartTimeBin < 0)
498 fStartTimeBin = 0;
499 fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
500 if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
501 fEndTimeBin = fData.GetValue()->size();
502
504 fNoOfFitBins = static_cast<UInt_t>(fEndTimeBin - fStartTimeBin);
505 else
506 fNoOfFitBins = 0;
507}
508
509//--------------------------------------------------------------------------
510// CalcTheory (protected)
511//--------------------------------------------------------------------------
528{
529 // feed the parameter vector
530 std::vector<Double_t> par;
531 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
532 for (UInt_t i=0; i<paramList->size(); i++)
533 par.push_back((*paramList)[i].fValue);
534
535 // calculate functions
536 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
537 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
538 }
539
540 // calculate asymmetry
541 Double_t asymFcnValue = 0.0;
542 Double_t a, b, f, alphaest;
543 Double_t time;
544
545 // Get estimated alpha
546 alphaest = fRunInfo->GetEstimatedAlpha();
547
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)) - (-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))-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))-(-f*(a*b+1.0)-(a-1.0))/((a+1.0)+f*(a*b-1.0));
597 break;
598 case 5: // alpha ?? , beta == 1
599 a = alphaest;
600 f = fTheory->Func(time, par, fFuncValues);
601 asymFcnValue = (f*(a+1.0)-(a-1.0))/((a+1.0)-f*(a-1.0)) - (-f*(a+1.0)-(a-1.0))/((a+1.0)+f*(a-1.0));
602 break;
603 case 6: // alpha ??, beta != 1
604 a = alphaest;
605 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
606 b = par[fRunInfo->GetBetaParamNo()-1];
607 } else { // beta is a function
608 // get function number
609 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
610 // evaluate function
611 b = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
612 }
613 f = fTheory->Func(time, par, fFuncValues);
614 asymFcnValue = (f*(a*b+1.0)-(a-1.0))/((a+1.0)-f*(a*b-1.0))-(-f*(a*b+1.0)-(a-1.0))/((a+1.0)+f*(a*b-1.0));
615 break;
616 default:
617 asymFcnValue = 0.0;
618 break;
619 }
620 fData.AppendTheoryValue(asymFcnValue);
621 }
622
623 // clean up
624 par.clear();
625}
626
627//--------------------------------------------------------------------------
628// PrepareData (protected)
629//--------------------------------------------------------------------------
652{
653 if (!fValid)
654 return false;
655
656 // keep the Global block info
657 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
658
659 // get forward/backward histo from PRunDataHandler object ------------------------
660 // get the correct run
661 PRawRunData *runData = fRawData->GetRunData(*(fRunInfo->GetRunName()));
662 if (!runData) { // run not found
663 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
664 std::cerr << std::endl;
665 return false;
666 }
667
668 // keep the field from the meta-data from the data-file
669 fMetaData.fField = runData->GetField();
670
671 // keep the energy from the meta-data from the data-file
672 fMetaData.fEnergy = runData->GetEnergy();
673
674 // keep the temperature(s) from the meta-data from the data-file
675 for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
676 fMetaData.fTemp.push_back(runData->GetTemperature(i));
677
678 // collect histogram numbers
679 PUIntVector forwardHistoNo;
680 PUIntVector backwardHistoNo;
681 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
682 forwardHistoNo.push_back(fRunInfo->GetForwardHistoNo(i));
683
684 if (!runData->IsPresent(forwardHistoNo[i])) {
685 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): **PANIC ERROR**:";
686 std::cerr << std::endl << ">> forwardHistoNo found = " << forwardHistoNo[i] << ", which is NOT present in the data file!?!?";
687 std::cerr << std::endl << ">> Will quit :-(";
688 std::cerr << std::endl;
689 // clean up
690 forwardHistoNo.clear();
691 backwardHistoNo.clear();
692 return false;
693 }
694 }
695 for (UInt_t i=0; i<fRunInfo->GetBackwardHistoNoSize(); i++) {
696 backwardHistoNo.push_back(fRunInfo->GetBackwardHistoNo(i));
697
698 if (!runData->IsPresent(backwardHistoNo[i])) {
699 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): **PANIC ERROR**:";
700 std::cerr << std::endl << ">> backwardHistoNo found = " << backwardHistoNo[i] << ", which is NOT present in the data file!?!?";
701 std::cerr << std::endl << ">> Will quit :-(";
702 std::cerr << std::endl;
703 // clean up
704 forwardHistoNo.clear();
705 backwardHistoNo.clear();
706 return false;
707 }
708 }
709/* //as35
710 if (forwardHistoNo.size() != backwardHistoNo.size()) {
711 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): **PANIC ERROR**:";
712 std::cerr << std::endl << ">> # of forward histograms different from # of backward histograms.";
713 std::cerr << std::endl << ">> Will quit :-(";
714 std::cerr << std::endl;
715 // clean up
716 forwardHistoNo.clear();
717 backwardHistoNo.clear();
718 return false;
719 }
720*/ //as35
721
722 // keep the time resolution in (s)
723 fTimeResolution = runData->GetTimeResolution()/1.0e3;
724 std::cout.precision(10);
725 std::cout << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ms)" << std::endl;
726
727 // get all the proper t0's and addt0's for the current RUN block
728 if (!GetProperT0(runData, globalBlock, forwardHistoNo, backwardHistoNo)) {
729 return false;
730 }
731
732 // keep the histo of each group at this point (addruns handled below)
733 std::vector<PDoubleVector> forward, backward;
734 forward.resize(forwardHistoNo.size()); // resize to number of groups
735 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
736 forward[i].resize(runData->GetDataBin(forwardHistoNo[i])->size());
737 forward[i] = *runData->GetDataBin(forwardHistoNo[i]);
738 }
739 backward.resize(backwardHistoNo.size()); // resize to number of groups
740 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
741 backward[i].resize(runData->GetDataBin(backwardHistoNo[i])->size());
742 backward[i] = *runData->GetDataBin(backwardHistoNo[i]);
743 }
744
745 // check if addrun's are present, and if yes add data
746 // check if there are runs to be added to the current one
747 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
748 PRawRunData *addRunData;
749 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
750 // get run to be added to the main one
751 addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
752 if (addRunData == nullptr) { // couldn't get run
753 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
754 std::cerr << std::endl;
755 return false;
756 }
757
758 // add forward run
759 UInt_t addRunSize;
760 for (UInt_t k=0; k<forwardHistoNo.size(); k++) { // fill each group
761 addRunSize = addRunData->GetDataBin(forwardHistoNo[k])->size();
762 for (UInt_t j=0; j<addRunData->GetDataBin(forwardHistoNo[k])->size(); j++) { // loop over the bin indices
763 // make sure that the index stays in the proper range
764 if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]) >= 0) &&
765 (j+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]) < addRunSize)) {
766 forward[k][j] += addRunData->GetDataBin(forwardHistoNo[k])->at(j+static_cast<Int_t>(fAddT0s[i-1][2*k])-static_cast<Int_t>(fT0s[2*k]));
767 }
768 }
769 }
770
771 // add backward run
772 for (UInt_t k=0; k<backwardHistoNo.size(); k++) { // fill each group
773 addRunSize = addRunData->GetDataBin(backwardHistoNo[k])->size();
774 for (UInt_t j=0; j<addRunData->GetDataBin(backwardHistoNo[k])->size(); j++) { // loop over the bin indices
775 // make sure that the index stays in the proper range
776 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) &&
777 (j+static_cast<Int_t>(fAddT0s[i-1][2*k+1])-static_cast<Int_t>(fT0s[2*k+1]) < addRunSize)) {
778 backward[k][j] += addRunData->GetDataBin(backwardHistoNo[k])->at(j+static_cast<Int_t>(fAddT0s[i-1][2*k+1])-static_cast<Int_t>(fT0s[2*k+1]));
779 }
780 }
781 }
782 }
783 }
784
785 // set forward histo data of the first group
786 fForwardp.resize(forward[0].size());
787 fForwardm.resize(forward[0].size());
788 for (UInt_t i=0; i<fForwardp.size(); i++) {
789 fForwardp[i] = forward[0][i];
790 fForwardm[i] = forward[1][i];
791 }
792 // set backward histo data of the first group
793 fBackwardp.resize(backward[0].size());
794 fBackwardm.resize(backward[0].size());
795 for (UInt_t i=0; i<fBackwardp.size(); i++) {
796 fBackwardp[i] = backward[0][i];
797 fBackwardm[i] = backward[1][i];
798 }
799
800 // subtract background from histogramms ------------------------------------------
801 if (fRunInfo->GetBkgFix(0) == PMUSR_UNDEFINED) { // no fixed background given
802 if (fRunInfo->GetBkgRange(0) >= 0) { // background range given
804 return false;
805 } else { // no background given to do the job, try to estimate it
806 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
807 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
808 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.1), 2);
809 fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[1]*0.6), 3);
810 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): **WARNING** Neither fix background nor background bins are given!";
811 std::cerr << std::endl << ">> Will try the following:";
812 std::cerr << std::endl << ">> forward: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
813 std::cerr << std::endl << ">> backward: bkg start = " << fRunInfo->GetBkgRange(2) << ", bkg end = " << fRunInfo->GetBkgRange(3);
814 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS MAKES ANY SENSE! Better check ...";
815 std::cerr << std::endl;
817 return false;
818 }
819 } else { // fixed background given
820 if (!SubtractFixBkg())
821 return false;
822 }
823
824 UInt_t histoNo[2] = {forwardHistoNo[0], backwardHistoNo[0]};
825
826 // get the data range (fgb/lgb) for the current RUN block
827 if (!GetProperDataRange(runData, histoNo)) {
828 return false;
829 }
830
831 // get the fit range for the current RUN block
832 GetProperFitRange(globalBlock);
833
834 // everything looks fine, hence fill data set
835 Bool_t status;
836 switch(fHandleTag) {
837 case kFit:
839 break;
840 case kView:
841 status = PrepareViewData(runData, histoNo);
842 break;
843 default:
844 status = false;
845 break;
846 }
847
848 // clean up
849 forwardHistoNo.clear();
850 backwardHistoNo.clear();
851
852 return status;
853}
854
855//--------------------------------------------------------------------------
856// SubtractFixBkg (private)
857//--------------------------------------------------------------------------
875{
876 Double_t dval;
877
878 // Order in RunInfo structure Fp, Fm, Bp, Bm
879 for (UInt_t i=0; i<fForwardp.size(); i++) {
880 // keep the error, and make sure that the bin is NOT empty
881 if (fForwardp[i] != 0.0)
882 dval = TMath::Sqrt(fForwardp[i]);
883 else
884 dval = 1.0;
885 fForwardpErr.push_back(dval);
886 fForwardp[i] -= fRunInfo->GetBkgFix(0);
887
888 // keep the error, and make sure that the bin is NOT empty
889 if (fForwardm[i] != 0.0)
890 dval = TMath::Sqrt(fForwardm[i]);
891 else
892 dval = 1.0;
893 fForwardmErr.push_back(dval);
894 fForwardm[i] -= fRunInfo->GetBkgFix(1);
895
896 // keep the error, and make sure that the bin is NOT empty
897 if (fBackwardp[i] != 0.0)
898 dval = TMath::Sqrt(fBackwardp[i]);
899 else
900 dval = 1.0;
901 fBackwardpErr.push_back(dval);
902 fBackwardp[i] -= fRunInfo->GetBkgFix(2);
903
904 // keep the error, and make sure that the bin is NOT empty
905 if (fBackwardm[i] != 0.0)
906 dval = TMath::Sqrt(fBackwardm[i]);
907 else
908 dval = 1.0;
909 fBackwardmErr.push_back(dval);
910 fBackwardm[i] -= fRunInfo->GetBkgFix(3);
911 }
912
913 return true;
914}
915
916//--------------------------------------------------------------------------
917// SubtractEstimatedBkg (private)
918//--------------------------------------------------------------------------
937{
938 Double_t beamPeriod = 0.0;
939
940 // check if data are from PSI, RAL, or TRIUMF
941 if (fRunInfo->GetInstitute()->Contains("psi"))
942 beamPeriod = ACCEL_PERIOD_PSI;
943 else if (fRunInfo->GetInstitute()->Contains("ral"))
944 beamPeriod = ACCEL_PERIOD_RAL;
945 else if (fRunInfo->GetInstitute()->Contains("triumf"))
946 beamPeriod = ACCEL_PERIOD_TRIUMF;
947 else
948 beamPeriod = 0.0;
949
950 // check if start and end are in proper order
951 Int_t start[2] = {fRunInfo->GetBkgRange(0), fRunInfo->GetBkgRange(2)};
952 Int_t end[2] = {fRunInfo->GetBkgRange(1), fRunInfo->GetBkgRange(3)};
953 for (UInt_t i=0; i<2; i++) {
954 if (end[i] < start[i]) {
955 std::cout << std::endl << "PRunAsymmetryBNMR::SubtractEstimatedBkg(): end = " << end[i] << " > start = " << start[i] << "! Will swap them!";
956 UInt_t keep = end[i];
957 end[i] = start[i];
958 start[i] = keep;
959 }
960 }
961
962 // calculate proper background range
963 for (UInt_t i=0; i<2; i++) {
964 if (beamPeriod != 0.0) {
965 Double_t timeBkg = static_cast<Double_t>(end[i]-start[i])*(fTimeResolution*fPacking); // length of the background intervall in time
966 UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
967 // correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
968 end[i] = start[i] + static_cast<UInt_t>((fullCycles*beamPeriod)/(fTimeResolution*fPacking));
969 std::cout << "PRunAsymmetryBNMR::SubtractEstimatedBkg(): Background " << start[i] << ", " << end[i] << std::endl;
970 if (end[i] == start[i])
971 end[i] = fRunInfo->GetBkgRange(2*i+1);
972 }
973 }
974
975 // check if start is within histogram bounds
976 if ((start[0] < 0) || (start[0] >= fForwardp.size()) ||
977 (start[1] < 0) || (start[1] >= fBackwardp.size())) {
978 std::cerr << std::endl << ">> PRunAsymmetryBNMR::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
979 std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForwardp.size() << "/" << fBackwardp.size() << ").";
980 std::cerr << std::endl << ">> background start (f/b) = (" << start[0] << "/" << start[1] << ").";
981 return false;
982 }
983
984 // check if end is within histogram bounds
985 if ((end[0] < 0) || (end[0] >= fForwardp.size()) ||
986 (end[1] < 0) || (end[1] >= fBackwardp.size())) {
987 std::cerr << std::endl << ">> PRunAsymmetryBNMR::SubtractEstimatedBkg(): **ERROR** background bin values out of bound!";
988 std::cerr << std::endl << ">> histo lengths (f/b) = (" << fForwardp.size() << "/" << fBackwardp.size() << ").";
989 std::cerr << std::endl << ">> background end (f/b) = (" << end[0] << "/" << end[1] << ").";
990 return false;
991 }
992
993 // calculate background
994 Double_t bkgp[2] = {0.0, 0.0};
995 Double_t errBkgp[2] = {0.0, 0.0};
996 Double_t bkgm[2] = {0.0, 0.0};
997 Double_t errBkgm[2] = {0.0, 0.0};
998
999 // forward
1000 for (UInt_t i=start[0]; i<=end[0]; i++) {
1001 bkgp[0] += fForwardp[i];
1002 bkgm[0] += fForwardm[i];
1003 }
1004 errBkgp[0] = TMath::Sqrt(bkgp[0])/(end[0] - start[0] + 1);
1005 bkgp[0] /= static_cast<Double_t>(end[0] - start[0] + 1);
1006 std::cout << std::endl << ">> estimated pos hel forward histo background: " << bkgp[0];
1007 errBkgm[0] = TMath::Sqrt(bkgp[0])/(end[0] - start[0] + 1);
1008 bkgm[0] /= static_cast<Double_t>(end[0] - start[0] + 1);
1009 std::cout << std::endl << ">> estimated neg hel forward histo background: " << bkgm[0];
1010
1011 // backward
1012 for (UInt_t i=start[1]; i<=end[1]; i++) {
1013 bkgp[1] += fBackwardp[i];
1014 bkgm[1] += fBackwardm[i];
1015 }
1016 errBkgp[1] = TMath::Sqrt(bkgp[1])/(end[1] - start[1] + 1);
1017 bkgp[1] /= static_cast<Double_t>(end[1] - start[1] + 1);
1018 std::cout << std::endl << ">> estimated pos hel backward histo background: " << bkgp[1];
1019 errBkgm[1] = TMath::Sqrt(bkgm[1])/(end[1] - start[1] + 1);
1020 bkgm[1] /= static_cast<Double_t>(end[1] - start[1] + 1);
1021 std::cout << std::endl << ">> estimated neg hel backward histo background: " << bkgm[1];
1022
1023 // correct error for forward, backward
1024 Double_t errVal = 0.0;
1025 for (UInt_t i=0; i<fForwardp.size(); i++) {
1026 if (fForwardp[i] > 0.0)
1027 errVal = TMath::Sqrt(fForwardp[i]+errBkgp[0]*errBkgp[0]);
1028 else
1029 errVal = 1.0;
1030 fForwardpErr.push_back(errVal);
1031 if (fBackwardp[i] > 0.0)
1032 errVal = TMath::Sqrt(fBackwardp[i]+errBkgp[1]*errBkgp[1]);
1033 else
1034 errVal = 1.0;
1035 fBackwardpErr.push_back(errVal);
1036 if (fForwardm[i] > 0.0)
1037 errVal = TMath::Sqrt(fForwardm[i]+errBkgm[0]*errBkgm[0]);
1038 else
1039 errVal = 1.0;
1040 fForwardmErr.push_back(errVal);
1041 if (fBackwardm[i] > 0.0)
1042 errVal = TMath::Sqrt(fBackwardm[i]+errBkgm[1]*errBkgm[1]);
1043 else
1044 errVal = 1.0;
1045 fBackwardmErr.push_back(errVal);
1046 }
1047
1048 // subtract background from data
1049 for (UInt_t i=0; i<fForwardp.size(); i++) {
1050 fForwardp[i] -= bkgp[0];
1051 fBackwardp[i] -= bkgp[1];
1052 fForwardm[i] -= bkgm[0];
1053 fBackwardm[i] -= bkgm[1];
1054 }
1055
1056 fRunInfo->SetBkgEstimated(bkgp[0], 0);
1057 fRunInfo->SetBkgEstimated(bkgp[1], 1);
1058 fRunInfo->SetBkgEstimated(bkgm[0], 3);
1059 fRunInfo->SetBkgEstimated(bkgm[1], 4);
1060
1061 // Get estimate for alpha once
1062 Double_t alpha;
1063 alpha = EstimateAlpha();
1064 fRunInfo->SetEstimatedAlpha(alpha);
1065
1066
1067 return true;
1068}
1069
1070//--------------------------------------------------------------------------
1071// PrepareFitData (protected)
1072//--------------------------------------------------------------------------
1085{
1086 // transform raw histo data. This is done the following way (for details see the manual):
1087 // first rebin the data, than calculate the asymmetry
1088
1089 // everything looks fine, hence fill packed forward and backward histo
1090 PRunData forwardpPacked;
1091 PRunData backwardpPacked;
1092 PRunData forwardmPacked;
1093 PRunData backwardmPacked;
1094 Double_t valuep = 0.0;
1095 Double_t errorp = 0.0;
1096 Double_t valuem = 0.0;
1097 Double_t errorm = 0.0;
1098
1099 // forward
1100 for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
1101 if (fPacking == 1) {
1102 forwardpPacked.AppendValue(fForwardp[i]);
1103 forwardpPacked.AppendErrorValue(fForwardpErr[i]);
1104 forwardmPacked.AppendValue(fForwardm[i]);
1105 forwardmPacked.AppendErrorValue(fForwardmErr[i]);
1106 } else { // packed data, i.e. fPacking > 1
1107 if (((i-fGoodBins[0]) % fPacking == 0) && (i != fGoodBins[0])) { // fill data
1108 // in order that after rebinning the fit does not need to be redone (important for plots)
1109 // the value is normalize to per bin
1110 valuep /= fPacking;
1111 valuem /= fPacking;
1112 forwardpPacked.AppendValue(valuep);
1113 forwardmPacked.AppendValue(valuem);
1114 if (valuep == 0.0) {
1115 forwardpPacked.AppendErrorValue(1.0);
1116 } else {
1117 forwardpPacked.AppendErrorValue(TMath::Sqrt(errorp)/fPacking);
1118 }
1119 if (valuem == 0.0) {
1120 forwardmPacked.AppendErrorValue(1.0);
1121 } else {
1122 forwardmPacked.AppendErrorValue(TMath::Sqrt(errorm)/fPacking);
1123 }
1124 valuep = 0.0;
1125 errorp = 0.0;
1126 valuem = 0.0;
1127 errorm = 0.0;
1128 }
1129 valuep += fForwardp[i];
1130 errorp += fForwardpErr[i]*fForwardpErr[i];
1131 valuem += fForwardm[i];
1132 errorm += fForwardmErr[i]*fForwardmErr[i];
1133 }
1134 }
1135
1136 // backward
1137 for (Int_t i=fGoodBins[2]; i<fGoodBins[3]; i++) {
1138 if (fPacking == 1) {
1139 backwardpPacked.AppendValue(fBackwardp[i]);
1140 backwardpPacked.AppendErrorValue(fBackwardpErr[i]);
1141 backwardmPacked.AppendValue(fBackwardm[i]);
1142 backwardmPacked.AppendErrorValue(fBackwardmErr[i]);
1143 } else { // packed data, i.e. fPacking > 1
1144 if (((i-fGoodBins[2]) % fPacking == 0) && (i != fGoodBins[2])) { // fill data
1145 // in order that after rebinning the fit does not need to be redone (important for plots)
1146 // the value is normalize to per bin
1147 valuep /= fPacking;
1148 valuem /= fPacking;
1149 backwardpPacked.AppendValue(valuep);
1150 backwardmPacked.AppendValue(valuem);
1151 if (valuep == 0.0) {
1152 backwardpPacked.AppendErrorValue(1.0);
1153 } else {
1154 backwardpPacked.AppendErrorValue(TMath::Sqrt(errorp)/fPacking);
1155 }
1156 if (valuem == 0.0) {
1157 backwardmPacked.AppendErrorValue(1.0);
1158 } else {
1159 backwardmPacked.AppendErrorValue(TMath::Sqrt(errorm)/fPacking);
1160 }
1161 valuep = 0.0;
1162 errorp = 0.0;
1163 valuem = 0.0;
1164 errorm = 0.0;
1165 }
1166 valuep += fBackwardp[i];
1167 errorp += fBackwardpErr[i]*fBackwardpErr[i];
1168 valuem += fBackwardm[i];
1169 errorm += fBackwardmErr[i]*fBackwardmErr[i];
1170 }
1171 }
1172
1173 // check if packed forward and backward hist have the same size, otherwise take the minimum size
1174 UInt_t noOfBins = forwardpPacked.GetValue()->size();
1175 if (forwardpPacked.GetValue()->size() != backwardpPacked.GetValue()->size()) {
1176 if (forwardpPacked.GetValue()->size() > backwardpPacked.GetValue()->size())
1177 noOfBins = backwardpPacked.GetValue()->size();
1178 }
1179
1180 // form asymmetry including error propagation
1181 Double_t asym,error;
1182 Double_t fp, bp, efp, ebp;
1183 Double_t fm, bm, efm, ebm;
1184 // fill data time start, and step
1185 // data start at data_start-t0 shifted by (pack-1)/2
1186 fData.SetDataTimeStart(fTimeResolution*(static_cast<Double_t>(fGoodBins[0])-fT0s[0]+static_cast<Double_t>(fPacking-1)/2.0));
1187 fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(fPacking));
1188 for (UInt_t i=0; i<noOfBins; i++) {
1189 // to make the formulae more readable
1190 fp = forwardpPacked.GetValue()->at(i);
1191 bp = backwardpPacked.GetValue()->at(i);
1192 efp = forwardpPacked.GetError()->at(i);
1193 ebp = backwardpPacked.GetError()->at(i);
1194 fm = forwardmPacked.GetValue()->at(i);
1195 bm = backwardmPacked.GetValue()->at(i);
1196 efm = forwardmPacked.GetError()->at(i);
1197 ebm = backwardmPacked.GetError()->at(i);
1198 // check that there are indeed bins
1199 if (fp+bp != 0.0)
1200 asym = (fp-bp) / (fp+bp) - (fm-bm) / (fm+bm);
1201 else
1202 asym = 0.0;
1203 fData.AppendValue(asym);
1204 // calculate the error
1205 if (fp+bp != 0.0)
1206 errorp = 2.0/((fp+bp)*(fp+bp))*TMath::Sqrt(bp*bp*efp*efp+ebp*ebp*fp*fp);
1207 else
1208 errorp = 1.0;
1209 if (fm+bm != 0.0)
1210 errorm = 2.0/((fm+bm)*(fm+bm))*TMath::Sqrt(bm*bm*efm*efm+ebm*ebm*fm*fm);
1211 else
1212 errorp = 1.0;
1213
1214 error = TMath::Sqrt(errorp*errorp+errorm*errorm);
1215 fData.AppendErrorValue(error);
1216 }
1217
1219
1220 // clean up
1221 fForwardp.clear();
1222 fForwardpErr.clear();
1223 fBackwardp.clear();
1224 fBackwardpErr.clear();
1225 fForwardm.clear();
1226 fForwardmErr.clear();
1227 fBackwardm.clear();
1228 fBackwardmErr.clear();
1229
1230 return true;
1231}
1232//--------------------------------------------------------------------------
1246Bool_t PRunAsymmetryBNMR::PrepareViewData(PRawRunData* runData, UInt_t histoNo[2])
1247{
1248 // check if view_packing is wished
1249 Int_t packing = fPacking;
1250 if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
1251 packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1252 }
1253
1254 // feed the parameter vector
1255 std::vector<Double_t> par;
1256 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
1257 for (UInt_t i=0; i<paramList->size(); i++)
1258 par.push_back((*paramList)[i].fValue);
1259
1260 // transform raw histo data. This is done the following way (for details see the manual):
1261 // first rebin the data, than calculate the asymmetry
1262
1263 // first get start data, end data, and t0
1264 Int_t start[2] = {fGoodBins[0], fGoodBins[2]};
1265 Int_t end[2] = {fGoodBins[1], fGoodBins[3]};
1266 Int_t t0[2] = {static_cast<Int_t>(fT0s[0]), static_cast<Int_t>(fT0s[1])};
1267
1268 // check if the data ranges and t0's between forward/backward are compatible
1269 Int_t fgb[2];
1270 if (start[0]-t0[0] != start[1]-t0[1]) { // wrong fgb aligning
1271 if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
1272 fgb[0] = start[0];
1273 fgb[1] = t0[1] + start[0]-t0[0];
1274 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareViewData(): **WARNING** needed to shift backward fgb from ";
1275 std::cerr << start[1] << " to " << fgb[1] << std::endl;
1276 } else {
1277 fgb[0] = t0[0] + start[1]-t0[1];
1278 fgb[1] = start[1];
1279 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareViewData(): **WARNING** needed to shift forward fgb from ";
1280 std::cerr << start[0] << " to " << fgb[0] << std::endl;
1281 }
1282 } else { // fgb aligning is correct
1283 fgb[0] = start[0];
1284 fgb[1] = start[1];
1285 }
1286
1287 Int_t val = fgb[0]-packing*(fgb[0]/packing);
1288 do {
1289 if (fgb[1] - fgb[0] < 0)
1290 val += packing;
1291 } while (val + fgb[1] - fgb[0] < 0);
1292
1293 start[0] = val;
1294 start[1] = val + fgb[1] - fgb[0];
1295
1296 // make sure that there are equal number of rebinned bins in forward and backward
1297 UInt_t noOfBins0 = (runData->GetDataBin(histoNo[0])->size()-start[0])/packing;
1298 UInt_t noOfBins1 = (runData->GetDataBin(histoNo[1])->size()-start[1])/packing;
1299 if (noOfBins0 > noOfBins1)
1300 noOfBins0 = noOfBins1;
1301 end[0] = start[0] + noOfBins0 * packing;
1302 end[1] = start[1] + noOfBins0 * packing;
1303
1304 // check if start, end, and t0 make any sense
1305 // 1st check if start and end are in proper order
1306 for (UInt_t i=0; i<2; i++) {
1307 if (end[i] < start[i]) { // need to swap them
1308 Int_t keep = end[i];
1309 end[i] = start[i];
1310 start[i] = keep;
1311 }
1312 // 2nd check if start is within proper bounds
1313 if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1314 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
1315 std::cerr << std::endl;
1316 return false;
1317 }
1318 // 3rd check if end is within proper bounds
1319 if ((end[i] < 0) || (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1320 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
1321 std::cerr << std::endl;
1322 return false;
1323 }
1324 // 4th check if t0 is within proper bounds
1325 if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1326 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareViewData(): **ERROR** t0 data bin doesn't make any sense!";
1327 std::cerr << std::endl;
1328 return false;
1329 }
1330 }
1331
1332 // everything looks fine, hence fill packed forward and backward histo
1333 PRunData forwardpPacked;
1334 PRunData backwardpPacked;
1335 PRunData forwardmPacked;
1336 PRunData backwardmPacked;
1337 Double_t valuep = 0.0;
1338 Double_t errorp = 0.0;
1339 Double_t valuem = 0.0;
1340 Double_t errorm = 0.0;
1341 Double_t value = 0.0;
1342 Double_t error = 0.0;
1343
1344 // forward
1345 for (Int_t i=start[0]; i<end[0]; i++) {
1346 if (packing == 1) {
1347 forwardpPacked.AppendValue(fForwardp[i]);
1348 forwardpPacked.AppendErrorValue(fForwardpErr[i]);
1349 forwardmPacked.AppendValue(fForwardm[i]);
1350 forwardmPacked.AppendErrorValue(fForwardmErr[i]);
1351 } else { // packed data, i.e. packing > 1
1352 if (((i-start[0]) % packing == 0) && (i != start[0])) { // fill data
1353 // in order that after rebinning the fit does not need to be redone (important for plots)
1354 // the value is normalize to per bin
1355 valuep /= packing;
1356 forwardpPacked.AppendValue(valuep);
1357 valuem /= packing;
1358 forwardmPacked.AppendValue(valuem);
1359 if (valuep == 0.0) {
1360 forwardpPacked.AppendErrorValue(1.0);
1361 } else {
1362 forwardpPacked.AppendErrorValue(TMath::Sqrt(errorp)/packing);
1363 }
1364 if (valuem == 0.0) {
1365 forwardmPacked.AppendErrorValue(1.0);
1366 } else {
1367 forwardmPacked.AppendErrorValue(TMath::Sqrt(errorm)/packing);
1368 }
1369 valuep = 0.0;
1370 errorp = 0.0;
1371 valuem = 0.0;
1372 errorm = 0.0;
1373 }
1374 valuep += fForwardp[i];
1375 errorp += fForwardpErr[i]*fForwardpErr[i];
1376 valuem += fForwardm[i];
1377 errorm += fForwardmErr[i]*fForwardmErr[i];
1378 }
1379 }
1380
1381 // backward
1382 for (Int_t i=start[1]; i<end[1]; i++) {
1383 if (packing == 1) {
1384 backwardpPacked.AppendValue(fBackwardp[i]);
1385 backwardpPacked.AppendErrorValue(fBackwardpErr[i]);
1386 backwardmPacked.AppendValue(fBackwardm[i]);
1387 backwardmPacked.AppendErrorValue(fBackwardmErr[i]);
1388 } else { // packed data, i.e. packing > 1
1389 if (((i-start[1]) % packing == 0) && (i != start[1])) { // fill data
1390 // in order that after rebinning the fit does not need to be redone (important for plots)
1391 // the value is normalize to per bin
1392 valuep /= packing;
1393 valuem /= packing;
1394 backwardpPacked.AppendValue(valuep);
1395 backwardmPacked.AppendValue(valuem);
1396 if (valuep == 0.0) {
1397 backwardpPacked.AppendErrorValue(1.0);
1398 } else {
1399 backwardpPacked.AppendErrorValue(TMath::Sqrt(errorp)/packing);
1400 }
1401 if (valuem == 0.0) {
1402 backwardmPacked.AppendErrorValue(1.0);
1403 } else {
1404 backwardmPacked.AppendErrorValue(TMath::Sqrt(errorm)/packing);
1405 }
1406 valuep = 0.0;
1407 errorp = 0.0;
1408 valuem = 0.0;
1409 errorm = 0.0;
1410 }
1411 valuep += fBackwardp[i];
1412 errorp += fBackwardpErr[i]*fBackwardpErr[i];
1413 valuem += fBackwardm[i];
1414 errorm += fBackwardmErr[i]*fBackwardmErr[i];
1415 }
1416 }
1417
1418 // check if packed forward and backward hist have the same size, otherwise take the minimum size
1419 UInt_t noOfBins = forwardpPacked.GetValue()->size();
1420 if (forwardpPacked.GetValue()->size() != backwardpPacked.GetValue()->size()) {
1421 if (forwardpPacked.GetValue()->size() > backwardpPacked.GetValue()->size())
1422 noOfBins = backwardpPacked.GetValue()->size();
1423 }
1424
1425 // form asymmetry including error propagation
1426 Double_t asym;
1427 Double_t fp, bp, efp, ebp, alpha, beta = 1.0;
1428 Double_t fm, bm, efm, ebm;
1429 // set data time start, and step
1430 // data start at data_start-t0
1431 fData.SetDataTimeStart(fTimeResolution*(static_cast<Double_t>(start[0])-t0[0]+static_cast<Double_t>(packing-1)/2.0));
1432 fData.SetDataTimeStep(fTimeResolution*static_cast<Double_t>(packing));
1433
1434 // Get estimate for alpha once
1435 alpha = EstimateAlpha();
1436
1437 // get the proper alpha and beta
1438 switch (fAlphaBetaTag) {
1439 case 1: // alpha == 1, beta == 1
1440 alpha = 1.0;
1441 beta = 1.0;
1442 break;
1443 case 2: // alpha != 1, beta == 1
1444 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1445 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1446 } else { // alpha is function
1447 // get function number
1448 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1449 // evaluate function
1450 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1451 }
1452 beta = 1.0;
1453 break;
1454 case 3: // alpha == 1, beta != 1
1455 alpha = 1.0;
1456 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1457 beta = par[fRunInfo->GetBetaParamNo()-1];
1458 } else { // beta is a function
1459 // get function number
1460 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1461 // evaluate function
1462 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1463 }
1464 break;
1465 case 4: // alpha != 1, beta != 1
1466 if (fRunInfo->GetAlphaParamNo() < MSR_PARAM_FUN_OFFSET) { // alpha is a parameter
1467 alpha = par[fRunInfo->GetAlphaParamNo()-1];
1468 } else { // alpha is function
1469 // get function number
1470 UInt_t funNo = fRunInfo->GetAlphaParamNo()-MSR_PARAM_FUN_OFFSET;
1471 // evaluate function
1472 alpha = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1473 }
1474 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1475 beta = par[fRunInfo->GetBetaParamNo()-1];
1476 } else { // beta is a function
1477 // get function number
1478 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1479 // evaluate function
1480 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1481 }
1482 break;
1483 case 5: // alpha ?? , beta == 1
1484 // use estimated value
1485 beta = 1.0;
1486 break;
1487 case 6: // alpha ??, beta != 1
1488 // use estimated value
1489 if (fRunInfo->GetBetaParamNo() < MSR_PARAM_FUN_OFFSET) { // beta is a parameter
1490 beta = par[fRunInfo->GetBetaParamNo()-1];
1491 } else { // beta is a function
1492 // get function number
1493 UInt_t funNo = fRunInfo->GetBetaParamNo()-MSR_PARAM_FUN_OFFSET;
1494 // evaluate function
1495 beta = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par, fMetaData);
1496 }
1497 break;
1498 default:
1499 break;
1500 }
1501
1502 for (UInt_t i=0; i<forwardpPacked.GetValue()->size(); i++) {
1503 // to make the formulae more readable
1504 fp = forwardpPacked.GetValue()->at(i);
1505 bp = backwardpPacked.GetValue()->at(i);
1506 efp = forwardpPacked.GetError()->at(i);
1507 ebp = backwardpPacked.GetError()->at(i);
1508 fm = forwardmPacked.GetValue()->at(i);
1509 bm = backwardmPacked.GetValue()->at(i);
1510 efm = forwardmPacked.GetError()->at(i);
1511 ebm = backwardmPacked.GetError()->at(i);
1512 // check that there are indeed bins
1513 if (fp+bp != 0.0)
1514 asym = (alpha*fp-bp) / (alpha*beta*fp+bp) - (alpha*fm-bm) / (alpha*beta*fm+bm);
1515 else
1516 asym = 0.0;
1517 fData.AppendValue(asym);
1518 // calculate the error
1519 if (fp+bp != 0.0)
1520 errorp = 2.0/((fp+bp)*(fp+bp))*TMath::Sqrt(bp*bp*efp*efp+ebp*ebp*fp*fp);
1521 else
1522 errorp = 1.0;
1523 if (fm+bm != 0.0)
1524 errorm = 2.0/((fm+bm)*(fm+bm))*TMath::Sqrt(bm*bm*efm*efm+ebm*ebm*fm*fm);
1525 else
1526 errorm = 1.0;
1527 error = TMath::Sqrt(errorp*errorp+errorm*errorm);
1528 fData.AppendErrorValue(error);
1529 }
1530
1532
1533 // clean up
1534 fForwardp.clear();
1535 fForwardpErr.clear();
1536 fBackwardp.clear();
1537 fBackwardpErr.clear();
1538 fForwardm.clear();
1539 fForwardmErr.clear();
1540 fBackwardm.clear();
1541 fBackwardmErr.clear();
1542
1543 // fill theory vector for kView
1544 // calculate functions
1545 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
1546 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
1547 }
1548
1549 // calculate theory
1550 UInt_t size = runData->GetDataBin(histoNo[0])->size();
1551
1552 Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
1553 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
1554 if (fTheoAsData) { // calculate theory only at the data points
1555 fData.SetTheoryTimeStep(fData.GetDataTimeStep());
1556 } else {
1557 // finer binning for the theory (8 times as many points = factor)
1558 size *= factor;
1559 fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
1560 }
1561
1562 Double_t time;
1563 for (UInt_t i=0; i<size; i++) {
1564 time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
1565 value = fTheory->Func(time, par, fFuncValues);
1566 if (fabs(value) > 10.0) { // dirty hack needs to be fixed!!
1567 value = 0.0;
1568 }
1569 fData.AppendTheoryValue(value);
1570 }
1571
1572 // clean up
1573 par.clear();
1574
1575 return true;
1576}
1577
1578
1579//--------------------------------------------------------------------------
1580// GetProperT0 (private)
1581//--------------------------------------------------------------------------
1600Bool_t PRunAsymmetryBNMR::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHistoNo, PUIntVector &backwardHistoNo)
1601{
1602 // feed all T0's
1603 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1604 fT0s.clear();
1605 // this strange fT0 size estimate is needed in case #forw histos != #back histos
1606 size_t size = 2*forwardHistoNo.size();
1607 if (backwardHistoNo.size() > forwardHistoNo.size())
1608 size = 2*backwardHistoNo.size();
1609 fT0s.resize(size);
1610 for (UInt_t i=0; i<fT0s.size(); i++) {
1611 fT0s[i] = -1.0;
1612 }
1613
1614 // fill in the T0's from the msr-file (if present)
1615 for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
1616 fT0s[i] = fRunInfo->GetT0Bin(i);
1617 }
1618
1619 // fill in the missing T0's from the GLOBAL block section (if present)
1620 for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
1621 if (fT0s[i] == -1) { // i.e. not given in the RUN block section
1622 fT0s[i] = globalBlock->GetT0Bin(i);
1623 }
1624 }
1625
1626 // fill in the missing T0's from the data file, if not already present in the msr-file
1627 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1628 if (fT0s[2*i] == -1.0) // i.e. not present in the msr-file, try the data file
1629 if (runData->GetT0Bin(forwardHistoNo[i]) > 0.0) {
1630 fT0s[2*i] = runData->GetT0Bin(forwardHistoNo[i]);
1631 fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1632 }
1633 }
1634 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1635 if (fT0s[2*i+1] == -1.0) // i.e. not present in the msr-file, try the data file
1636 if (runData->GetT0Bin(backwardHistoNo[i]) > 0.0) {
1637 fT0s[2*i+1] = runData->GetT0Bin(backwardHistoNo[i]);
1638 fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1639 }
1640 }
1641
1642 // 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
1643 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1644 if (fT0s[2*i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1645 fT0s[2*i] = runData->GetT0BinEstimated(forwardHistoNo[i]);
1646 fRunInfo->SetT0Bin(fT0s[2*i], 2*i);
1647
1648 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1649 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1650 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(forwardHistoNo[i]);
1651 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1652 std::cerr << std::endl;
1653 }
1654 }
1655 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1656 if (fT0s[2*i+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1657 fT0s[2*i+1] = runData->GetT0BinEstimated(backwardHistoNo[i]);
1658 fRunInfo->SetT0Bin(fT0s[2*i+1], 2*i+1);
1659
1660 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1661 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
1662 std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[i]);
1663 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1664 std::cerr << std::endl;
1665 }
1666 }
1667
1668 // check if t0 is within proper bounds
1669 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1670 if ((fT0s[2*i] < 0) || (fT0s[2*i] > static_cast<Int_t>(runData->GetDataBin(forwardHistoNo[i])->size()))) {
1671 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperT0(): **ERROR** t0 data bin (" << fT0s[2*i] << ") doesn't make any sense!";
1672 std::cerr << std::endl << ">> forwardHistoNo " << forwardHistoNo[i];
1673 std::cerr << std::endl;
1674 return false;
1675 }
1676 }
1677 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1678 if ((fT0s[2*i+1] < 0) || (fT0s[2*i+1] > static_cast<Int_t>(runData->GetDataBin(backwardHistoNo[i])->size()))) {
1679 std::cerr << std::endl << ">> PRunAsymmetryBNMR::PrepareData(): **ERROR** t0 data bin (" << fT0s[2*i+1] << ") doesn't make any sense!";
1680 std::cerr << std::endl << ">> backwardHistoNo " << backwardHistoNo[i];
1681 std::cerr << std::endl;
1682 return false;
1683 }
1684 }
1685
1686 // check if addrun's are present, and if yes add the necessary t0's
1687 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1688 PRawRunData *addRunData;
1689 fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1690 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1691 // get run to be added to the main one
1692 addRunData = fRawData->GetRunData(*(fRunInfo->GetRunName(i)));
1693 if (addRunData == 0) { // couldn't get run
1694 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1695 std::cerr << std::endl;
1696 return false;
1697 }
1698
1699 // feed all T0's
1700 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1701 fAddT0s[i-1].clear();
1702 fAddT0s[i-1].resize(2*forwardHistoNo.size());
1703 for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1704 fAddT0s[i-1][j] = -1.0;
1705 }
1706
1707 // fill in the T0's from the msr-file (if present)
1708 for (Int_t j=0; j<fRunInfo->GetAddT0BinSize(i-1); j++) {
1709 fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1, j);
1710 }
1711
1712 // fill in the T0's from the data file, if not already present in the msr-file
1713 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1714 if (fAddT0s[i-1][2*j] == -1.0) // i.e. not present in the msr-file, try the data file
1715 if (addRunData->GetT0Bin(forwardHistoNo[j]) > 0.0) {
1716 fAddT0s[i-1][2*j] = addRunData->GetT0Bin(forwardHistoNo[j]);
1717 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1718 }
1719 }
1720 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1721 if (fAddT0s[i-1][2*j+1] == -1.0) // i.e. not present in the msr-file, try the data file
1722 if (addRunData->GetT0Bin(backwardHistoNo[j]) > 0.0) {
1723 fAddT0s[i-1][2*j+1] = addRunData->GetT0Bin(backwardHistoNo[j]);
1724 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1725 }
1726 }
1727
1728 // 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
1729 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1730 if (fAddT0s[i-1][2*j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1731 fAddT0s[i-1][2*j] = addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1732 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j], i-1, 2*j);
1733
1734 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1735 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1736 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(forwardHistoNo[j]);
1737 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1738 std::cerr << std::endl;
1739 }
1740 }
1741 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1742 if (fAddT0s[i-1][2*j+1] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1743 fAddT0s[i-1][2*j+1] = addRunData->GetT0BinEstimated(backwardHistoNo[j]);
1744 fRunInfo->SetAddT0Bin(fAddT0s[i-1][2*j+1], i-1, 2*j+1);
1745
1746 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1747 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1748 std::cerr << std::endl << ">> will try the estimated one: backward t0 = " << runData->GetT0BinEstimated(backwardHistoNo[j]);
1749 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1750 std::cerr << std::endl;
1751 }
1752 }
1753 }
1754 }
1755
1756 return true;
1757}
1758
1759//--------------------------------------------------------------------------
1760// GetProperDataRange (private)
1761//--------------------------------------------------------------------------
1775Bool_t PRunAsymmetryBNMR::GetProperDataRange(PRawRunData* runData, UInt_t histoNo[2])
1776{
1777 // first get start/end data
1778 Int_t start[2] = {fRunInfo->GetDataRange(0), fRunInfo->GetDataRange(2)};
1779 Int_t end[2] = {fRunInfo->GetDataRange(1), fRunInfo->GetDataRange(3)};
1780 // check if data range has been provided in the RUN block. If not, try the GLOBAL block
1781 if (start[0] == -1) {
1782 start[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1783 }
1784 if (start[1] == -1) {
1785 start[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(2);
1786 }
1787 if (end[0] == -1) {
1788 end[0] = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1789 }
1790 if (end[1] == -1) {
1791 end[1] = fMsrInfo->GetMsrGlobal()->GetDataRange(3);
1792 }
1793
1794 Double_t t0[2] = {fT0s[0], fT0s[1]};
1795 Int_t offset = static_cast<Int_t>((10.0e-3/fTimeResolution)); // needed in case first good bin is not given, default = 10ns
1796
1797 // check if data range has been provided, and if not try to estimate them
1798 if (start[0] < 0) {
1799 start[0] = static_cast<Int_t>(t0[0])+offset;
1800 fRunInfo->SetDataRange(start[0], 0);
1801 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[0] << ".";
1802 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS DOES MAKE ANY SENSE.";
1803 std::cerr << std::endl;
1804 }
1805 if (start[1] < 0) {
1806 start[1] = static_cast<Int_t>(t0[1])+offset;
1807 fRunInfo->SetDataRange(start[1], 2);
1808 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start[1] << ".";
1809 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS DOES MAKE ANY SENSE.";
1810 std::cerr << std::endl;
1811 }
1812 if (end[0] < 0) {
1813 end[0] = runData->GetDataBin(histoNo[0])->size();
1814 fRunInfo->SetDataRange(end[0], 1);
1815 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **WARNING** data range (forward) was not provided, will try data range end = " << end[0] << ".";
1816 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS DOES MAKE ANY SENSE.";
1817 std::cerr << std::endl;
1818 }
1819 if (end[1] < 0) {
1820 end[1] = runData->GetDataBin(histoNo[1])->size();
1821 fRunInfo->SetDataRange(end[1], 3);
1822 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **WARNING** data range (backward) was not provided, will try data range end = " << end[1] << ".";
1823 std::cerr << std::endl << ">> NO GUARANTEE THAT THIS DOES MAKE ANY SENSE.";
1824 std::cerr << std::endl;
1825 }
1826
1827 // check if start, end, and t0 make any sense
1828 // 1st check if start and end are in proper order
1829 for (UInt_t i=0; i<2; i++) {
1830 if (end[i] < start[i]) { // need to swap them
1831 Int_t keep = end[i];
1832 end[i] = start[i];
1833 start[i] = keep;
1834 }
1835 // 2nd check if start is within proper bounds
1836 if ((start[i] < 0) || (start[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1837 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **ERROR** start data bin doesn't make any sense!";
1838 std::cerr << std::endl;
1839 return false;
1840 }
1841 // 3rd check if end is within proper bounds
1842 if (end[i] < 0) {
1843 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **ERROR** end data bin (" << end[i] << ") doesn't make any sense!";
1844 std::cerr << std::endl;
1845 return false;
1846 }
1847 if (end[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())) {
1848 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **WARNING** end data bin (" << end[i] << ") > histo length (" << static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()) << ").";
1849 std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
1850 std::cerr << std::endl;
1851 end[i] = static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size())-1;
1852 }
1853 // 4th check if t0 is within proper bounds
1854 if ((t0[i] < 0) || (t0[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
1855 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange(): **ERROR** t0 data bin doesn't make any sense!";
1856 std::cerr << std::endl;
1857 return false;
1858 }
1859 }
1860
1861 // check that start-t0 is the same for forward as for backward, otherwise take max(start[i]-t0[i])
1862 if (fabs(static_cast<Double_t>(start[0])-t0[0]) > fabs(static_cast<Double_t>(start[1])-t0[1])){
1863 start[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(start[0]) - t0[0]);
1864 end[1] = static_cast<Int_t>(t0[1] + static_cast<Double_t>(end[0]) - t0[0]);
1865 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange **WARNING** needed to shift backward data range.";
1866 std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(2) << ", " << fRunInfo->GetDataRange(3);
1867 std::cerr << std::endl << ">> used : " << start[1] << ", " << end[1];
1868 std::cerr << std::endl;
1869 }
1870 if (fabs(static_cast<Double_t>(start[0])-t0[0]) < fabs(static_cast<Double_t>(start[1])-t0[1])){
1871 start[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(start[1]) - t0[1]);
1872 end[0] = static_cast<Int_t>(t0[0] + static_cast<Double_t>(end[1]) - t0[1]);
1873 std::cerr << std::endl << ">> PRunAsymmetryBNMR::GetProperDataRange **WARNING** needed to shift forward data range.";
1874 std::cerr << std::endl << ">> given: " << fRunInfo->GetDataRange(0) << ", " << fRunInfo->GetDataRange(1);
1875 std::cerr << std::endl << ">> used : " << start[0] << ", " << end[0];
1876 std::cerr << std::endl;
1877 }
1878
1879 // keep good bins for potential latter use
1880 fGoodBins[0] = start[0];
1881 fGoodBins[1] = end[0];
1882 fGoodBins[2] = start[1];
1883 fGoodBins[3] = end[1];
1884
1885 return true;
1886}
1887
1888//--------------------------------------------------------------------------
1889// GetProperFitRange (private)
1890//--------------------------------------------------------------------------
1908{
1909 // set fit start/end time; first check RUN Block
1910 fFitStartTime = fRunInfo->GetFitRange(0);
1911 fFitEndTime = fRunInfo->GetFitRange(1);
1912 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1913 if (fRunInfo->IsFitRangeInBin()) {
1914 fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1915 fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1916 // write these times back into the data structure. This way it is available when writting the log-file
1917 fRunInfo->SetFitRange(fFitStartTime, 0);
1918 fRunInfo->SetFitRange(fFitEndTime, 1);
1919 }
1920 if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
1921 fFitStartTime = globalBlock->GetFitRange(0);
1922 fFitEndTime = globalBlock->GetFitRange(1);
1923 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1924 if (globalBlock->IsFitRangeInBin()) {
1925 fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1926 fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1927 // write these times back into the data structure. This way it is available when writting the log-file
1928 globalBlock->SetFitRange(fFitStartTime, 0);
1929 globalBlock->SetFitRange(fFitEndTime, 1);
1930 }
1931 }
1933 fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
1934 fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
1935 std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
1936 std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
1937 }
1938}
1939
1940
1941//--------------------------------------------------------------------------
1942// EstimateAlpha (private)
1943//--------------------------------------------------------------------------
1958{
1959
1960 Double_t SumF[2] = {0.0, 0.0};
1961 Double_t SumB[2] = {0.0, 0.0};
1962 Double_t alpha = 1;
1963
1964
1965 // forward
1966 for (Int_t i=0; i<fForwardp.size(); i++) {
1967 SumF[0] += fForwardp[i];
1968 SumF[1] += fForwardm[i];
1969 }
1970
1971 // backward
1972 for (Int_t i=0; i<fBackwardp.size(); i++) {
1973 SumB[0] += fBackwardp[i];
1974 SumB[1] += fBackwardm[i];
1975 }
1976
1977 // Spit out estimated alpha value from total counts (Bp+Bm)/(Fp+Fm)
1978 if ( (SumF[0]+SumF[1]) == 0) {
1979 alpha = 1;
1980 } else {
1981 alpha = (SumB[0]+SumB[1])/(SumF[0]+SumF[1]);
1982 }
1983 std::cout << std::endl << ">> PRunAsymmetryBNMR::EstimateAlpha(): alpha estimate=" << alpha << std::endl;
1984
1985 return alpha;
1986}
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
#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
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &forwardHisto, PUIntVector &backwardHistoNo)
Retrieves proper t0 values for all histograms.
UInt_t fAlphaBetaTag
Tag indicating α/β configuration: 1=both unity, 2=α free/β unity, 3=α unity/β free,...
Bool_t SubtractEstimatedBkg()
Estimates and subtracts background from histograms.
Bool_t SubtractFixBkg()
Subtracts fixed background from histograms.
virtual ~PRunAsymmetryBNMR()
Destructor.
virtual UInt_t GetNoOfFitBins()
Returns the number of bins used in the fit.
virtual Bool_t PrepareViewData(PRawRunData *runData, UInt_t histoNo[2])
Prepares data for viewing/plotting.
UInt_t fNoOfFitBins
Number of bins included in the fit.
PDoubleVector fBackwardmErr
Negative helicity backward histogram errors.
virtual Double_t EstimateAlpha()
Estimates α parameter from data.
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
Calculates expected chi-square (for statistical analysis).
virtual Bool_t PrepareFitData()
Prepares data specifically for fitting.
PDoubleVector fForwardm
Negative helicity forward histogram data.
virtual Bool_t PrepareData()
Prepares all data for fitting or viewing.
PDoubleVector fForwardpErr
Positive helicity forward histogram errors.
Bool_t fTheoAsData
If true, theory calculated only at data points; if false, extra points for nicer Fourier transforms.
Int_t fGoodBins[4]
Good bin boundaries: [0]=forward first, [1]=forward last, [2]=backward first, [3]=backward last.
Int_t fStartTimeBin
First bin index for fitting.
virtual void CalcNoOfFitBins()
Calculates the number of bins to be fitted.
PDoubleVector fBackwardpErr
Positive helicity backward histogram errors.
virtual void CalcTheory()
Calculates theoretical asymmetry function.
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
Calculates maximum likelihood estimator.
PDoubleVector fForwardmErr
Negative helicity forward histogram errors.
PDoubleVector fBackwardp
Positive helicity backward histogram data.
virtual Double_t CalcChiSquare(const std::vector< Double_t > &par)
Calculates chi-square for the current parameter set.
virtual void SetFitRangeBin(const TString fitRange)
Sets the fit range in bins.
PDoubleVector fBackwardm
Negative helicity backward histogram data.
Int_t fEndTimeBin
Last bin index for fitting.
PRunAsymmetryBNMR()
Default constructor.
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
Determines the proper fit range from global block.
PDoubleVector fForwardp
Positive helicity forward histogram data.
virtual Bool_t GetProperDataRange(PRawRunData *runData, UInt_t histoNo[2])
Retrieves proper data range for histograms.
Int_t fPacking
Bin packing factor from RUN or GLOBAL block.
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
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.
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