musrfit 1.10.0
PRunMuMinus.cpp
Go to the documentation of this file.
1/***************************************************************************
2
3 PRunMuMinus.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 <iostream>
39
40#include <TString.h>
41#include <TObjArray.h>
42#include <TObjString.h>
43
44#include "PRunMuMinus.h"
45
46//--------------------------------------------------------------------------
47// Constructor
48//--------------------------------------------------------------------------
65{
66 fNoOfFitBins = 0;
67 fPacking = -1;
68 fTheoAsData = false;
69
70 // the 2 following variables are need in case fit range is given in bins, and since
71 // the fit range can be changed in the command block, these variables need to be accessible
72 fGoodBins[0] = -1;
73 fGoodBins[1] = -1;
74
75 fStartTimeBin = -1;
76 fEndTimeBin = -1;
77
79}
80
81//--------------------------------------------------------------------------
82// Constructor
83//--------------------------------------------------------------------------
123PRunMuMinus::PRunMuMinus(PMsrHandler *msrInfo, PRunDataHandler *rawData, UInt_t runNo, EPMusrHandleTag tag, Bool_t theoAsData) :
124 PRunBase(msrInfo, rawData, runNo, tag), fTheoAsData(theoAsData)
125{
126 fNoOfFitBins = 0;
127
128 fPacking = fRunInfo->GetPacking();
129 if (fPacking == -1) { // i.e. packing is NOT given in the RUN-block, it must be given in the GLOBAL-block
130 fPacking = fMsrInfo->GetMsrGlobal()->GetPacking();
131 }
132 if (fPacking == -1) { // this should NOT happen, somethin is severely wrong
133 std::cerr << std::endl << ">> PRunMuMinus::PRunMuMinus: **SEVERE ERROR**: Couldn't find any packing information!";
134 std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
135 std::cerr << std::endl;
136 fValid = false;
137 return;
138 }
139
140 // the 2 following variables are need in case fit range is given in bins, and since
141 // the fit range can be changed in the command block, these variables need to be accessible
142 fGoodBins[0] = -1;
143 fGoodBins[1] = -1;
144
145 fStartTimeBin = -1;
146 fEndTimeBin = -1;
147
148 if (!PrepareData()) {
149 std::cerr << std::endl << ">> PRunMuMinus::PRunMuMinus: **SEVERE ERROR**: Couldn't prepare data for fitting!";
150 std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
151 std::cerr << std::endl;
152 fValid = false;
153 }
154}
155
156//--------------------------------------------------------------------------
157// Destructor
158//--------------------------------------------------------------------------
166{
167 fForward.clear();
168}
169
170//--------------------------------------------------------------------------
171// CalcChiSquare
172//--------------------------------------------------------------------------
207Double_t PRunMuMinus::CalcChiSquare(const std::vector<Double_t>& par)
208{
209 Double_t chisq = 0.0;
210 Double_t diff = 0.0;
211
212 // calculate functions
213 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
214 Int_t funcNo = fMsrInfo->GetFuncNo(i);
215 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
216 }
217
218 // calculate chi square
219 Double_t time(1.0);
220 Int_t i;
221
222 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
223 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
224 // for a given set of parameters---which should be done outside of the parallelized loop.
225 // For all other functions it means a tiny and acceptable overhead.
226 time = fTheory->Func(time, par, fFuncValues);
227
228 #ifdef HAVE_GOMP
229 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
230 if (chunk < 10)
231 chunk = 10;
232 #pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
233 #endif
234 for (i=fStartTimeBin; i < fEndTimeBin; ++i) {
235 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
236 diff = fData.GetValue()->at(i) - fTheory->Func(time, par, fFuncValues);
237 chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
238 }
239
240 return chisq;
241}
242
243//--------------------------------------------------------------------------
244// CalcChiSquareExpected (public)
245//--------------------------------------------------------------------------
275Double_t PRunMuMinus::CalcChiSquareExpected(const std::vector<Double_t>& par)
276{
277 Double_t chisq = 0.0;
278 Double_t diff = 0.0;
279 Double_t theo = 0.0;
280
281 // calculate functions
282 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
283 Int_t funcNo = fMsrInfo->GetFuncNo(i);
284 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
285 }
286
287 // calculate chi square
288 Double_t time(1.0);
289 Int_t i;
290
291 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
292 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
293 // for a given set of parameters---which should be done outside of the parallelized loop.
294 // For all other functions it means a tiny and acceptable overhead.
295 time = fTheory->Func(time, par, fFuncValues);
296
297 #ifdef HAVE_GOMP
298 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
299 if (chunk < 10)
300 chunk = 10;
301 #pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
302 #endif
303 for (i=fStartTimeBin; i < fEndTimeBin; ++i) {
304 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
305 theo = fTheory->Func(time, par, fFuncValues);
306 diff = fData.GetValue()->at(i) - theo;
307 chisq += diff*diff / theo;
308 }
309
310 return 0.0;
311}
312
313//--------------------------------------------------------------------------
314// CalcMaxLikelihood
315//--------------------------------------------------------------------------
360Double_t PRunMuMinus::CalcMaxLikelihood(const std::vector<Double_t>& par)
361{
362 Double_t mllh = 0.0; // maximum log likelihood assuming poisson distribution for the single bin
363
364 // calculate functions
365 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
366 Int_t funcNo = fMsrInfo->GetFuncNo(i);
367 fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par, fMetaData);
368 }
369
370 // calculate maximum log likelihood
371 Double_t theo;
372 Double_t data;
373 Double_t time(1.0);
374 Int_t i;
375
376 // Calculate the theory function once to ensure one function evaluation for the current set of parameters.
377 // This is needed for the LF and user functions where some non-thread-save calculations only need to be calculated once
378 // for a given set of parameters---which should be done outside of the parallelized loop.
379 // For all other functions it means a tiny and acceptable overhead.
380 time = fTheory->Func(time, par, fFuncValues);
381
382 #ifdef HAVE_GOMP
383 Int_t chunk = (fEndTimeBin - fStartTimeBin)/omp_get_num_procs();
384 if (chunk < 10)
385 chunk = 10;
386 #pragma omp parallel for default(shared) private(i,time,theo,data) schedule(dynamic,chunk) reduction(-:mllh)
387 #endif
388 for (i=fStartTimeBin; i < fEndTimeBin; ++i) {
389 time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
390 // calculate theory for the given parameter set
391 theo = fTheory->Func(time, par, fFuncValues);
392
393 data = fData.GetValue()->at(i);
394
395 if (theo <= 0.0) {
396 std::cerr << ">> PRunMuMinus::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << std::endl;
397 continue;
398 }
399
400 if (data > 1.0e-9) {
401 mllh += (theo-data) + data*log(data/theo);
402 } else {
403 mllh += (theo-data);
404 }
405 }
406
407 return 2.0*mllh;
408}
409
410//--------------------------------------------------------------------------
411// GetNoOfFitBins (public)
412//--------------------------------------------------------------------------
429{
431
432 return fNoOfFitBins;
433}
434
435//--------------------------------------------------------------------------
436// SetFitRangeBin (public)
437//--------------------------------------------------------------------------
449void PRunMuMinus::SetFitRangeBin(const TString fitRange)
450{
451 TObjArray *tok = nullptr;
452 TObjString *ostr = nullptr;
453 TString str;
454 Ssiz_t idx = -1;
455 Int_t offset = 0;
456
457 tok = fitRange.Tokenize(" \t");
458
459 if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
460 // handle fgb+n0 entry
461 ostr = dynamic_cast<TObjString*>(tok->At(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 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
471
472 // handle lgb-n1 entry
473 ostr = dynamic_cast<TObjString*>(tok->At(2));
474 str = ostr->GetString();
475 // check if there is an offset present
476 idx = str.First("-");
477 if (idx != -1) { // offset present
478 str.Remove(0, idx+1);
479 if (str.IsFloat()) // if str is a valid number, convert is to an integer
480 offset = str.Atoi();
481 }
482 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
483 } else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) { // structure FIT_RANGE fgb[+n00] lgb[-n01] [fgb[+n10] lgb[-n11] ... fgb[+nN0] lgb[-nN1]]
484 Int_t pos = 2*(fRunNo+1)-1;
485
486 if (pos + 1 >= tok->GetEntries()) {
487 std::cerr << std::endl << ">> PRunMuMinus::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
488 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
489 } else {
490 // handle fgb+n0 entry
491 ostr = dynamic_cast<TObjString*>(tok->At(pos));
492 str = ostr->GetString();
493 // check if there is an offset present
494 idx = str.First("+");
495 if (idx != -1) { // offset present
496 str.Remove(0, idx+1);
497 if (str.IsFloat()) // if str is a valid number, convert is to an integer
498 offset = str.Atoi();
499 }
500 fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
501
502 // handle lgb-n1 entry
503 ostr = dynamic_cast<TObjString*>(tok->At(pos+1));
504 str = ostr->GetString();
505 // check if there is an offset present
506 idx = str.First("-");
507 if (idx != -1) { // offset present
508 str.Remove(0, idx+1);
509 if (str.IsFloat()) // if str is a valid number, convert is to an integer
510 offset = str.Atoi();
511 }
512 fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
513 }
514 } else { // error
515 std::cerr << std::endl << ">> PRunMuMinus::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
516 std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
517 }
518
519 // clean up
520 if (tok) {
521 delete tok;
522 }
523}
524
525//--------------------------------------------------------------------------
526// CalcNoOfFitBins (public)
527//--------------------------------------------------------------------------
532{
533 // In order not having to loop over all bins and to stay consistent with the chisq method, calculate the start and end bins explicitly
534 fStartTimeBin = static_cast<Int_t>(ceil((fFitStartTime - fData.GetDataTimeStart())/fData.GetDataTimeStep()));
535 if (fStartTimeBin < 0)
536 fStartTimeBin = 0;
537 fEndTimeBin = static_cast<Int_t>(floor((fFitEndTime - fData.GetDataTimeStart())/fData.GetDataTimeStep())) + 1;
538 if (fEndTimeBin > static_cast<Int_t>(fData.GetValue()->size()))
539 fEndTimeBin = fData.GetValue()->size();
540
543 else
544 fNoOfFitBins = 0;
545}
546
547//--------------------------------------------------------------------------
548// CalcTheory
549//--------------------------------------------------------------------------
554{
555 // feed the parameter vector
556 std::vector<Double_t> par;
557 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
558 for (UInt_t i=0; i<paramList->size(); i++)
559 par.push_back((*paramList)[i].fValue);
560
561 // calculate functions
562 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
563 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
564 }
565
566 // calculate theory
567 UInt_t size = fData.GetValue()->size();
568 Double_t start = fData.GetDataTimeStart();
569 Double_t resolution = fData.GetDataTimeStep();
570 Double_t time;
571 for (UInt_t i=0; i<size; i++) {
572 time = start + static_cast<Double_t>(i)*resolution;
573 fData.AppendTheoryValue(fTheory->Func(time, par, fFuncValues));
574 }
575
576 // clean up
577 par.clear();
578}
579
580//--------------------------------------------------------------------------
581// PrepareData
582//--------------------------------------------------------------------------
597{
598 Bool_t success = true;
599
600 if (!fValid)
601 return false;
602
603 // keep the Global block info
604 PMsrGlobalBlock *globalBlock = fMsrInfo->GetMsrGlobal();
605
606 // get the proper run
607 PRawRunData* runData = fRawData->GetRunData(*fRunInfo->GetRunName());
608 if (!runData) { // couldn't get run
609 std::cerr << std::endl << ">> PRunMuMinus::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
610 std::cerr << std::endl;
611 return false;
612 }
613
614 // keep the field from the meta-data from the data-file
615 fMetaData.fField = runData->GetField();
616
617 // keep the energy from the meta-data from the data-file
618 fMetaData.fEnergy = runData->GetEnergy();
619
620 // keep the temperature(s) from the meta-data from the data-file
621 for (unsigned int i=0; i<runData->GetNoOfTemperatures(); i++)
622 fMetaData.fTemp.push_back(runData->GetTemperature(i));
623
624 // collect histogram numbers
625 PUIntVector histoNo; // histoNo = msr-file forward + redGreen_offset - 1
626 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
627 histoNo.push_back(fRunInfo->GetForwardHistoNo(i));
628
629 if (!runData->IsPresent(histoNo[i])) {
630 std::cerr << std::endl << ">> PRunMuMinus::PrepareData(): **PANIC ERROR**:";
631 std::cerr << std::endl << ">> histoNo found = " << histoNo[i] << ", which is NOT present in the data file!?!?";
632 std::cerr << std::endl << ">> Will quit :-(";
633 std::cerr << std::endl;
634 histoNo.clear();
635 return false;
636 }
637 }
638
639 // keep the time resolution in (us)
640 fTimeResolution = runData->GetTimeResolution()/1.0e3;
641 std::cout.precision(10);
642 std::cout << std::endl << ">> PRunMuMinus::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
643
644 // get all the proper t0's and addt0's for the current RUN block
645 if (!GetProperT0(runData, globalBlock, histoNo)) {
646 return false;
647 }
648
649 // keep the histo of each group at this point (addruns handled below)
650 std::vector<PDoubleVector> forward;
651 forward.resize(histoNo.size()); // resize to number of groups
652 for (UInt_t i=0; i<histoNo.size(); i++) {
653 forward[i].resize(runData->GetDataBin(histoNo[i])->size());
654 forward[i] = *runData->GetDataBin(histoNo[i]);
655 }
656
657 // check if there are runs to be added to the current one
658 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
659 PRawRunData *addRunData;
660 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
661
662 // get run to be added to the main one
663 addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
664 if (addRunData == nullptr) { // couldn't get run
665 std::cerr << std::endl << ">> PRunMuMinus::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
666 std::cerr << std::endl;
667 return false;
668 }
669
670 // add forward run
671 UInt_t addRunSize;
672 for (UInt_t k=0; k<histoNo.size(); k++) { // fill each group
673 addRunSize = addRunData->GetDataBin(histoNo[k])->size();
674 for (UInt_t j=0; j<addRunData->GetDataBin(histoNo[k])->size(); j++) { // loop over the bin indices
675 // make sure that the index stays in the proper range
676 if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) >= 0) &&
677 (j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) < addRunSize)) {
678 forward[k][j] += addRunData->GetDataBin(histoNo[k])->at(j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]));
679 }
680 }
681 }
682 }
683 }
684
685 // set forward/backward histo data of the first group
686 fForward.resize(forward[0].size());
687 for (UInt_t i=0; i<fForward.size(); i++) {
688 fForward[i] = forward[0][i];
689 }
690
691 // group histograms, add all the remaining forward histograms of the group
692 for (UInt_t i=1; i<histoNo.size(); i++) { // loop over the groupings
693 for (UInt_t j=0; j<runData->GetDataBin(histoNo[i])->size(); j++) { // loop over the bin indices
694 // make sure that the index stays within proper range
695 if ((static_cast<Int_t>(j)+fT0s[i]-fT0s[0] >= 0) && (j+fT0s[i]-fT0s[0] < runData->GetDataBin(histoNo[i])->size())) {
696 fForward[j] += forward[i][j+static_cast<Int_t>(fT0s[i])-static_cast<Int_t>(fT0s[0])];
697 }
698 }
699 }
700
701 // get the data range (fgb/lgb) for the current RUN block
702 if (!GetProperDataRange()) {
703 return false;
704 }
705
706 // get the fit range for the current RUN block
707 GetProperFitRange(globalBlock);
708
709 // do the more fit/view specific stuff
710 if (fHandleTag == kFit)
711 success = PrepareFitData(runData, histoNo[0]);
712 else if (fHandleTag == kView)
713 success = PrepareRawViewData(runData, histoNo[0]);
714 else
715 success = false;
716
717 // cleanup
718 histoNo.clear();
719
720 return success;
721}
722
723//--------------------------------------------------------------------------
724// PrepareFitData (private)
725//--------------------------------------------------------------------------
740Bool_t PRunMuMinus::PrepareFitData(PRawRunData* runData, const UInt_t histoNo)
741{
742 // transform raw histo data. This is done the following way (for details see the manual):
743 // for the single histo fit, just the rebinned raw data are copied
744
745 // fill data set
746 Int_t t0 = static_cast<Int_t>(fT0s[0]);
747 Double_t value = 0.0;
748 // data start at data_start-t0
749 // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
750 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(fGoodBins[0])-0.5) + static_cast<Double_t>(fPacking)/2.0 - static_cast<Double_t>(t0)));
751 fData.SetDataTimeStep(fTimeResolution*fPacking);
752 for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
753 if (fPacking == 1) {
754 value = fForward[i];
755 fData.AppendValue(value);
756 if (value == 0.0)
757 fData.AppendErrorValue(1.0);
758 else
759 fData.AppendErrorValue(TMath::Sqrt(value));
760 } else { // packed data, i.e. fPacking > 1
761 if (((i-fGoodBins[0]) % fPacking == 0) && (i != fGoodBins[0])) { // fill data
762 fData.AppendValue(value);
763 if (value == 0.0)
764 fData.AppendErrorValue(1.0);
765 else
766 fData.AppendErrorValue(TMath::Sqrt(value));
767 // reset values
768 value = 0.0;
769 }
770 value += fForward[i];
771 }
772 }
773
775
776 return true;
777}
778
779//--------------------------------------------------------------------------
780// PrepareRawViewData (private)
781//--------------------------------------------------------------------------
798Bool_t PRunMuMinus::PrepareRawViewData(PRawRunData* runData, const UInt_t histoNo)
799{
800 // check if view_packing is wished
801 Int_t packing = fPacking;
802 if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
803 packing = fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
804 }
805
806 // calculate necessary norms
807 Double_t theoryNorm = 1.0;
808 if (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
809 theoryNorm = static_cast<Double_t>(fMsrInfo->GetMsrPlotList()->at(0).fViewPacking)/static_cast<Double_t>(fPacking);
810 }
811
812 // raw data, since PMusrCanvas is doing ranging etc.
813 // start = the first bin which is a multiple of packing backward from first good data bin
814 Int_t start = fGoodBins[0] - (fGoodBins[0]/packing)*packing;
815 // end = last bin starting from start which is a multipl of packing and still within the data
816 Int_t end = start + ((fForward.size()-start)/packing)*packing;
817 // check if data range has been provided, and if not try to estimate them
818 if (start < 0) {
819 Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
820 start = (static_cast<Int_t>(fT0s[0])+offset) - ((static_cast<Int_t>(fT0s[0])+offset)/packing)*packing;
821 end = start + ((fForward.size()-start)/packing)*packing;
822 std::cerr << std::endl << ">> PRunMuMinus::PrepareData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
823 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
824 std::cerr << std::endl;
825 }
826 // check if start, end, and t0 make any sense
827 // 1st check if start and end are in proper order
828 if (end < start) { // need to swap them
829 Int_t keep = end;
830 end = start;
831 start = keep;
832 }
833 // 2nd check if start is within proper bounds
834 if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
835 std::cerr << std::endl << ">> PRunMuMinus::PrepareRawViewData(): **ERROR** start data bin doesn't make any sense!";
836 std::cerr << std::endl;
837 return false;
838 }
839 // 3rd check if end is within proper bounds
840 if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
841 std::cerr << std::endl << ">> PRunMuMinus::PrepareRawViewData(): **ERROR** end data bin doesn't make any sense!";
842 std::cerr << std::endl;
843 return false;
844 }
845
846 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
847 if (fRunInfo->IsFitRangeInBin()) {
848 fFitStartTime = (fRunInfo->GetDataRange(0) + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
849 fFitEndTime = (fRunInfo->GetDataRange(1) - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
850 }
851
852 // everything looks fine, hence fill data set
853 Int_t t0 = static_cast<Int_t>(fT0s[0]);
854 Double_t value = 0.0;
855 // data start time = (binStart - 0.5) + pack/2 - t0, with pack and binStart used as double
856 fData.SetDataTimeStart(fTimeResolution*((static_cast<Double_t>(start)-0.5) + static_cast<Double_t>(packing)/2.0 - static_cast<Double_t>(t0)));
857 fData.SetDataTimeStep(fTimeResolution*packing);
858
859 for (Int_t i=start; i<end; i++) {
860 if (((i-start) % packing == 0) && (i != start)) { // fill data
861 fData.AppendValue(value);
862 if (value == 0.0)
863 fData.AppendErrorValue(1.0);
864 else
865 fData.AppendErrorValue(TMath::Sqrt(value));
866 // reset values
867 value = 0.0;
868 }
869 value += fForward[i];
870 }
871
873
874 // fill theory vector for kView
875 // feed the parameter vector
876 std::vector<Double_t> par;
877 PMsrParamList *paramList = fMsrInfo->GetMsrParamList();
878 for (UInt_t i=0; i<paramList->size(); i++)
879 par.push_back((*paramList)[i].fValue);
880
881 // calculate functions
882 for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
883 fFuncValues[i] = fMsrInfo->EvalFunc(fMsrInfo->GetFuncNo(i), *fRunInfo->GetMap(), par, fMetaData);
884 }
885
886 // calculate theory
887 UInt_t size = fForward.size();
888/* //as35
889 Double_t factor = 1.0;
890 if (fData.GetValue()->size() * 10 > fForward.size()) {
891 size = fData.GetValue()->size() * 10;
892 factor = static_cast<Double_t>(fForward.size()) / static_cast<Double_t>(size);
893 }
894 Double_t time;
895 Double_t theoryValue;
896 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
897 fData.SetTheoryTimeStep(fTimeResolution*factor);
898*/ //as35
899
900 Int_t factor = 8; // 8 times more points for the theory (if fTheoAsData == false)
901 fData.SetTheoryTimeStart(fData.GetDataTimeStart());
902 if (fTheoAsData) { // calculate theory only at the data points
903 fData.SetTheoryTimeStep(fData.GetDataTimeStep());
904 } else {
905 // finer binning for the theory (8 times as many points = factor)
906 size *= factor;
907 fData.SetTheoryTimeStep(fData.GetDataTimeStep()/(Double_t)factor);
908 }
909
910 Double_t time;
911 Double_t theoryValue;
912 for (UInt_t i=0; i<size; i++) {
913 time = fData.GetTheoryTimeStart() + i*fData.GetTheoryTimeStep();
914 theoryValue = fTheory->Func(time, par, fFuncValues);
915 if (fabs(theoryValue) > 1.0e10) { // dirty hack needs to be fixed!!
916 theoryValue = 0.0;
917 }
918 fData.AppendTheoryValue(theoryNorm*theoryValue);
919 }
920
921 // clean up
922 par.clear();
923
924 return true;
925}
926
927//--------------------------------------------------------------------------
928// GetProperT0 (private)
929//--------------------------------------------------------------------------
945Bool_t PRunMuMinus::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globalBlock, PUIntVector &histoNo)
946{
947 // feed all T0's
948 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
949 fT0s.clear();
950 fT0s.resize(histoNo.size());
951 for (UInt_t i=0; i<fT0s.size(); i++) {
952 fT0s[i] = -1.0;
953 }
954
955 // fill in the T0's from the msr-file (if present)
956 for (UInt_t i=0; i<fRunInfo->GetT0BinSize(); i++) {
957 fT0s[i] = fRunInfo->GetT0Bin(i);
958 }
959
960 // fill in the T0's from the GLOBAL block section (if present)
961 for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
962 if (fT0s[i] == -1.0) { // i.e. not given in the RUN block section
963 fT0s[i] = globalBlock->GetT0Bin(i);
964 }
965 }
966
967 // fill in the T0's from the data file, if not already present in the msr-file
968 for (UInt_t i=0; i<histoNo.size(); i++) {
969 if (fT0s[i] == -1.0) { // i.e. not present in the msr-file, try the data file
970 if (runData->GetT0Bin(histoNo[i]) > 0.0) {
971 fT0s[i] = runData->GetT0Bin(histoNo[i]);
972 fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
973 }
974 }
975 }
976
977 // 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
978 for (UInt_t i=0; i<histoNo.size(); i++) {
979 if (fT0s[i] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
980 fT0s[i] = runData->GetT0BinEstimated(histoNo[i]);
981 fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
982
983 std::cerr << std::endl << ">> PRunMuMinus::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
984 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
985 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(histoNo[i]);
986 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
987 std::cerr << std::endl;
988 }
989 }
990
991 // check if t0 is within proper bounds
992 for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
993 if ((fT0s[i] < 0) || (fT0s[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
994 std::cerr << std::endl << ">> PRunMuMinus::GetProperT0(): **ERROR** t0 data bin (" << fT0s[i] << ") doesn't make any sense!";
995 std::cerr << std::endl;
996 return false;
997 }
998 }
999
1000 // check if there are runs to be added to the current one
1001 if (fRunInfo->GetRunNameSize() > 1) { // runs to be added present
1002 PRawRunData *addRunData;
1003 fAddT0s.resize(fRunInfo->GetRunNameSize()-1); // resize to the number of addruns
1004 for (UInt_t i=1; i<fRunInfo->GetRunNameSize(); i++) {
1005
1006 // get run to be added to the main one
1007 addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
1008 if (addRunData == nullptr) { // couldn't get run
1009 std::cerr << std::endl << ">> PRunMuMinus::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
1010 std::cerr << std::endl;
1011 return false;
1012 }
1013
1014 // feed all T0's
1015 // first init T0's, T0's are stored as (forward T0, backward T0, etc.)
1016 fAddT0s[i-1].resize(histoNo.size());
1017 for (UInt_t j=0; j<fAddT0s[i-1].size(); j++) {
1018 fAddT0s[i-1][j] = -1.0;
1019 }
1020
1021 // fill in the T0's from the msr-file (if present)
1022 for (UInt_t j=0; j<fRunInfo->GetT0BinSize(); j++) {
1023 fAddT0s[i-1][j] = fRunInfo->GetAddT0Bin(i-1,j); // addRunIdx starts at 0
1024 }
1025
1026 // fill in the T0's from the data file, if not already present in the msr-file
1027 for (UInt_t j=0; j<histoNo.size(); j++) {
1028 if (fAddT0s[i-1][j] == -1.0) // i.e. not present in the msr-file, try the data file
1029 if (addRunData->GetT0Bin(histoNo[j]) > 0.0) {
1030 fAddT0s[i-1][j] = addRunData->GetT0Bin(histoNo[j]);
1031 fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1032 }
1033 }
1034
1035 // 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
1036 for (UInt_t j=0; j<histoNo.size(); j++) {
1037 if (fAddT0s[i-1][j] == -1.0) { // i.e. not present in the msr-file and data file, use the estimated T0
1038 fAddT0s[i-1][j] = addRunData->GetT0BinEstimated(histoNo[j]);
1039 fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
1040
1041 std::cerr << std::endl << ">> PRunMuMinus::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
1042 std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
1043 std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(histoNo[j]);
1044 std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
1045 std::cerr << std::endl;
1046 }
1047 }
1048
1049 // check if t0 is within proper bounds
1050 for (UInt_t j=0; j<fRunInfo->GetForwardHistoNoSize(); j++) {
1051 if ((fAddT0s[i-1][j] < 0) || (fAddT0s[i-1][j] > static_cast<Int_t>(addRunData->GetDataBin(histoNo[j])->size()))) {
1052 std::cerr << std::endl << ">> PRunMuMinus::GetProperT0(): **ERROR** addt0 data bin (" << fAddT0s[i-1][j] << ") doesn't make any sense!";
1053 std::cerr << std::endl;
1054 return false;
1055 }
1056 }
1057 }
1058 }
1059
1060 return true;
1061}
1062
1063//--------------------------------------------------------------------------
1064// GetProperDataRange (private)
1065//--------------------------------------------------------------------------
1077{
1078 // get start/end data
1079 Int_t start;
1080 Int_t end;
1081 start = fRunInfo->GetDataRange(0);
1082 end = fRunInfo->GetDataRange(1);
1083
1084 // check if data range has been given in the RUN block, if not try to get it from the GLOBAL block
1085 if (start < 0) {
1086 start = fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1087 }
1088 if (end < 0) {
1089 end = fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1090 }
1091
1092 // check if data range has been provided, and if not try to estimate them
1093 if (start < 0) {
1094 Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
1095 start = static_cast<Int_t>(fT0s[0])+offset;
1096 fRunInfo->SetDataRange(start, 0);
1097 std::cerr << std::endl << ">> PRunMuMinus::GetProperDataRange(): **WARNING** data range was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start << ".";
1098 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1099 std::cerr << std::endl;
1100 }
1101 if (end < 0) {
1102 end = fForward.size();
1103 fRunInfo->SetDataRange(end, 1);
1104 std::cerr << std::endl << ">> PRunMuMinus::GetProperDataRange(): **WARNING** data range was not provided, will try data range end = " << end << ".";
1105 std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
1106 std::cerr << std::endl;
1107 }
1108
1109 // check if start and end make any sense
1110 // 1st check if start and end are in proper order
1111 if (end < start) { // need to swap them
1112 Int_t keep = end;
1113 end = start;
1114 start = keep;
1115 }
1116 // 2nd check if start is within proper bounds
1117 if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
1118 std::cerr << std::endl << ">> PRunMuMinus::GetProperDataRange(): **ERROR** start data bin (" << start << ") doesn't make any sense!";
1119 std::cerr << std::endl;
1120 return false;
1121 }
1122 // 3rd check if end is within proper bounds
1123 if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
1124 std::cerr << std::endl << ">> PRunMuMinus::GetProperDataRange(): **ERROR** end data bin (" << end << ") doesn't make any sense!";
1125 std::cerr << std::endl;
1126 return false;
1127 }
1128
1129 // keep good bins for potential later use
1130 fGoodBins[0] = start;
1131 fGoodBins[1] = end;
1132
1133 return true;
1134}
1135
1136//--------------------------------------------------------------------------
1137// GetProperFitRange (private)
1138//--------------------------------------------------------------------------
1152{
1153 // set fit start/end time; first check RUN Block
1154 fFitStartTime = fRunInfo->GetFitRange(0);
1155 fFitEndTime = fRunInfo->GetFitRange(1);
1156 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1157 if (fRunInfo->IsFitRangeInBin()) {
1158 fFitStartTime = (fGoodBins[0] + fRunInfo->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1159 fFitEndTime = (fGoodBins[1] - fRunInfo->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1160 // write these times back into the data structure. This way it is available when writting the log-file
1161 fRunInfo->SetFitRange(fFitStartTime, 0);
1162 fRunInfo->SetFitRange(fFitEndTime, 1);
1163 }
1164 if (fFitStartTime == PMUSR_UNDEFINED) { // fit start/end NOT found in the RUN block, check GLOBAL block
1165 fFitStartTime = globalBlock->GetFitRange(0);
1166 fFitEndTime = globalBlock->GetFitRange(1);
1167 // if fit range is given in bins (and not time), the fit start/end time can be calculated at this point now
1168 if (globalBlock->IsFitRangeInBin()) {
1169 fFitStartTime = (fGoodBins[0] + globalBlock->GetFitRangeOffset(0) - fT0s[0]) * fTimeResolution; // (fgb+n0-t0)*dt
1170 fFitEndTime = (fGoodBins[1] - globalBlock->GetFitRangeOffset(1) - fT0s[0]) * fTimeResolution; // (lgb-n1-t0)*dt
1171 // write these times back into the data structure. This way it is available when writting the log-file
1172 globalBlock->SetFitRange(fFitStartTime, 0);
1173 globalBlock->SetFitRange(fFitEndTime, 1);
1174 }
1175 }
1177 fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
1178 fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
1179 std::cerr << ">> PRunMuMinus::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
1180 std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
1181 }
1182}
std::vector< UInt_t > PUIntVector
Definition PMusr.h:361
EPMusrHandleTag
Definition PMusr.h:413
@ kEmpty
No operation active.
Definition PMusr.h:414
@ 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
std::vector< PMsrParamStructure > PMsrParamList
Definition PMusr.h:1022
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 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
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 ~PRunMuMinus()
Virtual destructor cleaning up allocated resources.
Bool_t fTheoAsData
Theory calculation mode flag.
PRunMuMinus()
Default constructor creating an empty, invalid μ⁻ run object.
Int_t fPacking
Bin packing factor (REQUIRED for μ⁻).
Int_t fGoodBins[2]
Good bin markers for bin-based fit range specification.
virtual Bool_t PrepareData()
Main data preparation routine for μ⁻ fitting and viewing.
virtual Double_t CalcChiSquareExpected(const std::vector< Double_t > &par)
Calculates expected χ² based on theory predictions (statistical diagnostic).
virtual Bool_t GetProperDataRange()
Determines data range (region of valid histogram data).
virtual void SetFitRangeBin(const TString fitRange)
Sets fit range using bin-offset specification (COMMANDS block syntax).
virtual void GetProperFitRange(PMsrGlobalBlock *globalBlock)
Determines fit range from MSR file settings.
virtual UInt_t GetNoOfFitBins()
Returns the number of bins included in the fit range.
virtual void CalcTheory()
Evaluates theory function at all data points (or high-resolution grid).
PDoubleVector fForward
Forward detector histogram data (background-corrected, packed).
virtual Bool_t GetProperT0(PRawRunData *runData, PMsrGlobalBlock *globalBlock, PUIntVector &histoNo)
Determines and validates t0 values for μ⁻ histogram.
Int_t fEndTimeBin
Last bin index in fit range (exclusive: loop as i < fEndTimeBin)
Int_t fStartTimeBin
First bin index in fit range (inclusive, 0-based after packing)
virtual Bool_t PrepareFitData(PRawRunData *runData, const UInt_t histoNo)
Prepares μ⁻ histogram data for fitting.
virtual Double_t CalcChiSquare(const std::vector< Double_t > &par)
Calculates χ² between μ⁻ data and theory (least-squares fit metric).
virtual void CalcNoOfFitBins()
Calculates start/end bin indices from fit time range.
virtual Double_t CalcMaxLikelihood(const std::vector< Double_t > &par)
Calculates negative log-likelihood for Poisson statistics (low-count fit metric).
virtual Bool_t PrepareRawViewData(PRawRunData *runData, const UInt_t histoNo)
Prepares μ⁻ histogram data for viewing/plotting (minimal processing).
UInt_t fNoOfFitBins
Number of bins within fit range (between fStartTimeBin and fEndTimeBin)