45#include <TObjString.h>
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;
130 Bool_t alphaFixedToOne =
false;
131 Bool_t alphaSetDefault =
false;
132 if (
fRunInfo->GetAlphaParamNo() == -1) {
137 alphaSetDefault = true;
138 }
else if ((
fRunInfo->GetAlphaParamNo() < 0) || (
fRunInfo->GetAlphaParamNo() > (Int_t)param->size())) {
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;
145 if (((*param)[fRunInfo->GetAlphaParamNo()-1].fStep == 0.0) &&
146 ((*param)[fRunInfo->GetAlphaParamNo()-1].fValue == 1.0))
147 alphaFixedToOne = true;
151 Bool_t betaFixedToOne =
false;
153 betaFixedToOne = true;
154 }
else if ((
fRunInfo->GetBetaParamNo() < 0) || (
fRunInfo->GetBetaParamNo() > (Int_t)param->size())) {
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;
161 if (((*param)[fRunInfo->GetBetaParamNo()-1].fStep == 0.0) &&
162 ((*param)[fRunInfo->GetBetaParamNo()-1].fValue == 1.0))
163 betaFixedToOne = true;
167 if (alphaFixedToOne && betaFixedToOne)
169 else if (!alphaFixedToOne && betaFixedToOne & !alphaSetDefault)
171 else if (alphaFixedToOne && !betaFixedToOne)
173 else if (!alphaFixedToOne && betaFixedToOne & alphaSetDefault)
175 else if (!alphaFixedToOne && !betaFixedToOne & alphaSetDefault)
228 Double_t chisq = 0.0;
230 Double_t asymFcnValue = 0.0;
234 for (Int_t i=0; i<
fMsrInfo->GetNoOfFuncs(); i++) {
239 Double_t time(1.0),alphaest;
243 alphaest =
fRunInfo->GetEstimatedAlpha();
251 a = par[
fRunInfo->GetAlphaParamNo()-1];
263 b = par[
fRunInfo->GetBetaParamNo()-1];
273 a = par[
fRunInfo->GetAlphaParamNo()-1];
281 b = par[
fRunInfo->GetBetaParamNo()-1];
296 b = par[
fRunInfo->GetBetaParamNo()-1];
320 #pragma omp parallel for default(shared) private(i,time,diff,asymFcnValue,f) schedule(dynamic,chunk) reduction(+:chisq)
323 time =
fData.GetDataTimeStart() +
static_cast<Double_t
>(i)*
fData.GetDataTimeStep();
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));
364 std::cout << std::endl <<
"PRunAsymmetryBNMR::CalcMaxLikelihood(): not implemented yet ..." << std::endl;
409 TObjArray *tok =
nullptr;
410 TObjString *ostr =
nullptr;
415 tok = fitRange.Tokenize(
" \t");
417 if (tok->GetEntries() == 3) {
419 ostr =
dynamic_cast<TObjString*
>(tok->At(1));
420 str = ostr->GetString();
422 idx = str.First(
"+");
424 str.Remove(0, idx+1);
431 ostr =
dynamic_cast<TObjString*
>(tok->At(2));
432 str = ostr->GetString();
434 idx = str.First(
"-");
436 str.Remove(0, idx+1);
441 }
else if ((tok->GetEntries() > 3) && (tok->GetEntries() % 2 == 1)) {
442 Int_t pos = 2*(
fRunNo+1)-1;
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;
449 ostr =
static_cast<TObjString*
>(tok->At(pos));
450 str = ostr->GetString();
452 idx = str.First(
"+");
454 str.Remove(0, idx+1);
461 ostr =
static_cast<TObjString*
>(tok->At(pos+1));
462 str = ostr->GetString();
464 idx = str.First(
"-");
466 str.Remove(0, idx+1);
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;
530 std::vector<Double_t> par;
532 for (UInt_t i=0; i<paramList->size(); i++)
533 par.push_back((*paramList)[i].fValue);
536 for (Int_t i=0; i<
fMsrInfo->GetNoOfFuncs(); i++) {
541 Double_t asymFcnValue = 0.0;
542 Double_t a, b, f, alphaest;
546 alphaest =
fRunInfo->GetEstimatedAlpha();
548 for (UInt_t i=0; i<
fData.GetValue()->size(); i++) {
549 time =
fData.GetDataTimeStart() +
static_cast<Double_t
>(i)*
fData.GetDataTimeStep();
556 a = par[
fRunInfo->GetAlphaParamNo()-1];
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));
568 b = par[
fRunInfo->GetBetaParamNo()-1];
576 asymFcnValue = f*(b+1.0)/(2.0-f*(b-1.0))-f*(b+1.0)/(2.0+f*(b-1.0));
580 a = par[
fRunInfo->GetAlphaParamNo()-1];
588 b = par[
fRunInfo->GetBetaParamNo()-1];
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));
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));
606 b = par[
fRunInfo->GetBetaParamNo()-1];
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));
620 fData.AppendTheoryValue(asymFcnValue);
663 std::cerr << std::endl <<
">> PRunAsymmetryBNMR::PrepareData(): **ERROR** Couldn't get run " <<
fRunInfo->GetRunName()->Data() <<
"!";
664 std::cerr << std::endl;
681 for (UInt_t i=0; i<
fRunInfo->GetForwardHistoNoSize(); i++) {
682 forwardHistoNo.push_back(
fRunInfo->GetForwardHistoNo(i));
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;
690 forwardHistoNo.clear();
691 backwardHistoNo.clear();
695 for (UInt_t i=0; i<
fRunInfo->GetBackwardHistoNoSize(); i++) {
696 backwardHistoNo.push_back(
fRunInfo->GetBackwardHistoNo(i));
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;
704 forwardHistoNo.clear();
705 backwardHistoNo.clear();
724 std::cout.precision(10);
725 std::cout << std::endl <<
">> PRunAsymmetryBNMR::PrepareData(): time resolution=" << std::fixed << runData->
GetTimeResolution() <<
"(ms)" << std::endl;
728 if (!
GetProperT0(runData, globalBlock, forwardHistoNo, backwardHistoNo)) {
733 std::vector<PDoubleVector> forward, backward;
734 forward.resize(forwardHistoNo.size());
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]);
739 backward.resize(backwardHistoNo.size());
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]);
747 if (
fRunInfo->GetRunNameSize() > 1) {
749 for (UInt_t i=1; i<
fRunInfo->GetRunNameSize(); i++) {
752 if (addRunData ==
nullptr) {
753 std::cerr << std::endl <<
">> PRunAsymmetryBNMR::PrepareData(): **ERROR** Couldn't get addrun " <<
fRunInfo->GetRunName(i)->Data() <<
"!";
754 std::cerr << std::endl;
760 for (UInt_t k=0; k<forwardHistoNo.size(); k++) {
761 addRunSize = addRunData->
GetDataBin(forwardHistoNo[k])->size();
762 for (UInt_t j=0; j<addRunData->
GetDataBin(forwardHistoNo[k])->size(); j++) {
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]));
772 for (UInt_t k=0; k<backwardHistoNo.size(); k++) {
773 addRunSize = addRunData->
GetDataBin(backwardHistoNo[k])->size();
774 for (UInt_t j=0; j<addRunData->
GetDataBin(backwardHistoNo[k])->size(); j++) {
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]));
788 for (UInt_t i=0; i<
fForwardp.size(); i++) {
802 if (
fRunInfo->GetBkgRange(0) >= 0) {
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;
824 UInt_t histoNo[2] = {forwardHistoNo[0], backwardHistoNo[0]};
849 forwardHistoNo.clear();
850 backwardHistoNo.clear();
879 for (UInt_t i=0; i<
fForwardp.size(); i++) {
938 Double_t beamPeriod = 0.0;
941 if (
fRunInfo->GetInstitute()->Contains(
"psi"))
943 else if (
fRunInfo->GetInstitute()->Contains(
"ral"))
945 else if (
fRunInfo->GetInstitute()->Contains(
"triumf"))
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];
963 for (UInt_t i=0; i<2; i++) {
964 if (beamPeriod != 0.0) {
966 UInt_t fullCycles =
static_cast<UInt_t
>(timeBkg/beamPeriod);
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);
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] <<
").";
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] <<
").";
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};
1000 for (UInt_t i=start[0]; i<=end[0]; i++) {
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];
1012 for (UInt_t i=start[1]; i<=end[1]; i++) {
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];
1024 Double_t errVal = 0.0;
1025 for (UInt_t i=0; i<
fForwardp.size(); i++) {
1027 errVal = TMath::Sqrt(
fForwardp[i]+errBkgp[0]*errBkgp[0]);
1032 errVal = TMath::Sqrt(
fBackwardp[i]+errBkgp[1]*errBkgp[1]);
1037 errVal = TMath::Sqrt(
fForwardm[i]+errBkgm[0]*errBkgm[0]);
1042 errVal = TMath::Sqrt(
fBackwardm[i]+errBkgm[1]*errBkgm[1]);
1049 for (UInt_t i=0; i<
fForwardp.size(); i++) {
1056 fRunInfo->SetBkgEstimated(bkgp[0], 0);
1057 fRunInfo->SetBkgEstimated(bkgp[1], 1);
1058 fRunInfo->SetBkgEstimated(bkgm[0], 3);
1059 fRunInfo->SetBkgEstimated(bkgm[1], 4);
1064 fRunInfo->SetEstimatedAlpha(alpha);
1094 Double_t valuep = 0.0;
1095 Double_t errorp = 0.0;
1096 Double_t valuem = 0.0;
1097 Double_t errorm = 0.0;
1114 if (valuep == 0.0) {
1119 if (valuem == 0.0) {
1151 if (valuep == 0.0) {
1156 if (valuem == 0.0) {
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();
1181 Double_t asym,error;
1182 Double_t fp, bp, efp, ebp;
1183 Double_t fm, bm, efm, ebm;
1188 for (UInt_t i=0; i<noOfBins; i++) {
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);
1200 asym = (fp-bp) / (fp+bp) - (fm-bm) / (fm+bm);
1203 fData.AppendValue(asym);
1206 errorp = 2.0/((fp+bp)*(fp+bp))*TMath::Sqrt(bp*bp*efp*efp+ebp*ebp*fp*fp);
1210 errorm = 2.0/((fm+bm)*(fm+bm))*TMath::Sqrt(bm*bm*efm*efm+ebm*ebm*fm*fm);
1214 error = TMath::Sqrt(errorp*errorp+errorm*errorm);
1215 fData.AppendErrorValue(error);
1250 if (
fMsrInfo->GetMsrPlotList()->at(0).fViewPacking > 0) {
1251 packing =
fMsrInfo->GetMsrPlotList()->at(0).fViewPacking;
1255 std::vector<Double_t> par;
1257 for (UInt_t i=0; i<paramList->size(); i++)
1258 par.push_back((*paramList)[i].fValue);
1266 Int_t t0[2] = {
static_cast<Int_t
>(
fT0s[0]),
static_cast<Int_t
>(
fT0s[1])};
1270 if (start[0]-t0[0] != start[1]-t0[1]) {
1271 if (abs(start[0]-t0[0]) > abs(start[1]-t0[1])) {
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;
1277 fgb[0] = t0[0] + start[1]-t0[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;
1287 Int_t val = fgb[0]-packing*(fgb[0]/packing);
1289 if (fgb[1] - fgb[0] < 0)
1291 }
while (val + fgb[1] - fgb[0] < 0);
1294 start[1] = val + fgb[1] - fgb[0];
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;
1306 for (UInt_t i=0; i<2; i++) {
1307 if (end[i] < start[i]) {
1308 Int_t keep = end[i];
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;
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;
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;
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;
1345 for (Int_t i=start[0]; i<end[0]; i++) {
1352 if (((i-start[0]) % packing == 0) && (i != start[0])) {
1359 if (valuep == 0.0) {
1364 if (valuem == 0.0) {
1382 for (Int_t i=start[1]; i<end[1]; i++) {
1389 if (((i-start[1]) % packing == 0) && (i != start[1])) {
1396 if (valuep == 0.0) {
1401 if (valuem == 0.0) {
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();
1427 Double_t fp, bp, efp, ebp, alpha, beta = 1.0;
1428 Double_t fm, bm, efm, ebm;
1431 fData.SetDataTimeStart(
fTimeResolution*(
static_cast<Double_t
>(start[0])-t0[0]+
static_cast<Double_t
>(packing-1)/2.0));
1445 alpha = par[
fRunInfo->GetAlphaParamNo()-1];
1457 beta = par[
fRunInfo->GetBetaParamNo()-1];
1467 alpha = par[
fRunInfo->GetAlphaParamNo()-1];
1475 beta = par[
fRunInfo->GetBetaParamNo()-1];
1490 beta = par[
fRunInfo->GetBetaParamNo()-1];
1502 for (UInt_t i=0; i<forwardpPacked.
GetValue()->size(); i++) {
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);
1514 asym = (alpha*fp-bp) / (alpha*beta*fp+bp) - (alpha*fm-bm) / (alpha*beta*fm+bm);
1517 fData.AppendValue(asym);
1520 errorp = 2.0/((fp+bp)*(fp+bp))*TMath::Sqrt(bp*bp*efp*efp+ebp*ebp*fp*fp);
1524 errorm = 2.0/((fm+bm)*(fm+bm))*TMath::Sqrt(bm*bm*efm*efm+ebm*ebm*fm*fm);
1527 error = TMath::Sqrt(errorp*errorp+errorm*errorm);
1528 fData.AppendErrorValue(error);
1545 for (Int_t i=0; i<
fMsrInfo->GetNoOfFuncs(); i++) {
1550 UInt_t size = runData->
GetDataBin(histoNo[0])->size();
1553 fData.SetTheoryTimeStart(
fData.GetDataTimeStart());
1555 fData.SetTheoryTimeStep(
fData.GetDataTimeStep());
1559 fData.SetTheoryTimeStep(
fData.GetDataTimeStep()/(Double_t)factor);
1563 for (UInt_t i=0; i<size; i++) {
1564 time =
fData.GetTheoryTimeStart() +
static_cast<Double_t
>(i)*
fData.GetTheoryTimeStep();
1566 if (fabs(value) > 10.0) {
1569 fData.AppendTheoryValue(value);
1606 size_t size = 2*forwardHistoNo.size();
1607 if (backwardHistoNo.size() > forwardHistoNo.size())
1608 size = 2*backwardHistoNo.size();
1610 for (UInt_t i=0; i<
fT0s.size(); i++) {
1615 for (UInt_t i=0; i<
fRunInfo->GetT0BinSize(); i++) {
1621 if (
fT0s[i] == -1) {
1627 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1628 if (
fT0s[2*i] == -1.0)
1629 if (runData->
GetT0Bin(forwardHistoNo[i]) > 0.0) {
1634 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1635 if (
fT0s[2*i+1] == -1.0)
1636 if (runData->
GetT0Bin(backwardHistoNo[i]) > 0.0) {
1643 for (UInt_t i=0; i<forwardHistoNo.size(); i++) {
1644 if (
fT0s[2*i] == -1.0) {
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;
1655 for (UInt_t i=0; i<backwardHistoNo.size(); i++) {
1656 if (
fT0s[2*i+1] == -1.0) {
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;
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;
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;
1687 if (
fRunInfo->GetRunNameSize() > 1) {
1690 for (UInt_t i=1; i<
fRunInfo->GetRunNameSize(); i++) {
1693 if (addRunData == 0) {
1694 std::cerr << std::endl <<
">> PRunAsymmetryBNMR::GetProperT0(): **ERROR** Couldn't get addrun " <<
fRunInfo->GetRunName(i)->Data() <<
"!";
1695 std::cerr << std::endl;
1702 fAddT0s[i-1].resize(2*forwardHistoNo.size());
1703 for (UInt_t j=0; j<
fAddT0s[i-1].size(); j++) {
1708 for (Int_t j=0; j<
fRunInfo->GetAddT0BinSize(i-1); j++) {
1713 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1714 if (
fAddT0s[i-1][2*j] == -1.0)
1715 if (addRunData->
GetT0Bin(forwardHistoNo[j]) > 0.0) {
1720 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1721 if (
fAddT0s[i-1][2*j+1] == -1.0)
1722 if (addRunData->
GetT0Bin(backwardHistoNo[j]) > 0.0) {
1729 for (UInt_t j=0; j<forwardHistoNo.size(); j++) {
1730 if (
fAddT0s[i-1][2*j] == -1.0) {
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;
1741 for (UInt_t j=0; j<backwardHistoNo.size(); j++) {
1742 if (
fAddT0s[i-1][2*j+1] == -1.0) {
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;
1781 if (start[0] == -1) {
1782 start[0] =
fMsrInfo->GetMsrGlobal()->GetDataRange(0);
1784 if (start[1] == -1) {
1785 start[1] =
fMsrInfo->GetMsrGlobal()->GetDataRange(2);
1788 end[0] =
fMsrInfo->GetMsrGlobal()->GetDataRange(1);
1791 end[1] =
fMsrInfo->GetMsrGlobal()->GetDataRange(3);
1794 Double_t t0[2] = {
fT0s[0],
fT0s[1]};
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;
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;
1813 end[0] = runData->
GetDataBin(histoNo[0])->size();
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;
1820 end[1] = runData->
GetDataBin(histoNo[1])->size();
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;
1829 for (UInt_t i=0; i<2; i++) {
1830 if (end[i] < start[i]) {
1831 Int_t keep = end[i];
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;
1843 std::cerr << std::endl <<
">> PRunAsymmetryBNMR::GetProperDataRange(): **ERROR** end data bin (" << end[i] <<
") doesn't make any sense!";
1844 std::cerr << std::endl;
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;
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;
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;
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;
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;
1960 Double_t SumF[2] = {0.0, 0.0};
1961 Double_t SumB[2] = {0.0, 0.0};
1966 for (Int_t i=0; i<
fForwardp.size(); i++) {
1978 if ( (SumF[0]+SumF[1]) == 0) {
1981 alpha = (SumB[0]+SumB[1])/(SumF[0]+SumF[1]);
1983 std::cout << std::endl <<
">> PRunAsymmetryBNMR::EstimateAlpha(): alpha estimate=" << alpha << std::endl;
std::vector< UInt_t > PUIntVector
#define ACCEL_PERIOD_TRIUMF
TRIUMF accelerator cycle: 43.37 ns.
@ kFit
Fitting mode - perform least-squares fit to data.
@ kView
Viewing mode - display data and theory without fitting.
#define MSR_PARAM_FUN_OFFSET
Offset added to function indices for parameter parsing.
std::vector< PMsrParamStructure > PMsrParamList
#define ACCEL_PERIOD_PSI
PSI (Paul Scherrer Institute) accelerator cycle: 19.75 ns.
#define ACCEL_PERIOD_RAL
RAL (Rutherford Appleton Lab) - pulsed beam.
virtual void SetFitRange(Double_t dval, UInt_t idx)
virtual Double_t GetFitRange(UInt_t idx)
virtual UInt_t GetT0BinSize()
virtual Bool_t IsFitRangeInBin()
virtual Int_t GetFitRangeOffset(UInt_t idx)
virtual Double_t GetT0Bin(UInt_t idx=0)
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)
virtual const Double_t GetT0Bin(const UInt_t histoNo)
virtual const Double_t GetTimeResolution()
virtual const Bool_t IsPresent(UInt_t histoNo)
virtual const UInt_t GetNoOfTemperatures()
virtual const Double_t GetField()
virtual const Double_t GetEnergy()
virtual const Double_t GetT0BinEstimated(const UInt_t histoNo)
virtual const PDoublePairVector * GetTemperature() const
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.
Bool_t fValid
Flag indicating if run object initialized successfully; false if any error occurred.
Double_t fFitEndTime
Fit range end time in microseconds (μs) relative to t0.
PDoubleVector fFuncValues
Cached values of user-defined functions from FUNCTIONS block, evaluated at current parameters.
PMsrHandler * fMsrInfo
Pointer to MSR file handler (owned externally, not deleted here)
PMetaData fMetaData
Experimental metadata extracted from data file header (magnetic field, temperature,...
std::unique_ptr< PTheory > fTheory
Theory function evaluator (smart pointer, automatically deleted)
std::vector< PDoubleVector > fAddT0s
Time-zero bin values for additional runs to be added to main run.
EPMusrHandleTag fHandleTag
Operation mode: kFit (fitting), kView (display only), kEmpty (uninitialized)
PRunData fData
Processed data container: background-corrected, packed, with theory values.
PRunDataHandler * fRawData
Pointer to raw data handler (owned externally, not deleted here)
PDoubleVector fT0s
Time-zero bin values for all histograms in this run (forward, backward, etc.)
PRunBase()
Default constructor.
Int_t fRunNo
Run number (0-based index in MSR file RUN blocks)
PMsrRunBlock * fRunInfo
Pointer to this run's RUN block settings within fMsrInfo.
Double_t fFitStartTime
Fit range start time in microseconds (μs) relative to t0.
Raw data file reader and format converter for μSR data.
virtual void AppendErrorValue(Double_t dval)
virtual void AppendValue(Double_t dval)
virtual const PDoubleVector * GetValue()
Returns pointer to data value vector (asymmetry, counts, or y-data)
virtual const PDoubleVector * GetError()
Returns pointer to data error vector (statistical uncertainties)