modernized code to C++11 and newer.

This allows to analyze the code by external code analyzers. Since a lot is adopted,
the version is changed to 1.4.3
This commit is contained in:
2019-04-16 15:34:49 +02:00
parent e6d424e900
commit 795cd75b1e
136 changed files with 6870 additions and 7085 deletions

View File

@ -8,7 +8,7 @@
***************************************************************************/
/***************************************************************************
* Copyright (C) 2007-2016 by Andreas Suter *
* Copyright (C) 2007-2019 by Andreas Suter *
* andreas.suter@psi.ch *
* *
* This program is free software; you can redistribute it and/or modify *
@ -38,7 +38,6 @@
#include <cmath>
#include <iostream>
#include <fstream>
using namespace std;
#include <TString.h>
#include <TObjArray.h>
@ -91,9 +90,9 @@ PRunSingleHisto::PRunSingleHisto(PMsrHandler *msrInfo, PRunDataHandler *rawData,
fPacking = fMsrInfo->GetMsrGlobal()->GetPacking();
}
if (fPacking == -1) { // this should NOT happen, somethin is severely wrong
cerr << endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't find any packing information!";
cerr << endl << ">> This is very bad :-(, will quit ...";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't find any packing information!";
std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
std::cerr << std::endl;
fValid = false;
return;
}
@ -107,9 +106,9 @@ PRunSingleHisto::PRunSingleHisto(PMsrHandler *msrInfo, PRunDataHandler *rawData,
fEndTimeBin = -1;
if (!PrepareData()) {
cerr << endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't prepare data for fitting!";
cerr << endl << ">> This is very bad :-(, will quit ...";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PRunSingleHisto: **SEVERE ERROR**: Couldn't prepare data for fitting!";
std::cerr << std::endl << ">> This is very bad :-(, will quit ...";
std::cerr << std::endl;
fValid = false;
}
}
@ -174,7 +173,7 @@ Double_t PRunSingleHisto::CalcChiSquare(const std::vector<Double_t>& par)
// calculate functions
for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
Int_t funcNo = fMsrInfo->GetFuncNo(i);
UInt_t funcNo = fMsrInfo->GetFuncNo(i);
fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par);
}
@ -195,7 +194,7 @@ Double_t PRunSingleHisto::CalcChiSquare(const std::vector<Double_t>& par)
#pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
#endif
for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
time = fData.GetDataTimeStart() + (Double_t)i*fData.GetDataTimeStep();
time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
diff = fData.GetValue()->at(i) -
(N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg);
chisq += diff*diff / (fData.GetError()->at(i)*fData.GetError()->at(i));
@ -280,7 +279,7 @@ Double_t PRunSingleHisto::CalcChiSquareExpected(const std::vector<Double_t>& par
#pragma omp parallel for default(shared) private(i,time,diff) schedule(dynamic,chunk) reduction(+:chisq)
#endif
for (i=fStartTimeBin; i < fEndTimeBin; ++i) {
time = fData.GetDataTimeStart() + (Double_t)i*fData.GetDataTimeStep();
time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
diff = fData.GetValue()->at(i) - theo;
chisq += diff*diff / theo;
@ -316,7 +315,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihood(const std::vector<Double_t>& par)
N0 = par[fRunInfo->GetNormParamNo()-1];
} else { // norm is a function
// get function number
UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
Int_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
// evaluate function
N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par);
}
@ -342,7 +341,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihood(const std::vector<Double_t>& par)
// calculate functions
for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
Int_t funcNo = fMsrInfo->GetFuncNo(i);
UInt_t funcNo = fMsrInfo->GetFuncNo(i);
fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par);
}
@ -371,7 +370,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihood(const std::vector<Double_t>& par)
#pragma omp parallel for default(shared) private(i,time,theo,data) schedule(dynamic,chunk) reduction(-:mllh)
#endif
for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
time = fData.GetDataTimeStart() + (Double_t)i*fData.GetDataTimeStep();
time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
// calculate theory for the given parameter set
theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
theo *= normalizer;
@ -379,7 +378,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihood(const std::vector<Double_t>& par)
data = normalizer*fData.GetValue()->at(i);
if (theo <= 0.0) {
cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << endl;
std::cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << std::endl;
continue;
}
@ -415,7 +414,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihoodExpected(const std::vector<Double_t>&
N0 = par[fRunInfo->GetNormParamNo()-1];
} else { // norm is a function
// get function number
UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
Int_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
// evaluate function
N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par);
}
@ -441,7 +440,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihoodExpected(const std::vector<Double_t>&
// calculate functions
for (Int_t i=0; i<fMsrInfo->GetNoOfFuncs(); i++) {
Int_t funcNo = fMsrInfo->GetFuncNo(i);
UInt_t funcNo = fMsrInfo->GetFuncNo(i);
fFuncValues[i] = fMsrInfo->EvalFunc(funcNo, *fRunInfo->GetMap(), par);
}
@ -470,7 +469,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihoodExpected(const std::vector<Double_t>&
#pragma omp parallel for default(shared) private(i,time,theo,data) schedule(dynamic,chunk) reduction(-:mllh)
#endif
for (i=fStartTimeBin; i<fEndTimeBin; ++i) {
time = fData.GetDataTimeStart() + (Double_t)i*fData.GetDataTimeStep();
time = fData.GetDataTimeStart() + static_cast<Double_t>(i)*fData.GetDataTimeStep();
// calculate theory for the given parameter set
theo = N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg;
theo *= normalizer;
@ -478,7 +477,7 @@ Double_t PRunSingleHisto::CalcMaxLikelihoodExpected(const std::vector<Double_t>&
data = normalizer*fData.GetValue()->at(i);
if (theo <= 0.0) {
cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << endl;
std::cerr << ">> PRunSingleHisto::CalcMaxLikelihood: **WARNING** NEGATIVE theory!!" << std::endl;
continue;
}
@ -511,7 +510,7 @@ void PRunSingleHisto::CalcTheory()
N0 = par[fRunInfo->GetNormParamNo()-1];
} else { // norm is a function
// get function number
UInt_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
Int_t funNo = fRunInfo->GetNormParamNo()-MSR_PARAM_FUN_OFFSET;
// evaluate function
N0 = fMsrInfo->EvalFunc(funNo, *fRunInfo->GetMap(), par);
}
@ -546,7 +545,7 @@ void PRunSingleHisto::CalcTheory()
Double_t resolution = fData.GetDataTimeStep();
Double_t time;
for (UInt_t i=0; i<size; i++) {
time = start + (Double_t)i*resolution;
time = start + static_cast<Double_t>(i)*resolution;
fData.AppendTheoryValue(N0*TMath::Exp(-time/tau)*(1.0+fTheory->Func(time, par, fFuncValues))+bkg);
}
@ -585,8 +584,8 @@ UInt_t PRunSingleHisto::GetNoOfFitBins()
*/
void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
{
TObjArray *tok = 0;
TObjString *ostr = 0;
TObjArray *tok = nullptr;
TObjString *ostr = nullptr;
TString str;
Ssiz_t idx = -1;
Int_t offset = 0;
@ -595,7 +594,7 @@ void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
if (tok->GetEntries() == 3) { // structure FIT_RANGE fgb+n0 lgb-n1
// handle fgb+n0 entry
ostr = (TObjString*) tok->At(1);
ostr = dynamic_cast<TObjString*>(tok->At(1));
str = ostr->GetString();
// check if there is an offset present
idx = str.First("+");
@ -607,7 +606,7 @@ void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
// handle lgb-n1 entry
ostr = (TObjString*) tok->At(2);
ostr = dynamic_cast<TObjString*>(tok->At(2));
str = ostr->GetString();
// check if there is an offset present
idx = str.First("-");
@ -621,11 +620,11 @@ void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
Int_t pos = 2*(fRunNo+1)-1;
if (pos + 1 >= tok->GetEntries()) {
cerr << endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
cerr << endl << ">> will ignore it. Sorry ..." << endl;
std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
} else {
// handle fgb+n0 entry
ostr = (TObjString*) tok->At(pos);
ostr = dynamic_cast<TObjString*>(tok->At(pos));
str = ostr->GetString();
// check if there is an offset present
idx = str.First("+");
@ -637,7 +636,7 @@ void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
fFitStartTime = (fGoodBins[0] + offset - fT0s[0]) * fTimeResolution;
// handle lgb-n1 entry
ostr = (TObjString*) tok->At(pos+1);
ostr = dynamic_cast<TObjString*>(tok->At(pos+1));
str = ostr->GetString();
// check if there is an offset present
idx = str.First("-");
@ -649,8 +648,8 @@ void PRunSingleHisto::SetFitRangeBin(const TString fitRange)
fFitEndTime = (fGoodBins[1] - offset - fT0s[0]) * fTimeResolution;
}
} else { // error
cerr << endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
cerr << endl << ">> will ignore it. Sorry ..." << endl;
std::cerr << std::endl << ">> PRunSingleHisto::SetFitRangeBin(): **ERROR** invalid FIT_RANGE command found: '" << fitRange << "'";
std::cerr << std::endl << ">> will ignore it. Sorry ..." << std::endl;
}
// clean up
@ -707,8 +706,8 @@ Bool_t PRunSingleHisto::PrepareData()
// get the proper run
PRawRunData* runData = fRawData->GetRunData(*fRunInfo->GetRunName());
if (!runData) { // couldn't get run
cerr << endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get run " << fRunInfo->GetRunName()->Data() << "!";
std::cerr << std::endl;
return false;
}
@ -718,10 +717,10 @@ Bool_t PRunSingleHisto::PrepareData()
histoNo.push_back(fRunInfo->GetForwardHistoNo(i));
if (!runData->IsPresent(histoNo[i])) {
cerr << endl << ">> PRunSingleHisto::PrepareData(): **PANIC ERROR**:";
cerr << endl << ">> histoNo found = " << histoNo[i] << ", which is NOT present in the data file!?!?";
cerr << endl << ">> Will quit :-(";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **PANIC ERROR**:";
std::cerr << std::endl << ">> histoNo found = " << histoNo[i] << ", which is NOT present in the data file!?!?";
std::cerr << std::endl << ">> Will quit :-(";
std::cerr << std::endl;
histoNo.clear();
return false;
}
@ -729,8 +728,8 @@ Bool_t PRunSingleHisto::PrepareData()
// keep the time resolution in (us)
fTimeResolution = runData->GetTimeResolution()/1.0e3;
cout.precision(10);
cout << endl << ">> PRunSingleHisto::PrepareData(): time resolution=" << fixed << runData->GetTimeResolution() << "(ns)" << endl;
std::cout.precision(10);
std::cout << std::endl << ">> PRunSingleHisto::PrepareData(): time resolution=" << std::fixed << runData->GetTimeResolution() << "(ns)" << std::endl;
// get all the proper t0's and addt0's for the current RUN block
if (!GetProperT0(runData, globalBlock, histoNo)) {
@ -738,7 +737,7 @@ Bool_t PRunSingleHisto::PrepareData()
}
// keep the histo of each group at this point (addruns handled below)
vector<PDoubleVector> forward;
std::vector<PDoubleVector> forward;
forward.resize(histoNo.size()); // resize to number of groups
for (UInt_t i=0; i<histoNo.size(); i++) {
forward[i].resize(runData->GetDataBin(histoNo[i])->size());
@ -752,9 +751,9 @@ Bool_t PRunSingleHisto::PrepareData()
// get run to be added to the main one
addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
if (addRunData == 0) { // couldn't get run
cerr << endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
cerr << endl;
if (addRunData == nullptr) { // couldn't get run
std::cerr << std::endl << ">> PRunSingleHisto::PrepareData(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
std::cerr << std::endl;
return false;
}
@ -764,8 +763,9 @@ Bool_t PRunSingleHisto::PrepareData()
addRunSize = addRunData->GetDataBin(histoNo[k])->size();
for (UInt_t j=0; j<addRunData->GetDataBin(histoNo[k])->size(); j++) { // loop over the bin indices
// make sure that the index stays in the proper range
if (((Int_t)j+(Int_t)fAddT0s[i-1][k]-(Int_t)fT0s[k] >= 0) && (j+(Int_t)fAddT0s[i-1][k]-(Int_t)fT0s[k] < addRunSize)) {
forward[k][j] += addRunData->GetDataBin(histoNo[k])->at(j+(Int_t)fAddT0s[i-1][k]-(Int_t)fT0s[k]);
if ((static_cast<Int_t>(j)+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) >= 0) &&
(j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]) < addRunSize)) {
forward[k][j] += addRunData->GetDataBin(histoNo[k])->at(j+static_cast<Int_t>(fAddT0s[i-1][k])-static_cast<Int_t>(fT0s[k]));
}
}
}
@ -782,8 +782,8 @@ Bool_t PRunSingleHisto::PrepareData()
for (UInt_t i=1; i<histoNo.size(); i++) { // loop over the groupings
for (UInt_t j=0; j<runData->GetDataBin(histoNo[i])->size(); j++) { // loop over the bin indices
// make sure that the index stays within proper range
if (((Int_t)j+fT0s[i]-fT0s[0] >= 0) && (j+fT0s[i]-fT0s[0] < runData->GetDataBin(histoNo[i])->size())) {
fForward[j] += forward[i][j+(Int_t)fT0s[i]-(Int_t)fT0s[0]];
if ((static_cast<Int_t>(j)+fT0s[i]-fT0s[0] >= 0) && (j+fT0s[i]-fT0s[0] < runData->GetDataBin(histoNo[i])->size())) {
fForward[j] += forward[i][j+static_cast<Int_t>(fT0s[i])-static_cast<Int_t>(fT0s[0])];
}
}
}
@ -854,10 +854,10 @@ Bool_t PRunSingleHisto::PrepareFitData(PRawRunData* runData, const UInt_t histoN
} else { // no background given to do the job, try estimate
fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
cerr << endl << ">> PRunSingleHisto::PrepareFitData(): **WARNING** Neither fix background nor background bins are given!";
cerr << endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
cerr << endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PrepareFitData(): **WARNING** Neither fix background nor background bins are given!";
std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
std::cerr << std::endl;
if (!EstimateBkg(histoNo))
return false;
}
@ -869,7 +869,7 @@ Bool_t PRunSingleHisto::PrepareFitData(PRawRunData* runData, const UInt_t histoN
}
// everything looks fine, hence fill data set
Int_t t0 = (Int_t)fT0s[0];
Int_t t0 = static_cast<Int_t>(fT0s[0]);
Double_t value = 0.0;
Double_t normalizer = 1.0;
// in order that after rebinning the fit does not need to be redone (important for plots)
@ -878,7 +878,7 @@ Bool_t PRunSingleHisto::PrepareFitData(PRawRunData* runData, const UInt_t histoN
normalizer = fPacking * (fTimeResolution * 1.0e3); // fTimeResolution us->ns
// data start at data_start-t0
// time shifted so that packing is included correctly, i.e. t0 == t0 after packing
fData.SetDataTimeStart(fTimeResolution*((Double_t)fGoodBins[0]-(Double_t)t0+(Double_t)(fPacking-1)/2.0));
fData.SetDataTimeStart(fTimeResolution*(static_cast<Double_t>(fGoodBins[0])-static_cast<Double_t>(t0)+static_cast<Double_t>(fPacking-1)/2.0));
fData.SetDataTimeStep(fTimeResolution*fPacking);
for (Int_t i=fGoodBins[0]; i<fGoodBins[1]; i++) {
if (fPacking == 1) {
@ -941,22 +941,22 @@ Bool_t PRunSingleHisto::PrepareRawViewData(PRawRunData* runData, const UInt_t hi
if (fScaleN0AndBkg) {
dataNorm = 1.0/ (packing * (fTimeResolution * 1.0e3)); // fTimeResolution us->ns
} else if (!fScaleN0AndBkg && (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0)) {
theoryNorm = (Double_t)fMsrInfo->GetMsrPlotList()->at(0).fViewPacking/(Double_t)fPacking;
theoryNorm = static_cast<Double_t>(fMsrInfo->GetMsrPlotList()->at(0).fViewPacking)/static_cast<Double_t>(fPacking);
}
// raw data, since PMusrCanvas is doing ranging etc.
// start = the first bin which is a multiple of packing backward from first good data bin
Int_t start = fGoodBins[0] - (fGoodBins[0]/packing)*packing;
// end = last bin starting from start which is a multipl of packing and still within the data
// end = last bin starting from start which is a multiple of packing and still within the data
Int_t end = start + ((fForward.size()-start)/packing)*packing;
// check if data range has been provided, and if not try to estimate them
if (start < 0) {
Int_t offset = (Int_t)(10.0e-3/fTimeResolution);
start = ((Int_t)fT0s[0]+offset) - (((Int_t)fT0s[0]+offset)/packing)*packing;
Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
start = (static_cast<Int_t>(fT0s[0])+offset) - ((static_cast<Int_t>(fT0s[0])+offset)/packing)*packing;
end = start + ((fForward.size()-start)/packing)*packing;
cerr << endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
cerr << endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
std::cerr << std::endl;
}
// check if start, end, and t0 make any sense
// 1st check if start and end are in proper order
@ -966,24 +966,24 @@ Bool_t PRunSingleHisto::PrepareRawViewData(PRawRunData* runData, const UInt_t hi
start = keep;
}
// 2nd check if start is within proper bounds
if ((start < 0) || (start > (Int_t)fForward.size())) {
cerr << endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** start data bin doesn't make any sense!";
cerr << endl;
if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** start data bin doesn't make any sense!";
std::cerr << std::endl;
return false;
}
// 3rd check if end is within proper bounds
if ((end < 0) || (end > (Int_t)fForward.size())) {
cerr << endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** end data bin doesn't make any sense!";
cerr << endl;
if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **ERROR** end data bin doesn't make any sense!";
std::cerr << std::endl;
return false;
}
// everything looks fine, hence fill data set
Int_t t0 = (Int_t)fT0s[0];
Int_t t0 = static_cast<Int_t>(fT0s[0]);
Double_t value = 0.0;
// data start at data_start-t0
// time shifted so that packing is included correctly, i.e. t0 == t0 after packing
fData.SetDataTimeStart(fTimeResolution*((Double_t)start-(Double_t)t0+(Double_t)(packing-1)/2.0));
fData.SetDataTimeStart(fTimeResolution*(static_cast<Double_t>(start)-static_cast<Double_t>(t0)+static_cast<Double_t>(packing-1)/2.0));
fData.SetDataTimeStep(fTimeResolution*packing);
for (Int_t i=start; i<end; i++) {
@ -1039,10 +1039,10 @@ Bool_t PRunSingleHisto::PrepareRawViewData(PRawRunData* runData, const UInt_t hi
} else { // no background given to do the job, try estimate
fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
cerr << endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** Neither fix background nor background bins are given!";
cerr << endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
cerr << endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PrepareRawViewData(): **WARNING** Neither fix background nor background bins are given!";
std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
std::cerr << std::endl;
if (!EstimateBkg(histoNo))
return false;
}
@ -1065,7 +1065,7 @@ Bool_t PRunSingleHisto::PrepareRawViewData(PRawRunData* runData, const UInt_t hi
Double_t factor = 1.0;
if (fData.GetValue()->size() * 10 > fForward.size()) {
size = fData.GetValue()->size() * 10;
factor = (Double_t)fForward.size() / (Double_t)size;
factor = static_cast<Double_t>(fForward.size()) / static_cast<Double_t>(size);
}
Double_t time;
Double_t theoryValue;
@ -1132,13 +1132,13 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
if (fScaleN0AndBkg) {
dataNorm = 1.0/ (packing * (fTimeResolution * 1.0e3)); // fTimeResolution us->ns
} else if (!fScaleN0AndBkg && (fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0)) {
theoryNorm = (Double_t)fMsrInfo->GetMsrPlotList()->at(0).fViewPacking/(Double_t)fPacking;
theoryNorm = static_cast<Double_t>(fMsrInfo->GetMsrPlotList()->at(0).fViewPacking)/static_cast<Double_t>(fPacking);
}
// transform raw histo data. This is done the following way (for details see the manual):
// for the single histo fit, just the rebinned raw data are copied
// first get start data, end data, and t0
Int_t t0 = (Int_t)fT0s[0];
Int_t t0 = static_cast<Int_t>(fT0s[0]);
// start = the first bin which is a multiple of packing backward from first good data bin
Int_t start = fGoodBins[0] - (fGoodBins[0]/packing)*packing;
@ -1147,12 +1147,12 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
// check if data range has been provided, and if not try to estimate them
if (start < 0) {
Int_t offset = (Int_t)(10.0e-3/fTimeResolution);
start = ((Int_t)fT0s[0]+offset) - (((Int_t)fT0s[0]+offset)/packing)*packing;
Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
start = (static_cast<Int_t>(fT0s[0])+offset) - ((static_cast<Int_t>(fT0s[0])+offset)/packing)*packing;
end = start + ((fForward.size()-start)/packing)*packing;
cerr << endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
cerr << endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** data range was not provided, will try data range start = " << start << ".";
std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
std::cerr << std::endl;
}
// check if start, end, and t0 make any sense
@ -1163,15 +1163,15 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
start = keep;
}
// 2nd check if start is within proper bounds
if ((start < 0) || (start > (Int_t)fForward.size())) {
cerr << endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
cerr << endl;
if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** start data bin doesn't make any sense!";
std::cerr << std::endl;
return false;
}
// 3rd check if end is within proper bounds
if ((end < 0) || (end > (Int_t)fForward.size())) {
cerr << endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
cerr << endl;
if ((end < 0) || (end > static_cast<Int_t>(fForward.size()))) {
std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **ERROR** end data bin doesn't make any sense!";
std::cerr << std::endl;
return false;
}
@ -1213,10 +1213,10 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
} else { // no background given to do the job, try estimate
fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.1), 0);
fRunInfo->SetBkgRange(static_cast<Int_t>(fT0s[0]*0.6), 1);
cerr << endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** Neither fix background nor background bins are given!";
cerr << endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
cerr << endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::PrepareViewData(): **WARNING** Neither fix background nor background bins are given!";
std::cerr << std::endl << ">> Will try the following: bkg start = " << fRunInfo->GetBkgRange(0) << ", bkg end = " << fRunInfo->GetBkgRange(1);
std::cerr << std::endl << ">> NO WARRANTY THAT THIS MAKES ANY SENSE! Better check ...";
std::cerr << std::endl;
if (!EstimateBkg(histoNo))
return false;
}
@ -1235,7 +1235,7 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
Double_t time = 0.0;
// data start at data_start-t0 shifted by (pack-1)/2
fData.SetDataTimeStart(fTimeResolution*((Double_t)start-(Double_t)t0+(Double_t)(packing-1)/2.0));
fData.SetDataTimeStart(fTimeResolution*(static_cast<Double_t>(start)-static_cast<Double_t>(t0)+static_cast<Double_t>(packing-1)/2.0));
fData.SetDataTimeStep(fTimeResolution*packing);
// data is always normalized to (per nsec!!)
@ -1244,7 +1244,7 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
for (Int_t i=start; i<end; i++) {
if (((i-start) % packing == 0) && (i != start)) { // fill data
value *= dataNorm;
time = (((Double_t)i-(Double_t)(packing-1)/2.0)-t0)*fTimeResolution;
time = ((static_cast<Double_t>(i)-static_cast<Double_t>(packing-1)/2.0)-t0)*fTimeResolution;
expval = TMath::Exp(+time/tau)/N0;
fData.AppendValue(-1.0+expval*(value-bkg));
fData.AppendErrorValue(expval*TMath::Sqrt(value*dataNorm));
@ -1285,7 +1285,7 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
value = 0.0;
error = 0.0;
}
time = ((Double_t)i-t0)*fTimeResolution;
time = (static_cast<Double_t>(i)-t0)*fTimeResolution;
expval = TMath::Exp(+time/tau)/N0;
rrf_val = (-1.0+expval*(fForward[i]/(fTimeResolution*1.0e3)-bkg))*TMath::Cos(wRRF * time + phaseRRF);
value += rrf_val;
@ -1310,7 +1310,7 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
// check if a finer binning for the theory is needed
if (fData.GetValue()->size() * 10 > fForward.size()) {
size = fData.GetValue()->size() * 10;
factor = (Double_t)fForward.size() / (Double_t)size;
factor = static_cast<Double_t>(fForward.size()) / static_cast<Double_t>(size);
}
fData.SetTheoryTimeStart(fData.GetDataTimeStart());
fData.SetTheoryTimeStep(fTimeResolution*factor);
@ -1321,7 +1321,7 @@ Bool_t PRunSingleHisto::PrepareViewData(PRawRunData* runData, const UInt_t histo
}
for (UInt_t i=0; i<size; i++) {
time = fData.GetTheoryTimeStart() + (Double_t)i*fData.GetTheoryTimeStep();
time = fData.GetTheoryTimeStart() + static_cast<Double_t>(i)*fData.GetTheoryTimeStep();
theoryValue = fTheory->Func(time, par, fFuncValues);
if (wRRF != 0.0) {
theoryValue *= 2.0*TMath::Cos(wRRF * time + phaseRRF);
@ -1399,7 +1399,7 @@ Bool_t PRunSingleHisto::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globa
// fill in the T0's from the GLOBAL block section (if present)
for (UInt_t i=0; i<globalBlock->GetT0BinSize(); i++) {
if (fT0s[i] == -1) { // i.e. not given in the RUN block section
if (fT0s[i] == -1.0) { // i.e. not given in the RUN block section
fT0s[i] = globalBlock->GetT0Bin(i);
}
}
@ -1420,19 +1420,19 @@ Bool_t PRunSingleHisto::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globa
fT0s[i] = runData->GetT0BinEstimated(histoNo[i]);
fRunInfo->SetT0Bin(fT0s[i], i); // keep value for the msr-file
cerr << endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
cerr << endl << ">> run: " << fRunInfo->GetRunName()->Data();
cerr << endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(histoNo[i]);
cerr << endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName()->Data();
std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << runData->GetT0BinEstimated(histoNo[i]);
std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
std::cerr << std::endl;
}
}
// check if t0 is within proper bounds
for (UInt_t i=0; i<fRunInfo->GetForwardHistoNoSize(); i++) {
if ((fT0s[i] < 0) || (fT0s[i] > (Int_t)runData->GetDataBin(histoNo[i])->size())) {
cerr << endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** t0 data bin (" << fT0s[i] << ") doesn't make any sense!";
cerr << endl;
if ((fT0s[i] < 0.0) || (fT0s[i] > static_cast<Int_t>(runData->GetDataBin(histoNo[i])->size()))) {
std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** t0 data bin (" << fT0s[i] << ") doesn't make any sense!";
std::cerr << std::endl;
return false;
}
}
@ -1445,9 +1445,9 @@ Bool_t PRunSingleHisto::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globa
// get run to be added to the main one
addRunData = fRawData->GetRunData(*fRunInfo->GetRunName(i));
if (addRunData == 0) { // couldn't get run
cerr << endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
cerr << endl;
if (addRunData == nullptr) { // couldn't get run
std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** Couldn't get addrun " << fRunInfo->GetRunName(i)->Data() << "!";
std::cerr << std::endl;
return false;
}
@ -1478,19 +1478,19 @@ Bool_t PRunSingleHisto::GetProperT0(PRawRunData* runData, PMsrGlobalBlock *globa
fAddT0s[i-1][j] = addRunData->GetT0BinEstimated(histoNo[j]);
fRunInfo->SetAddT0Bin(fAddT0s[i-1][j], i-1, j); // keep value for the msr-file
cerr << endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
cerr << endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
cerr << endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(histoNo[j]);
cerr << endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **WARRNING** NO t0's found, neither in the run data nor in the msr-file!";
std::cerr << std::endl << ">> run: " << fRunInfo->GetRunName(i)->Data();
std::cerr << std::endl << ">> will try the estimated one: forward t0 = " << addRunData->GetT0BinEstimated(histoNo[j]);
std::cerr << std::endl << ">> NO WARRANTY THAT THIS OK!! For instance for LEM this is almost for sure rubbish!";
std::cerr << std::endl;
}
}
// check if t0 is within proper bounds
for (UInt_t j=0; j<fRunInfo->GetForwardHistoNoSize(); j++) {
if ((fAddT0s[i-1][j] < 0) || (fAddT0s[i-1][j] > (Int_t)addRunData->GetDataBin(histoNo[j])->size())) {
cerr << endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** addt0 data bin (" << fAddT0s[i-1][j] << ") doesn't make any sense!";
cerr << endl;
if ((fAddT0s[i-1][j] < 0.0) || (fAddT0s[i-1][j] > static_cast<Int_t>(addRunData->GetDataBin(histoNo[j])->size()))) {
std::cerr << std::endl << ">> PRunSingleHisto::GetProperT0(): **ERROR** addt0 data bin (" << fAddT0s[i-1][j] << ") doesn't make any sense!";
std::cerr << std::endl;
return false;
}
}
@ -1531,19 +1531,19 @@ Bool_t PRunSingleHisto::GetProperDataRange()
// check if data range has been provided, and if not try to estimate them
if (start < 0) {
Int_t offset = (Int_t)(10.0e-3/fTimeResolution);
start = (Int_t)fT0s[0]+offset;
Int_t offset = static_cast<Int_t>(10.0e-3/fTimeResolution);
start = static_cast<Int_t>(fT0s[0])+offset;
fRunInfo->SetDataRange(start, 0);
cerr << endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start << ".";
cerr << endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range start = t0+" << offset << "(=10ns) = " << start << ".";
std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
std::cerr << std::endl;
}
if (end < 0) {
end = fForward.size();
fRunInfo->SetDataRange(end, 1);
cerr << endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range end = " << end << ".";
cerr << endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** data range was not provided, will try data range end = " << end << ".";
std::cerr << std::endl << ">> NO WARRANTY THAT THIS DOES MAKE ANY SENSE.";
std::cerr << std::endl;
}
// check if start and end make any sense
@ -1554,22 +1554,22 @@ Bool_t PRunSingleHisto::GetProperDataRange()
start = keep;
}
// 2nd check if start is within proper bounds
if ((start < 0) || (start > (Int_t)fForward.size())) {
cerr << endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** start data bin (" << start << ") doesn't make any sense!";
cerr << endl;
if ((start < 0) || (start > static_cast<Int_t>(fForward.size()))) {
std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** start data bin (" << start << ") doesn't make any sense!";
std::cerr << std::endl;
return false;
}
// 3rd check if end is within proper bounds
if (end < 0) {
cerr << endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** end data bin (" << end << ") doesn't make any sense!";
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **ERROR** end data bin (" << end << ") doesn't make any sense!";
std::cerr << std::endl;
return false;
}
if (end > (Int_t)fForward.size()) {
cerr << endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** end data bin (" << end << ") > histo length (" << (Int_t)fForward.size() << ").";
cerr << endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << endl;
cerr << endl;
end = (Int_t)fForward.size()-1;
if (end > static_cast<Int_t>(fForward.size())) {
std::cerr << std::endl << ">> PRunSingleHisto::GetProperDataRange(): **WARNING** end data bin (" << end << ") > histo length (" << fForward.size() << ").";
std::cerr << std::endl << ">> Will set end = (histo length - 1). Consider to change it in the msr-file." << std::endl;
std::cerr << std::endl;
end = static_cast<Int_t>(fForward.size())-1;
}
// keep good bins for potential later use
@ -1622,8 +1622,8 @@ void PRunSingleHisto::GetProperFitRange(PMsrGlobalBlock *globalBlock)
if ((fFitStartTime == PMUSR_UNDEFINED) || (fFitEndTime == PMUSR_UNDEFINED)) {
fFitStartTime = (fGoodBins[0] - fT0s[0]) * fTimeResolution; // (fgb-t0)*dt
fFitEndTime = (fGoodBins[1] - fT0s[0]) * fTimeResolution; // (lgb-t0)*dt
cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << endl;
cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << endl;
std::cerr << ">> PRunSingleHisto::GetProperFitRange(): **WARNING** Couldn't get fit start/end time!" << std::endl;
std::cerr << ">> Will set it to fgb/lgb which given in time is: " << fFitStartTime << "..." << fFitEndTime << " (usec)" << std::endl;
}
}
@ -1646,7 +1646,7 @@ void PRunSingleHisto::EstimateN0()
assert(param);
if (paramNo > param->size()) {
cerr << endl << ">> PRunSingleHisto::EstimateN0: **ERROR** found parameter number " << paramNo << ", which is larger than the number of parameters = " << param->size() << endl;
std::cerr << std::endl << ">> PRunSingleHisto::EstimateN0: **ERROR** found parameter number " << paramNo << ", which is larger than the number of parameters = " << param->size() << std::endl;
return;
}
@ -1663,7 +1663,7 @@ void PRunSingleHisto::EstimateN0()
if ((paramNoBkg > 10000) || (paramNoBkg == -1)) { // i.e. fun or map
scaleBkg = false;
} else {
if (paramNoBkg-1 < (Int_t)param->size()) {
if (paramNoBkg-1 < static_cast<Int_t>(param->size())) {
bkg = param->at(paramNoBkg-1).fValue;
errBkg = param->at(paramNoBkg-1).fStep;
}
@ -1673,7 +1673,7 @@ void PRunSingleHisto::EstimateN0()
Double_t dt = fTimeResolution;
Double_t tau = PMUON_LIFETIME;
UInt_t t0 = (UInt_t)round(fT0s[0]);
UInt_t t0 = static_cast<UInt_t>(round(fT0s[0]));
Double_t dval = 0.0;
Double_t nom = 0.0;
Double_t denom = 0.0;
@ -1681,13 +1681,13 @@ void PRunSingleHisto::EstimateN0()
// calc nominator
for (UInt_t i=t0; i<fForward.size(); i++) {
xx = exp(-dt*(Double_t)(i-t0)/tau);
xx = exp(-dt*static_cast<Double_t>(i-t0)/tau);
nom += xx;
}
// calc denominator
for (UInt_t i=t0; i<fForward.size(); i++) {
xx = exp(-dt*(Double_t)(i-t0)/tau);
xx = exp(-dt*static_cast<Double_t>(i-t0)/tau);
dval = fForward[i];
if (dval > 0)
denom += xx*xx/dval;
@ -1707,9 +1707,9 @@ void PRunSingleHisto::EstimateN0()
errBkg *= rescale;
}
cout << ">> PRunSingleHisto::EstimateN0: found N0=" << param->at(paramNo-1).fValue << ", will set it to N0=" << N0 << endl;
std::cout << ">> PRunSingleHisto::EstimateN0: found N0=" << param->at(paramNo-1).fValue << ", will set it to N0=" << N0 << std::endl;
if (scaleBkg)
cout << ">> PRunSingleHisto::EstimateN0: found Bkg=" << param->at(paramNoBkg-1).fValue << ", will set it to Bkg=" << bkg << endl;
std::cout << ">> PRunSingleHisto::EstimateN0: found Bkg=" << param->at(paramNoBkg-1).fValue << ", will set it to Bkg=" << bkg << std::endl;
fMsrInfo->SetMsrParamValue(paramNo-1, N0);
fMsrInfo->SetMsrParamStep(paramNo-1, sqrt(fabs(N0)));
if (scaleBkg) {
@ -1745,41 +1745,41 @@ Bool_t PRunSingleHisto::EstimateBkg(UInt_t histoNo)
beamPeriod = 0.0;
// check if start and end are in proper order
UInt_t start = fRunInfo->GetBkgRange(0);
UInt_t end = fRunInfo->GetBkgRange(1);
Int_t start = fRunInfo->GetBkgRange(0);
Int_t end = fRunInfo->GetBkgRange(1);
if (end < start) {
cout << endl << "PRunSingleHisto::EstimatBkg(): end = " << end << " > start = " << start << "! Will swap them!";
UInt_t keep = end;
std::cout << std::endl << "PRunSingleHisto::EstimatBkg(): end = " << end << " > start = " << start << "! Will swap them!";
Int_t keep = end;
end = start;
start = keep;
}
// calculate proper background range
if (beamPeriod != 0.0) {
Double_t timeBkg = (Double_t)(end-start)*(fTimeResolution*fPacking); // length of the background intervall in time
UInt_t fullCycles = (UInt_t)(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
Double_t timeBkg = static_cast<Double_t>(end-start)*(fTimeResolution*fPacking); // length of the background intervall in time
UInt_t fullCycles = static_cast<UInt_t>(timeBkg/beamPeriod); // how many proton beam cylces can be placed within the proposed background intervall
// correct the end of the background intervall such that the background is as close as possible to a multiple of the proton cylce
end = start + (UInt_t) ((fullCycles*beamPeriod)/(fTimeResolution*fPacking));
cout << endl << "PRunSingleHisto::EstimatBkg(): Background " << start << ", " << end;
end = start + static_cast<UInt_t>((fullCycles*beamPeriod)/(fTimeResolution*fPacking));
std::cout << std::endl << "PRunSingleHisto::EstimatBkg(): Background " << start << ", " << end;
if (end == start)
end = fRunInfo->GetBkgRange(1);
}
// check if start is within histogram bounds
if (start >= fForward.size()) {
cerr << endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
cerr << endl << ">> histo lengths = " << fForward.size();
cerr << endl << ">> background start = " << start;
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
std::cerr << std::endl << ">> histo lengths = " << fForward.size();
std::cerr << std::endl << ">> background start = " << start;
std::cerr << std::endl;
return false;
}
// check if end is within histogram bounds
if (end >= fForward.size()) {
cerr << endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
cerr << endl << ">> histo lengths = " << fForward.size();
cerr << endl << ">> background end = " << end;
cerr << endl;
std::cerr << std::endl << ">> PRunSingleHisto::EstimatBkg(): **ERROR** background bin values out of bound!";
std::cerr << std::endl << ">> histo lengths = " << fForward.size();
std::cerr << std::endl << ">> background end = " << end;
std::cerr << std::endl;
return false;
}
@ -1821,13 +1821,13 @@ Bool_t PRunSingleHisto::IsScaleN0AndBkg()
PMsrLines *cmd = fMsrInfo->GetMsrCommands();
for (UInt_t i=0; i<cmd->size(); i++) {
if (cmd->at(i).fLine.Contains("SCALE_N0_BKG", TString::kIgnoreCase)) {
TObjArray *tokens = 0;
TObjString *ostr = 0;
TObjArray *tokens = nullptr;
TObjString *ostr = nullptr;
TString str;
tokens = cmd->at(i).fLine.Tokenize(" \t");
if (tokens->GetEntries() != 2) {
cerr << endl << ">> PRunSingleHisto::IsScaleN0AndBkg(): **WARNING** Found uncorrect 'SCALE_N0_BKG' command, will ignore it.";
cerr << endl << ">> Allowed commands: SCALE_N0_BKG TRUE | FALSE" << endl;
std::cerr << std::endl << ">> PRunSingleHisto::IsScaleN0AndBkg(): **WARNING** Found uncorrect 'SCALE_N0_BKG' command, will ignore it.";
std::cerr << std::endl << ">> Allowed commands: SCALE_N0_BKG TRUE | FALSE" << std::endl;
return willScale;
}
ostr = dynamic_cast<TObjString*>(tokens->At(1));