57int main(
int argc,
char** argv)
65 double projGateLow = -1.;
66 double projGateHigh = -1.;
67 double projBgLow = -1.;
68 double projBgHigh = -1.;
69 double timeRandomNorm = -1.;
72 double peakHigh = -1.;
73 std::string baseName =
"AngularCorrelation";
74 std::string bgName = baseName +
"BG";
75 std::string mixedName = baseName +
"Mixed";
76 std::string inputFile;
77 std::string theoryFile;
78 std::string outputFile =
"AngularCorrelations.root";
79 std::string settingsFile;
83 parser.
option(
"h help ?", &help,
true).description(
"Show this help message.");
84 parser.
option(
"projection-low proj-low", &projGateLow,
true).description(
"Low edge of projection gate.");
85 parser.
option(
"projection-high proj-high", &projGateHigh,
true).description(
"High edge of projection gate.");
86 parser.
option(
"background-low bg-low", &projBgLow,
true).description(
"Low edge of background gate (for multiple gates use settings file).");
87 parser.
option(
"background-high bg-high", &projBgHigh,
true).description(
"High edge of background gate (for multiple gates use settings file).");
88 parser.
option(
"time-random-normalization", &timeRandomNorm,
true).description(
"Normalization factor for subtraction of time-random matrix. If negative (default) it will be calculated automatically.");
89 parser.
option(
"base-name", &baseName,
true).description(
"Base name of matrices.").default_value(
"AngularCorrelation");
90 parser.
option(
"background-name", &bgName,
true).description(
"Name of backround matrices.").default_value(
"AngularCorrelationBG");
91 parser.
option(
"mixed-name", &mixedName,
true).description(
"Name of mixed matrices.").default_value(
"AngularCorrelationMixed");
92 parser.
option(
"peak-pos peak", &peakPos,
true).description(
"Peak position (for multiple peaks use settings file).");
93 parser.
option(
"peak-low", &peakLow,
true).description(
"Low edge of peak fit.");
94 parser.
option(
"peak-high", &peakHigh,
true).description(
"High edge of peak fit.");
95 parser.
option(
"settings", &settingsFile,
true).description(
"Settings file with user settings, these do not overwrite anything provided on command line!");
96 parser.
option(
"input", &inputFile,
true).description(
"Input file with gamma-gamma matrices for each angle (coincident, time-random, and event mixed).");
97 parser.
option(
"theory", &theoryFile,
true).description(
"File with simulated z0-, z2-, and z4-graphs.");
98 parser.
option(
"output", &outputFile,
true).description(
"Name of output file, default is \"AngularCorrelations.root\".");
100 parser.
parse(argc, argv,
true);
103 std::cout << parser << std::endl;
108 if(inputFile.empty()) {
109 std::cerr <<
"Need an input file!" << std::endl;
110 std::cout << parser << std::endl;
114 TFile input(inputFile.c_str());
116 if(!input.IsOpen()) {
117 std::cerr <<
"Failed to open input file " << inputFile << std::endl;
118 std::cout << parser << std::endl;
122 auto* settings =
static_cast<TUserSettings*
>(input.Get(
"UserSettings"));
126 if(!settingsFile.empty()) {
127 if(settings ==
nullptr) {
130 settings->ReadSettings(settingsFile);
135 if(settings ==
nullptr || settings->empty()) {
136 std::cerr <<
"Failed to get user settings from input file." << std::endl;
137 std::cout << parser << std::endl;
142 if(projGateLow == -1.) { projGateLow = settings->GetDouble(
"Projection.Low"); }
143 if(projGateHigh == -1.) { projGateHigh = settings->GetDouble(
"Projection.High"); }
144 if(timeRandomNorm == -1.) { timeRandomNorm = settings->GetDouble(
"TimeRandomNormalization"); }
145 if(peakPos == -1.) { peakPos = settings->GetDouble(
"Peak.Position"); }
146 if(peakLow == -1.) { peakLow = settings->GetDouble(
"Peak.Low"); }
147 if(peakHigh == -1.) { peakHigh = settings->GetDouble(
"Peak.High"); }
148 if(baseName ==
"AngularCorrelation") { baseName = settings->GetString(
"Histograms.BaseName",
"AngularCorrelation"); }
149 if(bgName ==
"AngularCorrelationBG") { bgName = settings->GetString(
"Histograms.BackgroundName",
"AngularCorrelationBG"); }
150 if(mixedName ==
"AngularCorrelationMixed") { mixedName = settings->GetString(
"Histograms.MixedName",
"AngularCorrelationMixed"); }
151 if(theoryFile.empty()) {
153 theoryFile = settings->GetString(
"Theory",
true);
154 }
catch(std::out_of_range&) {}
156 if(outputFile ==
"AngularCorrelations.root") { outputFile = settings->GetString(
"Output",
"AngularCorrelations.root"); }
159 std::vector<std::tuple<double, double, double>> bgPeakPos;
160 std::vector<double>
bgLow;
161 std::vector<double>
bgHigh;
163 if(projBgLow == -1. || projBgHigh == -1. || projBgLow >= projBgHigh) {
164 for(
int i = 1;; ++i) {
166 auto low = settings->GetDouble(Form(
"Background.%d.Low", i),
true);
167 auto high = settings->GetDouble(Form(
"Background.%d.High", i),
true);
169 std::cout << i <<
". background gate, low edge not lower than the high edge: " << low <<
" >= " << high << std::endl;
172 bgLow.push_back(low);
174 }
catch(std::out_of_range&) {
179 bgLow.push_back(projBgLow);
180 bgHigh.push_back(projBgHigh);
184 for(
int i = 1;; ++i) {
186 auto pos = settings->GetDouble(Form(
"Background.Peak.%d.Position", i),
true);
187 if(pos <= peakLow || pos >= peakHigh) {
188 std::cout << i <<
". background peak outside of fit range: " << pos <<
" <= " << peakLow <<
" or " << pos <<
" >= " << peakHigh << std::endl;
192 auto low = settings->GetDouble(Form(
"Background.Peak.%d.Low", i), -1.);
193 auto high = settings->GetDouble(Form(
"Background.Peak.%d.Low", i), -1.);
194 bgPeakPos.emplace_back(pos, low, high);
195 }
catch(std::out_of_range&) {
203 std::vector<double> peakParameter(6, -2.);
204 std::vector<double> peakParameterLow(6, 0.);
205 std::vector<double> peakParameterHigh(6, -1.);
206 std::vector<double> bgPeakParameter(6, -2.);
207 std::vector<double> bgPeakParameterLow(6, 0.);
208 std::vector<double> bgPeakParameterHigh(6, -1.);
209 for(
size_t i = 0; i < peakParameterLow.size(); ++i) {
210 peakParameter[i] = settings->GetDouble(Form(
"Peak.Parameter.%d",
static_cast<int>(i)), -2.);
211 peakParameterLow[i] = settings->GetDouble(Form(
"Peak.Parameter.%d.Low",
static_cast<int>(i)), 0.);
212 peakParameterHigh[i] = settings->GetDouble(Form(
"Peak.Parameter.%d.High",
static_cast<int>(i)), -1.);
213 bgPeakParameter[i] = settings->GetDouble(Form(
"Background.Peak.Parameter.%d",
static_cast<int>(i)), -2.);
214 bgPeakParameterLow[i] = settings->GetDouble(Form(
"Background.Peak.Parameter.%d.Low",
static_cast<int>(i)), 0.);
215 bgPeakParameterHigh[i] = settings->GetDouble(Form(
"Background.Peak.Parameter.%d.High",
static_cast<int>(i)), -1.);
218 if(peakParameterLow[i] <= peakParameterHigh[i] && (peakParameter[i] < peakParameterLow[i] || peakParameterHigh[i] < peakParameter[i])) {
219 bool output = peakParameter[i] != -2.;
220 if(output) { std::cout <<
"Warning, " << i <<
". peak parameter (" << peakParameter[i] <<
") is out of range " << peakParameterLow[i] <<
" - " << peakParameterHigh[i] <<
", resetting it to "; }
221 peakParameter[i] = (peakParameterHigh[i] + peakParameterLow[i]) / 2.;
222 if(output) { std::cout << peakParameter[i] << std::endl; }
224 if(bgPeakParameterLow[i] <= bgPeakParameterHigh[i] && (bgPeakParameter[i] < bgPeakParameterLow[i] || bgPeakParameterHigh[i] < bgPeakParameter[i])) {
225 bool output = bgPeakParameter[i] != -2.;
226 if(output) { std::cout <<
"Warning, " << i <<
". background peak parameter (" << bgPeakParameter[i] <<
") is out of range " << bgPeakParameterLow[i] <<
" - " << bgPeakParameterHigh[i] <<
", resetting it to "; }
227 bgPeakParameter[i] = (bgPeakParameterHigh[i] + bgPeakParameterLow[i]) / 2.;
228 if(output) { std::cout << bgPeakParameter[i] << std::endl; }
231 if(i == 1 && peakParameterLow[i] == 0. && peakParameterHigh[i] == -1.) {
232 peakParameter[i] = peakPos;
233 peakParameterLow[i] = peakPos;
234 peakParameterHigh[i] = peakPos;
239 std::vector<double> backgroundParameter(4, -2.);
240 std::vector<double> backgroundParameterLow(4, 0.);
241 std::vector<double> backgroundParameterHigh(4, -1.);
242 backgroundParameter[0] = settings->GetDouble(
"Background.Offset", -2.);
243 backgroundParameterLow[0] = settings->GetDouble(
"Background.Offset.Low", 0.);
244 backgroundParameterHigh[0] = settings->GetDouble(
"Background.Offset.High", -1.);
245 backgroundParameter[1] = settings->GetDouble(
"Background.Linear", -2.);
246 backgroundParameterLow[1] = settings->GetDouble(
"Background.Linear.Low", 0.);
247 backgroundParameterHigh[1] = settings->GetDouble(
"Background.Linear.High", -1.);
248 backgroundParameter[2] = settings->GetDouble(
"Background.Quadratic", -2.);
249 backgroundParameterLow[2] = settings->GetDouble(
"Background.Quadratic.Low", 0.);
250 backgroundParameterHigh[2] = settings->GetDouble(
"Background.Quadratic.High", -1.);
251 backgroundParameter[3] = settings->GetDouble(
"Background.Xoffset", -2.);
252 backgroundParameterLow[3] = settings->GetDouble(
"Background.Xoffset.Low", 0.);
253 backgroundParameterHigh[3] = settings->GetDouble(
"Background.Xoffset.High", -1.);
254 for(
size_t i = 0; i < backgroundParameter.size(); ++i) {
255 if(backgroundParameterLow[i] <= backgroundParameterHigh[i] && (backgroundParameter[i] < backgroundParameterLow[i] || backgroundParameterHigh[i] < backgroundParameter[i])) {
256 bool output = backgroundParameter[i] != -2.;
257 if(output) { std::cout <<
"Warning, " << i <<
". background parameter (" << backgroundParameter[i] <<
") is out of range " << backgroundParameterLow[i] <<
" - " << backgroundParameterHigh[i] <<
", resetting it to "; }
258 backgroundParameter[i] = (backgroundParameterHigh[i] + backgroundParameterLow[i]) / 2.;
259 if(output) { std::cout << backgroundParameter[i] << std::endl; }
265 std::vector<int> twoJLow = settings->GetIntVector(
"TwoJ.Low");
266 std::vector<int> twoJMiddle = settings->GetIntVector(
"TwoJ.Middle");
267 std::vector<int> twoJHigh = settings->GetIntVector(
"TwoJ.High");
270 double confidenceLevel = settings->GetDouble(
"ConfidenceLevel", 1.535);
273 if(projGateLow >= projGateHigh) {
274 std::cerr <<
"Need a projection gate with a low edge that is smaller than the high edge, " << projGateLow <<
" >= " << projGateHigh << std::endl;
275 std::cout << parser << std::endl;
279 if(peakLow >= peakHigh) {
280 std::cerr <<
"Need a fit range with a low edge that is smaller than the high edge, " << peakLow <<
" >= " << peakHigh << std::endl;
281 std::cout << parser << std::endl;
285 if(peakPos >= peakHigh || peakPos <= peakLow) {
286 std::cerr <<
"Need a peak within the fit range, " << peakPos <<
" not within " << peakLow <<
" - " << peakHigh << std::endl;
287 std::cout << parser << std::endl;
291 if(timeRandomNorm <= 0.) {
292 std::cerr <<
"Need a positive normalization factor for time random subtraction" << std::endl;
293 std::cout << parser << std::endl;
300 std::cerr <<
"Background gate information missing, either no low/high edges or a mismatching amount of low and high edges: " <<
bgLow.size() <<
" low edges, and " <<
bgHigh.size() <<
" high edges" << std::endl;
301 std::cout << parser << std::endl;
306 auto* angles =
static_cast<TGriffinAngles*
>(input.Get(
"GriffinAngles"));
308 if(angles ==
nullptr) {
309 std::cerr <<
"Failed to find 'GriffinAngles' in '" << inputFile <<
"'" << std::endl;
310 std::cout << parser << std::endl;
315 auto logFileName = outputFile.substr(0, outputFile.find_last_of(
'.')) +
".log";
316 std::ofstream logFile(logFileName.c_str());
323 TFile output(outputFile.c_str(),
"recreate");
325 auto* rawAngularDistribution =
new TGraphErrors(angles->NumberOfAngles());
326 rawAngularDistribution->SetName(
"RawAngularDistribution");
327 auto* angularDistribution =
new TGraphErrors(angles->NumberOfAngles());
328 angularDistribution->SetName(
"AngularDistribution");
329 auto* mixedAngularDistribution =
new TGraphErrors(angles->NumberOfAngles());
330 mixedAngularDistribution->SetName(
"MixedAngularDistribution");
331 auto* rawChiSquares =
new TGraph(angles->NumberOfAngles());
332 rawChiSquares->SetName(
"RawChiSquares");
333 auto* mixedChiSquares =
new TGraph(angles->NumberOfAngles());
334 mixedChiSquares->SetName(
"MixedChiSquares");
339#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
340 auto* fitDir = output.mkdir(
"fits",
"Projections with fits",
true);
342 auto* fitDir = output.mkdir(
"fits",
"Projections with fits");
345 logFile <<
"# columns are:" << std::endl
346 <<
"#ID p/m centroid +- uncertainty area +- uncertainty FWHM +- uncertainty red. chi^2" << std::endl;
349 for(
int i = 0; i < angles->NumberOfAngles(); ++i) {
351 auto* prompt =
static_cast<TH2*
>(input.Get(Form(
"%s%d", baseName.c_str(), i)));
352 if(prompt ==
nullptr) {
353 std::cerr <<
"Failed to find histogram '" << Form(
"%s%d", baseName.c_str(), i) <<
"', should have " << angles->NumberOfAngles() <<
" angles in total!" << std::endl;
354 std::cout << parser << std::endl;
357 auto* bg =
static_cast<TH2*
>(input.Get(Form(
"%s%d", bgName.c_str(), i)));
359 std::cerr <<
"Failed to find histogram '" << Form(
"%s%d", bgName.c_str(), i) <<
"', should have " << angles->NumberOfAngles() <<
" angles in total!" << std::endl;
360 std::cout << parser << std::endl;
363 auto* mixed =
static_cast<TH2*
>(input.Get(Form(
"%s%d", mixedName.c_str(), i)));
364 if(mixed ==
nullptr) {
365 std::cerr <<
"Failed to find histogram '" << Form(
"%s%d", mixedName.c_str(), i) <<
"', should have " << angles->NumberOfAngles() <<
" angles in total!" << std::endl;
366 std::cout << parser << std::endl;
372 prompt->Add(bg, -timeRandomNorm);
376 auto* proj = prompt->ProjectionX(Form(
"proj%d", i), prompt->GetYaxis()->FindBin(projGateLow), prompt->GetYaxis()->FindBin(projGateHigh));
377 auto* projMixed = mixed->ProjectionX(Form(
"projMixed%d", i), mixed->GetYaxis()->FindBin(projGateLow), mixed->GetYaxis()->FindBin(projGateHigh));
379 auto* projBg = prompt->ProjectionX(Form(
"projBg%d", i), prompt->GetYaxis()->FindBin(
bgLow[0]), prompt->GetYaxis()->FindBin(
bgHigh[0]));
380 auto* projMixedBg = mixed->ProjectionX(Form(
"projMixedBg%d", i), mixed->GetYaxis()->FindBin(
bgLow[0]), mixed->GetYaxis()->FindBin(
bgHigh[0]));
381 double bgGateWidth = prompt->GetYaxis()->FindBin(
bgHigh[0]) - prompt->GetYaxis()->FindBin(
bgLow[0]) + 1;
382 double mixedBgGateWidth = mixed->GetYaxis()->FindBin(
bgHigh[0]) - mixed->GetYaxis()->FindBin(
bgLow[0]) + 1;
383 for(
size_t g = 1; g <
bgLow.size(); ++g) {
384 projBg->Add(prompt->ProjectionX(Form(
"projBg%d", i), prompt->GetYaxis()->FindBin(
bgLow[g]), prompt->GetYaxis()->FindBin(
bgHigh[g])));
385 projMixedBg->Add(mixed->ProjectionX(Form(
"projMixedBg%d", i), mixed->GetYaxis()->FindBin(
bgLow[g]), mixed->GetYaxis()->FindBin(
bgHigh[g])));
386 bgGateWidth += prompt->GetYaxis()->FindBin(
bgHigh[g]) - prompt->GetYaxis()->FindBin(
bgLow[g]) + 1;
387 mixedBgGateWidth += mixed->GetYaxis()->FindBin(
bgHigh[g]) - mixed->GetYaxis()->FindBin(
bgLow[g]) + 1;
391 proj->Add(projBg, -(prompt->GetYaxis()->FindBin(projGateHigh) - prompt->GetYaxis()->FindBin(projGateLow) + 1) / bgGateWidth);
392 projMixed->Add(projMixedBg, -(mixed->GetYaxis()->FindBin(projGateHigh) - mixed->GetYaxis()->FindBin(projGateLow) + 1) / mixedBgGateWidth);
398 for(
size_t p = 0; p < peakParameterLow.size(); ++p) {
399 if(peakParameterLow[p] == peakParameterHigh[p]) {
401 }
else if(peakParameterLow[p] < peakParameterHigh[p]) {
403 peak.
GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
407 for(
auto bgPeak : bgPeakPos) {
408 auto* bgP =
new TRWPeak(std::get<0>(bgPeak));
410 if(std::get<1>(bgPeak) != -1. && std::get<2>(bgPeak) != -1. && std::get<1>(bgPeak) < std::get<2>(bgPeak)) {
411 bgP->GetFitFunction()->SetParLimits(1, std::get<1>(bgPeak), std::get<2>(bgPeak));
413 for(
size_t p = 0; p < bgPeakParameterLow.size(); ++p) {
414 if(bgPeakParameterLow[p] == bgPeakParameterHigh[p]) {
415 bgP->GetFitFunction()->FixParameter(p, bgPeakParameter[p]);
416 }
else if(bgPeakParameterLow[p] < bgPeakParameterHigh[p]) {
417 bgP->GetFitFunction()->SetParameter(p, bgPeakParameter[p]);
418 bgP->GetFitFunction()->SetParLimits(p, bgPeakParameterLow[p], bgPeakParameterHigh[p]);
423 for(
size_t p = 0; p < backgroundParameterLow.size(); ++p) {
424 if(backgroundParameterLow[p] == backgroundParameterHigh[p]) {
426 }
else if(backgroundParameterLow[p] < backgroundParameterHigh[p]) {
428 pf.
GetBackground()->SetParLimits(p, backgroundParameterLow[p], backgroundParameterHigh[p]);
433 pf.
Fit(proj,
"qretryfit");
436 logFile << std::setw(2) << i <<
" p "
437 << std::setw(10) << peak.
Centroid() <<
" +- " << std::setw(10) << peak.
CentroidErr() <<
" "
438 << std::setw(10) << peak.
Area() <<
" +- " << std::setw(10) << peak.
AreaErr() <<
" "
439 << std::setw(10) << peak.
FWHM() <<
" +- " << std::setw(10) << peak.
FWHMErr() <<
" "
444 for(
size_t p = 0; p < peakParameterLow.size(); ++p) {
445 if(peakParameterLow[p] == peakParameterHigh[p]) {
447 }
else if(peakParameterLow[p] < peakParameterHigh[p]) {
449 peakMixed.
GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
453 for(
auto bgPeak : bgPeakPos) {
454 auto* bgP =
new TRWPeak(std::get<0>(bgPeak));
456 if(std::get<1>(bgPeak) != -1. && std::get<2>(bgPeak) != -1. && std::get<1>(bgPeak) < std::get<2>(bgPeak)) {
457 bgP->GetFitFunction()->SetParLimits(1, std::get<1>(bgPeak), std::get<2>(bgPeak));
459 for(
size_t p = 0; p < bgPeakParameterLow.size(); ++p) {
460 if(bgPeakParameterLow[p] == bgPeakParameterHigh[p]) {
461 bgP->GetFitFunction()->FixParameter(p, bgPeakParameter[p]);
462 }
else if(bgPeakParameterLow[p] < bgPeakParameterHigh[p]) {
463 bgP->GetFitFunction()->SetParameter(p, bgPeakParameter[p]);
464 bgP->GetFitFunction()->SetParLimits(p, bgPeakParameterLow[p], bgPeakParameterHigh[p]);
469 for(
size_t p = 0; p < backgroundParameterLow.size(); ++p) {
470 if(backgroundParameterLow[p] == backgroundParameterHigh[p]) {
471 pfMixed.
GetBackground()->FixParameter(p, backgroundParameter[p]);
472 }
else if(backgroundParameterLow[p] < backgroundParameterHigh[p]) {
473 pfMixed.
GetBackground()->SetParameter(p, backgroundParameter[p]);
474 pfMixed.
GetBackground()->SetParLimits(p, backgroundParameterLow[p], backgroundParameterHigh[p]);
479 pfMixed.
Fit(projMixed,
"qretryfit");
482 logFile << std::setw(2) << i <<
" m "
483 << std::setw(8) << peakMixed.
Centroid() <<
" +- " << std::setw(8) << peakMixed.
CentroidErr() <<
" "
484 << std::setw(8) << peakMixed.
Area() <<
" +- " << std::setw(8) << peakMixed.
AreaErr() <<
" "
485 << std::setw(8) << peakMixed.
FWHM() <<
" +- " << std::setw(8) << peakMixed.
FWHMErr() <<
" "
494 rawAngularDistribution->SetPoint(i, angles->AverageAngle(i), peak.
Area());
495 rawAngularDistribution->SetPointError(i, 0., peak.
AreaErr());
496 mixedAngularDistribution->SetPoint(i, angles->AverageAngle(i), peakMixed.
Area());
497 mixedAngularDistribution->SetPointError(i, 0., peakMixed.
AreaErr());
498 angularDistribution->SetPoint(i, angles->AverageAngle(i), peak.
Area() / peakMixed.
Area());
499 angularDistribution->SetPointError(i, 0., peak.
Area() / peakMixed.
Area() * TMath::Sqrt(TMath::Power(peak.
AreaErr() / peak.
Area(), 2) + TMath::Power(peakMixed.
AreaErr() / peakMixed.
Area(), 2)));
501 rawChiSquares->SetPoint(i, angles->AverageAngle(i), peak.
GetReducedChi2());
502 mixedChiSquares->SetPoint(i, angles->AverageAngle(i), peakMixed.
GetReducedChi2());
504 std::cout <<
"Angle " << std::setw(3) << i <<
" of " << angles->NumberOfAngles() <<
" done\r" << std::flush;
506 std::cout <<
"Fitting of projections done." << std::endl;
509 auto* rawAngularDistributionCorr =
static_cast<TGraphErrors*
>(rawAngularDistribution->Clone(
"RawAngularDistributionCorrected"));
510 for(
int i = 0; i < rawAngularDistributionCorr->GetN(); ++i) {
511#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
512 if(angles->Count(rawAngularDistributionCorr->GetPointX(i)) != 0) {
513 rawAngularDistributionCorr->SetPointY(i, rawAngularDistributionCorr->GetPointY(i) / angles->Count(rawAngularDistributionCorr->GetPointX(i)));
514 rawAngularDistributionCorr->SetPointError(i, rawAngularDistributionCorr->GetErrorX(i), rawAngularDistributionCorr->GetErrorY(i) / angles->Count(rawAngularDistributionCorr->GetPointX(i)));
516 rawAngularDistributionCorr->SetPointY(i, 0);
520 rawAngularDistributionCorr->GetPoint(i, px, py);
521 if(angles->Count(px) != 0) {
522 rawAngularDistributionCorr->SetPoint(i, px, py / angles->Count(px));
523 rawAngularDistributionCorr->SetPointError(i, rawAngularDistributionCorr->GetErrorX(i), rawAngularDistributionCorr->GetErrorY(i) / angles->Count(px));
525 rawAngularDistributionCorr->SetPoint(i, px, 0);
529 auto* mixedAngularDistributionCorr =
static_cast<TGraphErrors*
>(mixedAngularDistribution->Clone(
"MixedAngularDistributionCorrected"));
530 for(
int i = 0; i < mixedAngularDistributionCorr->GetN(); ++i) {
531#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
532 if(angles->Count(mixedAngularDistributionCorr->GetPointX(i)) != 0) {
533 mixedAngularDistributionCorr->SetPointY(i, mixedAngularDistributionCorr->GetPointY(i) / angles->Count(mixedAngularDistributionCorr->GetPointX(i)));
534 mixedAngularDistributionCorr->SetPointError(i, mixedAngularDistributionCorr->GetErrorX(i), mixedAngularDistributionCorr->GetErrorY(i) / angles->Count(mixedAngularDistributionCorr->GetPointX(i)));
536 mixedAngularDistributionCorr->SetPointY(i, 0);
540 mixedAngularDistributionCorr->GetPoint(i, px, py);
541 if(angles->Count(px) != 0) {
542 mixedAngularDistributionCorr->SetPoint(i, px, py / angles->Count(px));
543 mixedAngularDistributionCorr->SetPointError(i, mixedAngularDistributionCorr->GetErrorX(i), mixedAngularDistributionCorr->GetErrorY(i) / angles->Count(px));
545 mixedAngularDistributionCorr->SetPoint(i, px, 0);
559 if(!theoryFile.empty()) {
560 TFile theory(theoryFile.c_str());
562 if(theory.IsOpen()) {
564 auto* z0 =
static_cast<TGraphErrors*
>(theory.Get(
"graph000"));
565 auto* z2 =
static_cast<TGraphErrors*
>(theory.Get(
"graph100"));
566 auto* z4 =
static_cast<TGraphErrors*
>(theory.Get(
"graph010"));
568 if(z0 !=
nullptr && z2 !=
nullptr && z4 !=
nullptr && z0->GetN() == z2->GetN() && z0->GetN() == z4->GetN()) {
572 if(z0->GetN() == angularDistribution->GetN() || ((angles->Grouping() || angles->Folding()) && z0->GetN() == (angles->Addback() ? 49 : 51))) {
574 if(z0->GetN() != angularDistribution->GetN()) {
575 angles->FoldOrGroup(z0, z2, z4);
578 std::vector<TGraph*> spin;
579 std::vector<double> spinLabel;
580 std::vector<std::vector<double>> parameters;
581 logFile << std::endl;
583 if(twoJLow.size() > 1 && twoJMiddle.size() == 1 && twoJHigh.size() == 1) {
584 logFile <<
"# Mixing method, high 2J = " << twoJHigh.at(0) <<
", middle 2J = " << twoJMiddle.at(0) <<
", low 2J = " << twoJLow.at(0) <<
" - " << twoJLow.back() << std::endl;
585 for(
auto twoJ : twoJLow) {
586 parameters.emplace_back();
588 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJMiddle.at(0), twoJ, parameters.back(), logFile));
589 spinLabel.push_back(twoJ / 2.);
591 }
else if(twoJLow.size() == 1 && twoJMiddle.size() > 1 && twoJHigh.size() == 1) {
592 logFile <<
"# Mixing method, high 2J = " << twoJHigh.at(0) <<
", middle 2J = " << twoJMiddle.at(0) <<
" - " << twoJMiddle.back() <<
", low 2J = " << twoJLow.at(0) << std::endl;
593 for(
auto twoJ : twoJMiddle) {
594 parameters.emplace_back();
596 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJ, twoJLow.at(0), parameters.back(), logFile));
597 spinLabel.push_back(twoJ / 2.);
599 }
else if(twoJLow.size() == 1 && twoJMiddle.size() == 1 && twoJHigh.size() > 1) {
600 logFile <<
"# Mixing method, high 2J = " << twoJHigh.at(0) <<
" - " << twoJHigh.back() <<
", middle 2J = " << twoJMiddle.at(0) <<
", low 2J = " << twoJLow.at(0) << std::endl;
601 for(
auto twoJ : twoJHigh) {
602 parameters.emplace_back();
604 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJ, twoJMiddle.at(0), twoJLow.at(0), parameters.back(), logFile));
605 spinLabel.push_back(twoJ / 2.);
608 logFile <<
"# Mixing method, high 2J = " << twoJHigh.at(0) <<
", middle 2J = " << twoJMiddle.at(0) <<
", low 2J = " << twoJLow.at(0) << std::endl;
609 parameters.emplace_back();
611 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJMiddle.at(0), twoJLow.at(0), parameters.back(), logFile));
612 spinLabel.push_back(twoJHigh.at(0) / 2.);
616 auto* canvas =
new TCanvas;
619 double min = TMath::MinElement(spin.at(0)->GetN(), spin.at(0)->GetY());
620 double max = TMath::MaxElement(spin.at(0)->GetN(), spin.at(0)->GetY());
621 for(
size_t i = 1; i < spin.size(); ++i) {
622 min = TMath::Min(min, TMath::MinElement(spin.at(i)->GetN(), spin.at(i)->GetY()));
623 max = TMath::Max(max, TMath::MaxElement(spin.at(i)->GetN(), spin.at(i)->GetY()));
625 min = TMath::Min(min, confidenceLevel);
629 for(first = 0; first < spin.size(); ++first) {
630 if(spin[first]->GetN() > 1) {
break; }
633 spin[first]->SetTitle(
"");
634 spin[first]->SetMinimum(0.9 * min);
635 spin[first]->SetMaximum(1.1 * max);
637 for(
size_t i = 0; i < spin.size(); ++i) {
638 spin[i]->SetLineColor(i + 1);
639 spin[i]->SetMarkerColor(i + 1);
640 spin[i]->SetLineWidth(2);
643 spin[first]->Draw(
"ac");
644 for(
size_t i = 0; i < spin.size(); ++i) {
645 if(i == first) {
continue; }
646 if(spin[i]->GetN() > 1) {
653 auto* confidenceLevelLine =
new TLine(-1.5, confidenceLevel, 1.5, confidenceLevel);
655 confidenceLevelLine->Draw();
657#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
658 auto* legend =
new TLegend(0.1, 0.3);
660 auto* legend =
new TLegend(0.7, 0.6, 0.8, 0.9);
662 for(
size_t i = 0; i < spin.size(); ++i) {
663 if(spin[i]->GetN() == 1) {
664 legend->AddEntry(spin[i], Form(
"J = %.1f", spinLabel[i]),
"p");
666 legend->AddEntry(spin[i], Form(
"J = %.1f", spinLabel[i]),
"l");
674 spin[first]->GetHistogram()->GetXaxis()->SetRangeUser(-1.5, 1.5);
675 spin[first]->GetHistogram()->GetXaxis()->SetTitle(
"atan(#delta) [rad]");
676 spin[first]->GetHistogram()->GetXaxis()->CenterTitle();
677 spin[first]->GetHistogram()->GetYaxis()->SetTitle(
"red. #chi^{2}");
678 spin[first]->GetHistogram()->GetYaxis()->CenterTitle();
682 z0->Write(
"graph000");
683 z2->Write(
"graph010");
684 z4->Write(
"graph100");
685 for(
size_t i = 0; i < spin.size(); ++i) {
686 spin[i]->Write(Form(
"spin%d",
static_cast<int>(i)));
688 canvas->Write(
"MixingCanvas");
691 std::vector<TGraphErrors*> spinFit(spin.size(),
new TGraphErrors(angularDistribution->GetN()));
692 auto* x = angularDistribution->GetX();
693 for(
int p = 0; p < angularDistribution->GetN(); ++p) {
694 for(
size_t i = 0; i < spinFit.size(); ++i) {
695 spinFit[i]->SetPoint(p, x[p], parameters[i][0] * ((1. - parameters[i][1] - parameters[i][2]) * z0->Eval(x[p]) + parameters[i][1] * z2->Eval(x[p]) + parameters[i][2] * z4->Eval(x[p])));
696 spinFit[i]->SetPointError(p, 0., TMath::Sqrt(parameters[i][0] * (TMath::Power((1. - parameters[i][1] - parameters[i][2]) *
GetYError(z0, x[p]), 2) + TMath::Power(parameters[i][1] *
GetYError(z2, x[p]), 2) + TMath::Power(parameters[i][2] *
GetYError(z4, x[p]), 2))));
699 for(
size_t i = 0; i < spin.size(); ++i) {
700 spinFit[i]->SetLineColor(i + 1);
701 spinFit[i]->SetMarkerColor(i + 1);
702 spinFit[i]->SetLineWidth(2);
703 spinFit[i]->Write(Form(
"SpinFit%d",
static_cast<int>(i)));
704 output.WriteObject(¶meters[i], Form(
"Parameters%d",
static_cast<int>(i)));
708 auto a2a4Parameters =
A2a4Method(angularDistribution, z0, z2, z4);
709 output.WriteObject(&a2a4Parameters,
"ParametersA2a4Fit");
711 std::cerr <<
"Mismatch between theory and data (" << z0->GetN() <<
" != " << angularDistribution->GetN() <<
"?) or neither grouping and folding are enabled (" << std::boolalpha << angles->Grouping() <<
", " << angles->Folding() <<
") or the ungrouped/folded theory doesn't match the required size for " << (angles->Addback() ?
"addback" :
"singles") <<
" data (" << (angles->Addback() ? 49 : 51) <<
")" << std::endl;
714 std::cerr <<
"Failed to find z0 (" << z0 <<
", \"graph000\"), z2 (" << z2 <<
", \"graph010\"), or z4 (" << z4 <<
", \"graph100\") in " << theoryFile <<
", or they had mismatched sizes" << std::endl;
717 std::cerr <<
"Failed to open " << theoryFile << std::endl;
720 std::cout <<
"No file with simulation results (--theory flag), so we won't produce chi2 vs mixing angle plot." << std::endl;
723 rawAngularDistribution->Write();
724 angularDistribution->Write();
725 mixedAngularDistribution->Write();
726 rawChiSquares->Write();
727 mixedChiSquares->Write();
728 rawAngularDistributionCorr->Write();
729 mixedAngularDistributionCorr->Write();
873TGraph*
MixingMethod(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4,
int twoJhigh,
int twoJmid,
int twoJlow, std::vector<double>& bestParameters, std::ofstream& logFile)
875 logFile <<
"# high 2J " << twoJhigh <<
", middle 2J " << twoJmid <<
", low 2J " << twoJlow << std::endl;
876 logFile <<
"# a0 a2 a4 red.chi^2" << std::endl;
877 TGraph* result =
nullptr;
878 Ac ac(data, z0, z2, z4);
879 ROOT::Fit::Fitter fitter;
881 fitter.SetFCN(nPar, ac);
882 for(
int i = 0; i < nPar; ++i) {
884 fitter.Config().ParSettings(i) = ROOT::Fit::ParameterSettings(Form(
"a_{%d}", 2 * i), 0.5, 0.0001, -10., 10.);
886 fitter.Config().MinimizerOptions().SetPrintLevel(0);
887 fitter.Config().SetMinimizer(
"Minuit2",
"Migrad");
892 double j1 = 0.5 * twoJhigh;
893 double j2 = 0.5 * twoJmid;
894 double j3 = 0.5 * twoJlow;
898 int l1a = TMath::Abs(twoJhigh - twoJmid) / 2;
899 if(l1a == 0) { l1a = 1; }
904 int l2a = TMath::Abs(twoJmid - twoJlow) / 2;
905 if(l2a == 0) { l2a = 1; }
909 if((twoJhigh == 0 && twoJmid == 0) || (twoJmid == 0 && twoJlow == 0)) {
910 std::cout <<
"!!!!!!!!!!!!!!! ERROR !!!!!!!!!!!!!!!" << std::endl
911 <<
"Can't have gamma transition between J=0 states (high " << twoJhigh <<
", mid " << twoJmid <<
", low " << twoJlow <<
")." << std::endl
912 <<
"Aborting..." << std::endl;
915 if(l1a == TMath::Abs(twoJhigh + twoJmid) / 2) {
922 if(l2a == TMath::Abs(twoJmid + twoJlow) / 2) {
943 double mixingAngle1Minimum = -TMath::Pi() / 2;
944 double mixingAngle1Maximum = TMath::Pi() / 2;
946 double stepSize1 = (mixingAngle1Maximum - mixingAngle1Minimum) / steps1;
948 double mixingAngle2Minimum = -TMath::Pi() / 2;
949 double mixingAngle2Maximum = TMath::Pi() / 2;
951 double stepSize2 = (mixingAngle2Maximum - mixingAngle2Minimum) / steps2;
955 mixingAngle1Minimum = 0;
959 mixingAngle2Minimum = 0;
963 result =
new TGraph(steps1);
964 double minChi2 = 1e6;
965 std::vector<double> bestErrors;
966 double bestMixingAngle1 = 0.;
967 double bestMixingAngle2 = 0.;
968 for(
int i = 0; i < steps1; i++) {
969 double mixangle1 = mixingAngle1Minimum + i * stepSize1;
970 double delta1 = TMath::Tan(mixangle1);
971 for(
int j = 0; j < steps2; j++) {
972 double mixangle2 = mixingAngle2Minimum + j * stepSize2;
973 double delta2 = TMath::Tan(mixangle2);
977 fitter.Config().ParSettings(1).Set(
"a_{2}", a2);
978 fitter.Config().ParSettings(1).Fix();
982 fitter.Config().ParSettings(2).Set(
"a_{4}", a4);
983 fitter.Config().ParSettings(2).Fix();
984 if(!fitter.FitFCN()) {
985 std::cerr << i <<
", " << j <<
": Fit failed using a2 " << a2 <<
", a4 " << a4 << std::endl;
988 auto fitResult = fitter.Result();
990 double chi2 = fitResult.MinFcnValue() / (ac.
Np() - fitResult.NFreeParameters());
993 result->SetPoint(i, mixangle1, chi2);
996 bestParameters = fitResult.Parameters();
997 bestErrors = fitResult.Errors();
998 bestMixingAngle1 = mixangle1;
999 bestMixingAngle2 = mixangle2;
1001 logFile << std::setw(12) << fitResult.Parameter(0) <<
" " << std::setw(10) << a2 <<
" " << std::setw(10) << a4 <<
" " << std::setw(10) << fitResult.MinFcnValue() / (ac.
Np() - fitResult.NFreeParameters()) << std::endl;
1006 auto minIndex = TMath::LocMin(result->GetN(), result->GetY());
1007 double x1 = std::numeric_limits<double>::quiet_NaN();
1008 double x2 = std::numeric_limits<double>::quiet_NaN();
1009 double y1 = std::numeric_limits<double>::quiet_NaN();
1010 double y2 = std::numeric_limits<double>::quiet_NaN();
1011 for(
int i = minIndex; i < result->GetN(); ++i) {
1012 if(result->GetPointY(i) > minChi2 + 1.) {
1013 x1 = result->GetPointX(i);
1014 x2 = result->GetPointX(i - 1);
1015 y1 = result->GetPointY(i);
1016 y2 = result->GetPointY(i - 1);
1020 double upperLimit = std::numeric_limits<double>::quiet_NaN();
1021 if(!std::isnan(x1)) {
1022 upperLimit = x1 - (x2 - x1) / (y2 - y1) * y1 + (x2 - x1) / (y2 - y1) * (minChi2 + 1.);
1024 x1 = std::numeric_limits<double>::quiet_NaN();
1025 for(
int i = minIndex; i >= 0; --i) {
1026 if(result->GetPointY(i) > minChi2 + 1.) {
1027 x1 = result->GetPointX(i);
1028 x2 = result->GetPointX(i + 1);
1029 y1 = result->GetPointY(i);
1030 y2 = result->GetPointY(i + 1);
1034 double lowerLimit = std::numeric_limits<double>::quiet_NaN();
1035 if(!std::isnan(x1)) {
1036 lowerLimit = x1 - (x2 - x1) / (y2 - y1) * y1 + (x2 - x1) / (y2 - y1) * (minChi2 + 1.);
1040 auto* fit =
static_cast<TGraphErrors*
>(z0->Clone(
"fit"));
1041 for(
int i = 0; i < fit->GetN(); ++i) {
1042#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
1043 fit->SetPointY(i, bestParameters[0] * ((1. - bestParameters[1] - bestParameters[2]) * z0->GetPointY(i) + bestParameters[1] * z2->GetPointY(i) + bestParameters[2] * z4->GetPointY(i)));
1045 fit->SetPoint(i, fit->GetX()[i], bestParameters[0] * ((1. - bestParameters[1] - bestParameters[2]) * z0->GetY()[i] + bestParameters[1] * z2->GetY()[i] + bestParameters[2] * z4->GetY()[i]));
1047 fit->SetPointError(i, 0., std::abs(bestParameters[0]) * TMath::Sqrt(TMath::Power((1. - bestParameters[1] - bestParameters[2]) * z0->GetErrorY(i), 2) + TMath::Power(bestParameters[1] * z2->GetErrorY(i), 2) + TMath::Power(bestParameters[2] * z4->GetErrorY(i), 2)));
1050 auto* residual =
static_cast<TGraphErrors*
>(data->Clone(
"residual"));
1051 for(
int i = 0; i < residual->GetN(); ++i) {
1052#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
1053 residual->SetPointY(i, data->GetPointY(i) - fit->GetPointY(i));
1055 residual->SetPoint(i, residual->GetX()[i], data->GetY()[i] - fit->GetY()[i]);
1057 residual->SetPointError(i, 0., TMath::Sqrt(TMath::Power(data->GetErrorY(i), 2) + TMath::Power(fit->GetErrorY(i), 2)));
1060 auto* canvas =
new TCanvas;
1061 TMultiGraph* multiGraph =
nullptr;
1063 if(steps1 == 1 && steps2 == 1) {
1064 std::cout <<
"Analyzed cascade " << twoJhigh / 2. <<
" -> " << twoJmid / 2. <<
" -> " << twoJlow / 2. <<
": red. chi^2 " << std::setw(12) << minChi2 <<
", a0 " << std::setw(12) << bestParameters[0] <<
" +- " << std::setw(12) << bestErrors[0] <<
", a2 " << std::setw(12) << bestParameters[1] <<
", a4 " << std::setw(12) << bestParameters[2] << std::endl;
1065 multiGraph =
PlotCanvas(data, fit, residual, bestParameters, bestErrors, minChi2);
1068 std::cout <<
"Varied mixing angle for cascade " << twoJhigh / 2. <<
" -> " << twoJmid / 2. <<
" -> " << twoJlow / 2. <<
": best red. chi^2 " << std::setw(12) << minChi2 <<
", at mixing angle " << std::setw(12) << bestMixingAngle2 <<
", a0 " << std::setw(12) << bestParameters[0] <<
" +- " << std::setw(12) << bestErrors[0] <<
", a2 " << std::setw(12) << bestParameters[1] <<
", a4 " << std::setw(12) << bestParameters[2] << std::endl;
1069 multiGraph =
PlotCanvas(data, fit, residual, bestParameters, bestErrors, minChi2, Form(
"mixing angle = %f", bestMixingAngle2));
1070 }
else if(steps2 == 1) {
1071 std::cout <<
"Varied mixing angle for cascade " << twoJhigh / 2. <<
" -> " << twoJmid / 2. <<
" -> " << twoJlow / 2. <<
": best red. chi^2 " << std::setw(12) << minChi2 <<
", at mixing angle " << std::setw(12) << bestMixingAngle1 <<
", a0 " << std::setw(12) << bestParameters[0] <<
" +- " << std::setw(12) << bestErrors[0] <<
", a2 " << std::setw(12) << bestParameters[1] <<
", a4 " << std::setw(12) << bestParameters[2] << std::endl;
1072 multiGraph =
PlotCanvas(data, fit, residual, bestParameters, bestErrors, minChi2, Form(
"mixing angle = %f (-%f/+%f)", bestMixingAngle1, bestMixingAngle1 - lowerLimit, upperLimit - bestMixingAngle1));
1074 std::cout <<
"Varied mixing angle for cascade " << twoJhigh / 2. <<
" -> " << twoJmid / 2. <<
" -> " << twoJlow / 2. <<
": best red. chi^2 " << std::setw(12) << minChi2 <<
", at mixing angle " << std::setw(12) << bestMixingAngle1 <<
" / " << std::setw(12) << bestMixingAngle2 <<
", a0 " << std::setw(12) << bestParameters[0] <<
" +- " << std::setw(12) << bestErrors[0] <<
", a2 " << std::setw(12) << bestParameters[1] <<
", a4 " << std::setw(12) << bestParameters[2] << std::endl;
1075 multiGraph =
PlotCanvas(data, fit, residual, bestParameters, bestErrors, minChi2, Form(
"mixing angle = %f/%f", bestMixingAngle1, bestMixingAngle2));
1079 fit->Write(Form(
"BestMixingFit%d_%d_%d", twoJhigh, twoJmid, twoJlow));
1080 residual->Write(Form(
"Residual%d_%d_%d", twoJhigh, twoJmid, twoJlow));
1081 multiGraph->Write(Form(
"FitComparison%d_%d_%d", twoJhigh, twoJmid, twoJlow));
1082 canvas->Write(Form(
"Canvas%d_%d_%d", twoJhigh, twoJmid, twoJlow));