54int main(
int argc,
char** argv)
62 double projGateLow = -1.;
63 double projGateHigh = -1.;
64 double projBgLow = -1.;
65 double projBgHigh = -1.;
66 double timeRandomNorm = -1.;
69 double peakHigh = -1.;
70 std::string baseName =
"AngularCorrelation";
71 std::string inputFile;
72 std::string theoryFile;
73 std::string outputFile =
"AngularCorrelations.root";
74 std::string settingsFile;
78 parser.
option(
"h help ?", &help,
true).description(
"Show this help message.");
79 parser.
option(
"projection-low proj-low", &projGateLow,
true).description(
"Low edge of projection gate.");
80 parser.
option(
"projection-high proj-high", &projGateHigh,
true).description(
"High edge of projection gate.");
81 parser.
option(
"background-low bg-low", &projBgLow,
true).description(
"Low edge of background gate (for multiple gates use settings file).");
82 parser.
option(
"background-high bg-high", &projBgHigh,
true).description(
"High edge of background gate (for multiple gates use settings file).");
83 parser.
option(
"time-random-normalization", &timeRandomNorm,
true).description(
"Normalization factor for subtraction of time-random matrix. If negative (default) it will be calculated automatically.");
84 parser.
option(
"base-name", &baseName,
true).description(
"Base name of matrices.").default_value(
"AngularCorrelation");
85 parser.
option(
"peak-pos peak", &peakPos,
true).description(
"Peak position (for multiple peaks use settings file).");
86 parser.
option(
"peak-low", &peakLow,
true).description(
"Low edge of peak fit.");
87 parser.
option(
"peak-high", &peakHigh,
true).description(
"High edge of peak fit.");
88 parser.
option(
"settings", &settingsFile,
true).description(
"Settings file with user settings, these do not overwrite anything provided on command line!");
89 parser.
option(
"input", &inputFile,
true).description(
"Input file with gamma-gamma matrices for each angle (coincident, time-random, and event mixed).");
90 parser.
option(
"theory", &theoryFile,
true).description(
"File with simulated z0-, z2-, and z4-graphs.");
91 parser.
option(
"output", &outputFile,
true).description(
"Name of output file, default is \"AngularCorrelations.root\".");
93 parser.
parse(argc, argv,
true);
96 std::cout << parser << std::endl;
101 if(inputFile.empty()) {
102 std::cerr <<
"Need an input file!" << std::endl;
106 TFile input(inputFile.c_str());
108 if(!input.IsOpen()) {
109 std::cerr <<
"Failed to open input file " << inputFile << std::endl;
113 auto* settings =
static_cast<TUserSettings*
>(input.Get(
"UserSettings"));
117 if(!settingsFile.empty()) {
118 if(settings ==
nullptr) {
121 settings->ReadSettings(settingsFile);
126 if(settings ==
nullptr || settings->empty()) {
127 std::cerr <<
"Failed to get user settings from input file." << std::endl;
132 if(projGateLow == -1.) { projGateLow = settings->GetDouble(
"Projection.Low"); }
133 if(projGateHigh == -1.) { projGateHigh = settings->GetDouble(
"Projection.High"); }
134 if(timeRandomNorm == -1.) { timeRandomNorm = settings->GetDouble(
"TimeRandomNormalization"); }
135 if(peakPos == -1.) { peakPos = settings->GetDouble(
"Peak.Position"); }
136 if(peakLow == -1.) { peakLow = settings->GetDouble(
"Peak.Low"); }
137 if(peakHigh == -1.) { peakHigh = settings->GetDouble(
"Peak.High"); }
139 std::vector<double> bgPeakPos;
140 std::vector<double>
bgLow;
141 std::vector<double>
bgHigh;
143 if(projBgLow == -1. || projBgHigh == -1. || projBgLow >= projBgHigh) {
144 for(
int i = 1;; ++i) {
146 auto low = settings->GetDouble(Form(
"Background.%d.Low", i),
true);
147 auto high = settings->GetDouble(Form(
"Background.%d.High", i),
true);
149 std::cout << i <<
". background gate, low edge not lower than the high edge: " << low <<
" >= " << high << std::endl;
152 bgLow.push_back(low);
154 }
catch(std::out_of_range& e) {
159 bgLow.push_back(projBgLow);
160 bgHigh.push_back(projBgHigh);
164 for(
int i = 1;; ++i) {
166 auto pos = settings->GetDouble(Form(
"BackgroundPeak.Position.%d", i),
true);
167 if(pos <= peakLow || pos >= peakHigh) {
168 std::cout << i <<
". background peak outside of fit range: " << pos <<
" <= " << peakLow <<
" or " << pos <<
" >= " << peakHigh << std::endl;
171 bgPeakPos.push_back(pos);
172 }
catch(std::out_of_range& e) {
180 std::vector<double> peakParameter(6, -2.);
181 std::vector<double> peakParameterLow(6, 0.);
182 std::vector<double> peakParameterHigh(6, -1.);
183 std::vector<double> bgPeakParameter(6, -2.);
184 std::vector<double> bgPeakParameterLow(6, 0.);
185 std::vector<double> bgPeakParameterHigh(6, -1.);
186 for(
size_t i = 0; i < peakParameterLow.size(); ++i) {
187 peakParameter[i] = settings->GetDouble(Form(
"Peak.Parameter.%d",
static_cast<int>(i)), -2.);
188 peakParameterLow[i] = settings->GetDouble(Form(
"Peak.Parameter.%d.Low",
static_cast<int>(i)), 0.);
189 peakParameterHigh[i] = settings->GetDouble(Form(
"Peak.Parameter.%d.High",
static_cast<int>(i)), -1.);
190 bgPeakParameter[i] = settings->GetDouble(Form(
"BackgroundPeak.Parameter.%d",
static_cast<int>(i)), -2.);
191 bgPeakParameterLow[i] = settings->GetDouble(Form(
"BackgroundPeak.Parameter.%d.Low",
static_cast<int>(i)), 0.);
192 bgPeakParameterHigh[i] = settings->GetDouble(Form(
"BackgroundPeak.Parameter.%d.High",
static_cast<int>(i)), -1.);
195 if(peakParameterLow[i] <= peakParameterHigh[i] && (peakParameter[i] < peakParameterLow[i] || peakParameterHigh[i] < peakParameter[i])) {
196 bool output = peakParameter[i] != -2.;
197 if(output) { std::cout <<
"Warning, " << i <<
". peak parameter (" << peakParameter[i] <<
") is out of range " << peakParameterLow[i] <<
" - " << peakParameterHigh[i] <<
", resetting it to "; }
198 peakParameter[i] = (peakParameterHigh[i] + peakParameterLow[i]) / 2.;
199 if(output) { std::cout << peakParameter[i] << std::endl; }
201 if(bgPeakParameterLow[i] <= bgPeakParameterHigh[i] && (bgPeakParameter[i] < bgPeakParameterLow[i] || bgPeakParameterHigh[i] < bgPeakParameter[i])) {
202 bool output = bgPeakParameter[i] != -2.;
203 if(output) { std::cout <<
"Warning, " << i <<
". background peak parameter (" << bgPeakParameter[i] <<
") is out of range " << bgPeakParameterLow[i] <<
" - " << bgPeakParameterHigh[i] <<
", resetting it to "; }
204 bgPeakParameter[i] = (bgPeakParameterHigh[i] + bgPeakParameterLow[i]) / 2.;
205 if(output) { std::cout << bgPeakParameter[i] << std::endl; }
210 std::vector<double> backgroundParameter(4, -2.);
211 std::vector<double> backgroundParameterLow(4, 0.);
212 std::vector<double> backgroundParameterHigh(4, -1.);
213 backgroundParameter[0] = settings->GetDouble(
"Background.Offset", -2.);
214 backgroundParameterLow[0] = settings->GetDouble(
"Background.Offset.Low", 0.);
215 backgroundParameterHigh[0] = settings->GetDouble(
"Background.Offset.High", -1.);
216 backgroundParameter[1] = settings->GetDouble(
"Background.Linear", -2.);
217 backgroundParameterLow[1] = settings->GetDouble(
"Background.Linear.Low", 0.);
218 backgroundParameterHigh[1] = settings->GetDouble(
"Background.Linear.High", -1.);
219 backgroundParameter[2] = settings->GetDouble(
"Background.Quadratic", -2.);
220 backgroundParameterLow[2] = settings->GetDouble(
"Background.Quadratic.Low", 0.);
221 backgroundParameterHigh[2] = settings->GetDouble(
"Background.Quadratic.High", -1.);
222 backgroundParameter[3] = settings->GetDouble(
"Background.Xoffset", -2.);
223 backgroundParameterLow[3] = settings->GetDouble(
"Background.Xoffset.Low", 0.);
224 backgroundParameterHigh[3] = settings->GetDouble(
"Background.Xoffset.High", -1.);
225 for(
size_t i = 0; i < backgroundParameter.size(); ++i) {
226 if(backgroundParameterLow[i] <= backgroundParameterHigh[i] && (backgroundParameter[i] < backgroundParameterLow[i] || backgroundParameterHigh[i] < backgroundParameter[i])) {
227 bool output = backgroundParameter[i] != -2.;
228 if(output) { std::cout <<
"Warning, " << i <<
". background parameter (" << backgroundParameter[i] <<
") is out of range " << backgroundParameterLow[i] <<
" - " << backgroundParameterHigh[i] <<
", resetting it to "; }
229 backgroundParameter[i] = (backgroundParameterHigh[i] + backgroundParameterLow[i]) / 2.;
230 if(output) { std::cout << backgroundParameter[i] << std::endl; }
236 std::vector<int> twoJLow = settings->GetIntVector(
"TwoJ.Low");
237 std::vector<int> twoJMiddle = settings->GetIntVector(
"TwoJ.Middle");
238 std::vector<int> twoJHigh = settings->GetIntVector(
"TwoJ.High");
241 double confidenceLevel = settings->GetDouble(
"ConfidenceLevel", 1.535);
244 if(projGateLow >= projGateHigh) {
245 std::cerr <<
"Need a projection gate with a low edge that is smaller than the high edge, " << projGateLow <<
" >= " << projGateHigh << std::endl;
249 if(peakLow >= peakHigh) {
250 std::cerr <<
"Need a fit range with a low edge that is smaller than the high edge, " << peakLow <<
" >= " << peakHigh << std::endl;
254 if(peakPos >= peakHigh || peakPos <= peakLow) {
255 std::cerr <<
"Need a peak within the fit range, " << peakPos <<
" not within " << peakLow <<
" - " << peakHigh << std::endl;
259 if(timeRandomNorm <= 0.) {
260 std::cerr <<
"Need a positive normalization factor for time random subtraction" << std::endl;
267 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;
272 auto* angles =
static_cast<TGriffinAngles*
>(input.Get(
"GriffinAngles"));
274 if(angles ==
nullptr) {
275 std::cerr <<
"Failed to find 'GriffinAngles' in '" << inputFile <<
"'" << std::endl;
284 TFile output(outputFile.c_str(),
"recreate");
286 auto* rawAngularDistribution =
new TGraphErrors(angles->NumberOfAngles());
287 rawAngularDistribution->SetName(
"RawAngularDistribution");
288 auto* angularDistribution =
new TGraphErrors(angles->NumberOfAngles());
289 angularDistribution->SetName(
"AngularDistribution");
290 auto* mixedAngularDistribution =
new TGraphErrors(angles->NumberOfAngles());
291 mixedAngularDistribution->SetName(
"MixedAngularDistribution");
292 auto* rawChiSquares =
new TGraph(angles->NumberOfAngles());
293 rawChiSquares->SetName(
"RawChiSquares");
294 auto* mixedChiSquares =
new TGraph(angles->NumberOfAngles());
295 mixedChiSquares->SetName(
"MixedChiSquares");
300#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
301 auto* fitDir = output.mkdir(
"fits",
"Projections with fits",
true);
303 auto* fitDir = output.mkdir(
"fits",
"Projections with fits");
307 for(
int i = 0; i < angles->NumberOfAngles(); ++i) {
309 auto* prompt =
static_cast<TH2*
>(input.Get(Form(
"AngularCorrelation%d", i)));
310 if(prompt ==
nullptr) {
311 std::cerr <<
"Failed to find histogram '" << Form(
"AngularCorrelation%d", i) <<
"', should have " << angles->NumberOfAngles() <<
" angles in total!" << std::endl;
314 auto* bg =
static_cast<TH2*
>(input.Get(Form(
"AngularCorrelationBG%d", i)));
316 std::cerr <<
"Failed to find histogram '" << Form(
"AngularCorrelationBG%d", i) <<
"', should have " << angles->NumberOfAngles() <<
" angles in total!" << std::endl;
319 auto* mixed =
static_cast<TH2*
>(input.Get(Form(
"AngularCorrelationMixed%d", i)));
320 if(mixed ==
nullptr) {
321 std::cerr <<
"Failed to find histogram '" << Form(
"AngularCorrelationMixed%d", i) <<
"', should have " << angles->NumberOfAngles() <<
" angles in total!" << std::endl;
327 prompt->Add(bg, -timeRandomNorm);
331 auto* proj = prompt->ProjectionX(Form(
"proj%d", i), prompt->GetYaxis()->FindBin(projGateLow), prompt->GetYaxis()->FindBin(projGateHigh));
332 auto* projMixed = mixed->ProjectionX(Form(
"projMixed%d", i), mixed->GetYaxis()->FindBin(projGateLow), mixed->GetYaxis()->FindBin(projGateHigh));
334 auto* projBg = prompt->ProjectionX(Form(
"projBg%d", i), prompt->GetYaxis()->FindBin(
bgLow[0]), prompt->GetYaxis()->FindBin(
bgHigh[0]));
335 auto* projMixedBg = mixed->ProjectionX(Form(
"projMixedBg%d", i), mixed->GetYaxis()->FindBin(
bgLow[0]), mixed->GetYaxis()->FindBin(
bgHigh[0]));
336 double bgGateWidth = prompt->GetYaxis()->FindBin(
bgHigh[0]) - prompt->GetYaxis()->FindBin(
bgLow[0]) + 1;
337 double mixedBgGateWidth = mixed->GetYaxis()->FindBin(
bgHigh[0]) - mixed->GetYaxis()->FindBin(
bgLow[0]) + 1;
338 for(
size_t g = 1; g <
bgLow.size(); ++g) {
339 projBg->Add(prompt->ProjectionX(Form(
"projBg%d", i), prompt->GetYaxis()->FindBin(
bgLow[g]), prompt->GetYaxis()->FindBin(
bgHigh[g])));
340 projMixedBg->Add(mixed->ProjectionX(Form(
"projMixedBg%d", i), mixed->GetYaxis()->FindBin(
bgLow[g]), mixed->GetYaxis()->FindBin(
bgHigh[g])));
341 bgGateWidth += prompt->GetYaxis()->FindBin(
bgHigh[g]) - prompt->GetYaxis()->FindBin(
bgLow[g]) + 1;
342 mixedBgGateWidth += mixed->GetYaxis()->FindBin(
bgHigh[g]) - mixed->GetYaxis()->FindBin(
bgLow[g]) + 1;
346 proj->Add(projBg, -(prompt->GetYaxis()->FindBin(projGateHigh) - prompt->GetYaxis()->FindBin(projGateLow) + 1) / bgGateWidth);
347 projMixed->Add(projMixedBg, -(mixed->GetYaxis()->FindBin(projGateHigh) - mixed->GetYaxis()->FindBin(projGateLow) + 1) / mixedBgGateWidth);
352 for(
size_t p = 0; p < peakParameterLow.size(); ++p) {
353 if(peakParameterLow[p] == peakParameterHigh[p]) {
355 }
else if(peakParameterLow[p] < peakParameterHigh[p]) {
357 peak.
GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
361 for(
auto bgPeak : bgPeakPos) {
362 auto* bgP =
new TRWPeak(bgPeak);
363 for(
size_t p = 0; p < bgPeakParameterLow.size(); ++p) {
364 if(bgPeakParameterLow[p] == bgPeakParameterHigh[p]) {
365 bgP->GetFitFunction()->FixParameter(p, bgPeakParameter[p]);
366 }
else if(bgPeakParameterLow[p] < bgPeakParameterHigh[p]) {
367 bgP->GetFitFunction()->SetParameter(p, bgPeakParameter[p]);
368 bgP->GetFitFunction()->SetParLimits(p, bgPeakParameterLow[p], bgPeakParameterHigh[p]);
373 for(
size_t p = 0; p < backgroundParameterLow.size(); ++p) {
374 if(backgroundParameterLow[p] == backgroundParameterHigh[p]) {
376 }
else if(backgroundParameterLow[p] < backgroundParameterHigh[p]) {
378 pf.
GetBackground()->SetParLimits(p, backgroundParameterLow[p], backgroundParameterHigh[p]);
381 pf.
Fit(proj,
"qretryfit");
385 for(
size_t p = 0; p < peakParameterLow.size(); ++p) {
386 if(peakParameterLow[p] == peakParameterHigh[p]) {
388 }
else if(peakParameterLow[p] < peakParameterHigh[p]) {
390 peakMixed.
GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
394 for(
auto bgPeak : bgPeakPos) {
395 auto* bgP =
new TRWPeak(bgPeak);
396 for(
size_t p = 0; p < bgPeakParameterLow.size(); ++p) {
397 if(bgPeakParameterLow[p] == bgPeakParameterHigh[p]) {
398 bgP->GetFitFunction()->FixParameter(p, bgPeakParameter[p]);
399 }
else if(bgPeakParameterLow[p] < bgPeakParameterHigh[p]) {
400 bgP->GetFitFunction()->SetParameter(p, bgPeakParameter[p]);
401 bgP->GetFitFunction()->SetParLimits(p, bgPeakParameterLow[p], bgPeakParameterHigh[p]);
406 for(
size_t p = 0; p < backgroundParameterLow.size(); ++p) {
407 if(backgroundParameterLow[p] == backgroundParameterHigh[p]) {
408 pfMixed.
GetBackground()->FixParameter(p, backgroundParameter[p]);
409 }
else if(backgroundParameterLow[p] < backgroundParameterHigh[p]) {
410 pfMixed.
GetBackground()->SetParameter(p, backgroundParameter[p]);
411 pfMixed.
GetBackground()->SetParLimits(p, backgroundParameterLow[p], backgroundParameterHigh[p]);
414 pfMixed.
Fit(projMixed,
"qretryfit");
422 rawAngularDistribution->SetPoint(i, angles->AverageAngle(i), peak.
Area());
423 rawAngularDistribution->SetPointError(i, 0., peak.
AreaErr());
424 mixedAngularDistribution->SetPoint(i, angles->AverageAngle(i), peakMixed.
Area());
425 mixedAngularDistribution->SetPointError(i, 0., peakMixed.
AreaErr());
426 angularDistribution->SetPoint(i, angles->AverageAngle(i), peak.
Area() / peakMixed.
Area());
427 angularDistribution->SetPointError(i, 0., peak.
Area() / peakMixed.
Area() * TMath::Sqrt(TMath::Power(peak.
AreaErr() / peak.
Area(), 2) + TMath::Power(peakMixed.
AreaErr() / peakMixed.
Area(), 2)));
429 rawChiSquares->SetPoint(i, angles->AverageAngle(i), peak.
GetReducedChi2());
430 mixedChiSquares->SetPoint(i, angles->AverageAngle(i), peakMixed.
GetReducedChi2());
432 std::cout << std::setw(3) << i <<
" of " << angles->NumberOfAngles() <<
" done\r" << std::flush;
434 std::cout <<
"Fitting of projections done." << std::endl;
437 auto* rawAngularDistributionCorr =
static_cast<TGraphErrors*
>(rawAngularDistribution->Clone(
"RawAngularDistributionCorrected"));
438 for(
int i = 0; i < rawAngularDistributionCorr->GetN(); ++i) {
439#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
440 if(angles->Count(rawAngularDistributionCorr->GetPointX(i)) != 0) {
441 rawAngularDistributionCorr->SetPointY(i, rawAngularDistributionCorr->GetPointY(i) / angles->Count(rawAngularDistributionCorr->GetPointX(i)));
442 rawAngularDistributionCorr->SetPointError(i, rawAngularDistributionCorr->GetErrorX(i), rawAngularDistributionCorr->GetErrorY(i) / angles->Count(rawAngularDistributionCorr->GetPointX(i)));
444 rawAngularDistributionCorr->SetPointY(i, 0);
448 rawAngularDistributionCorr->GetPoint(i, px, py);
449 if(angles->Count(px) != 0) {
450 rawAngularDistributionCorr->SetPoint(i, px, py / angles->Count(px));
451 rawAngularDistributionCorr->SetPointError(i, rawAngularDistributionCorr->GetErrorX(i), rawAngularDistributionCorr->GetErrorY(i) / angles->Count(px));
453 rawAngularDistributionCorr->SetPoint(i, px, 0);
457 auto* mixedAngularDistributionCorr =
static_cast<TGraphErrors*
>(mixedAngularDistribution->Clone(
"MixedAngularDistributionCorrected"));
458 for(
int i = 0; i < mixedAngularDistributionCorr->GetN(); ++i) {
459#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
460 if(angles->Count(mixedAngularDistributionCorr->GetPointX(i)) != 0) {
461 mixedAngularDistributionCorr->SetPointY(i, mixedAngularDistributionCorr->GetPointY(i) / angles->Count(mixedAngularDistributionCorr->GetPointX(i)));
462 mixedAngularDistributionCorr->SetPointError(i, mixedAngularDistributionCorr->GetErrorX(i), mixedAngularDistributionCorr->GetErrorY(i) / angles->Count(mixedAngularDistributionCorr->GetPointX(i)));
464 mixedAngularDistributionCorr->SetPointY(i, 0);
468 mixedAngularDistributionCorr->GetPoint(i, px, py);
469 if(angles->Count(px) != 0) {
470 mixedAngularDistributionCorr->SetPoint(i, px, py / angles->Count(px));
471 mixedAngularDistributionCorr->SetPointError(i, mixedAngularDistributionCorr->GetErrorX(i), mixedAngularDistributionCorr->GetErrorY(i) / angles->Count(px));
473 mixedAngularDistributionCorr->SetPoint(i, px, 0);
484 if(!theoryFile.empty()) {
485 TFile theory(theoryFile.c_str());
487 if(theory.IsOpen()) {
489 auto* z0 =
static_cast<TGraphErrors*
>(theory.Get(
"graph000"));
490 auto* z2 =
static_cast<TGraphErrors*
>(theory.Get(
"graph010"));
491 auto* z4 =
static_cast<TGraphErrors*
>(theory.Get(
"graph100"));
493 if(z0 !=
nullptr && z2 !=
nullptr && z4 !=
nullptr && z0->GetN() == z2->GetN() && z0->GetN() == z4->GetN()) {
497 if(z0->GetN() == angularDistribution->GetN() || ((angles->Grouping() || angles->Folding()) && z0->GetN() == (angles->Addback() ? 49 : 51))) {
499 if(z0->GetN() != angularDistribution->GetN()) {
500 angles->FoldOrGroup(z0, z2, z4);
503 std::vector<TGraph*> spin;
504 std::vector<double> spinLabel;
505 std::vector<std::vector<double>> parameters;
507 if(twoJLow.size() > 1 && twoJMiddle.size() == 1 && twoJHigh.size() == 1) {
508 for(
auto twoJ : twoJLow) {
509 parameters.emplace_back();
510 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJMiddle.at(0), twoJ, parameters.back()));
511 spinLabel.push_back(twoJ / 2.);
513 }
else if(twoJLow.size() == 1 && twoJMiddle.size() > 1 && twoJHigh.size() == 1) {
514 for(
auto twoJ : twoJMiddle) {
515 parameters.emplace_back();
516 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJ, twoJLow.at(0), parameters.back()));
517 spinLabel.push_back(twoJ / 2.);
519 }
else if(twoJLow.size() == 1 && twoJMiddle.size() == 1 && twoJHigh.size() > 1) {
520 for(
auto twoJ : twoJHigh) {
521 parameters.emplace_back();
522 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJ, twoJMiddle.at(0), twoJLow.at(0), parameters.back()));
523 spinLabel.push_back(twoJ / 2.);
526 parameters.emplace_back();
527 spin.push_back(
MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJMiddle.at(0), twoJLow.at(0), parameters.back()));
528 spinLabel.push_back(twoJHigh.at(0) / 2.);
532 auto* canvas =
new TCanvas;
535 double min = TMath::MinElement(spin.at(0)->GetN(), spin.at(0)->GetY());
536 double max = TMath::MaxElement(spin.at(0)->GetN(), spin.at(0)->GetY());
537 for(
size_t i = 1; i < spin.size(); ++i) {
538 min = TMath::Min(min, TMath::MinElement(spin.at(i)->GetN(), spin.at(i)->GetY()));
539 max = TMath::Max(max, TMath::MaxElement(spin.at(i)->GetN(), spin.at(i)->GetY()));
541 min = TMath::Min(min, confidenceLevel);
545 for(first = 0; first < spin.size(); ++first) {
546 if(spin[first]->GetN() > 1) {
break; }
549 spin[first]->SetTitle(
"");
550 spin[first]->SetMinimum(0.9 * min);
551 spin[first]->SetMaximum(1.1 * max);
553 for(
size_t i = 0; i < spin.size(); ++i) {
554 spin[i]->SetLineColor(i + 1);
555 spin[i]->SetMarkerColor(i + 1);
556 spin[i]->SetLineWidth(2);
559 spin[first]->Draw(
"ac");
560 for(
size_t i = 0; i < spin.size(); ++i) {
561 if(i == first) {
continue; }
562 if(spin[i]->GetN() > 1) {
569 auto* confidenceLevelLine =
new TLine(-1.5, confidenceLevel, 1.5, confidenceLevel);
571 confidenceLevelLine->Draw();
573#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
574 auto* legend =
new TLegend(0.1, 0.3);
576 auto* legend =
new TLegend(0.7, 0.6, 0.8, 0.9);
578 for(
size_t i = 0; i < spin.size(); ++i) {
579 if(spin[i]->GetN() == 1) {
580 legend->AddEntry(spin[i], Form(
"J = %.1f", spinLabel[i]),
"p");
582 legend->AddEntry(spin[i], Form(
"J = %.1f", spinLabel[i]),
"l");
590 spin[first]->GetHistogram()->GetXaxis()->SetRangeUser(-1.5, 1.5);
591 spin[first]->GetHistogram()->GetXaxis()->SetTitle(
"atan(#delta) [rad]");
592 spin[first]->GetHistogram()->GetXaxis()->CenterTitle();
593 spin[first]->GetHistogram()->GetYaxis()->SetTitle(
"red. #chi^{2}");
594 spin[first]->GetHistogram()->GetYaxis()->CenterTitle();
598 z0->Write(
"graph000");
599 z2->Write(
"graph010");
600 z4->Write(
"graph100");
601 for(
size_t i = 0; i < spin.size(); ++i) {
602 spin[i]->Write(Form(
"spin%d",
static_cast<int>(i)));
604 canvas->Write(
"MixingCanvas");
607 std::vector<TGraphErrors*> spinFit(spin.size(),
new TGraphErrors(angularDistribution->GetN()));
608 auto* x = angularDistribution->GetX();
609 for(
int p = 0; p < angularDistribution->GetN(); ++p) {
610 for(
size_t i = 0; i < spinFit.size(); ++i) {
611 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])));
612 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))));
615 for(
size_t i = 0; i < spin.size(); ++i) {
616 spinFit[i]->SetLineColor(i + 1);
617 spinFit[i]->SetMarkerColor(i + 1);
618 spinFit[i]->SetLineWidth(2);
619 spinFit[i]->Write(Form(
"SpinFit%d",
static_cast<int>(i)));
620 output.WriteObject(¶meters[i], Form(
"Parameters%d",
static_cast<int>(i)));
623 auto a2a4Parameters =
A2a4Method(angularDistribution, z0, z2, z4);
624 output.WriteObject(&a2a4Parameters,
"ParametersA2a4Fit");
626 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;
629 std::cerr <<
"Failed to find z0 (" << z0 <<
", \"graph000\"), z2 (" << z2 <<
", \"graph010\"), or z4 (" << z4 <<
", \"graph100\") in " << theoryFile <<
", or they had mismatched sizes" << std::endl;
632 std::cerr <<
"Failed to open " << theoryFile << std::endl;
635 std::cout <<
"No file with simulation results (--theory flag), so we won't produce chi2 vs mixing angle plot." << std::endl;
639 rawAngularDistribution->Write();
640 angularDistribution->Write();
641 mixedAngularDistribution->Write();
642 rawChiSquares->Write();
643 mixedChiSquares->Write();
644 rawAngularDistributionCorr->Write();
645 mixedAngularDistributionCorr->Write();
824std::vector<double>
A2a4Method(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4)
826 assert(data->GetN() == z0->GetN());
827 assert(data->GetN() == z2->GetN());
828 assert(data->GetN() == z4->GetN());
831 auto* cosTheta =
new TGraphErrors(*data);
832 auto* x = cosTheta->GetX();
833 for(
int i = 0; i < cosTheta->GetN(); ++i) {
834 x[i] = TMath::Cos(x[i] / 180. * TMath::Pi());
839 legendre->SetParNames(
"a_{0}",
"a_{2}",
"a_{4}");
840 legendre->SetParameters(1., 0.5, 0.5);
841 cosTheta->Fit(legendre,
"QN0");
844 Ac ac(data, z0, z2, z4);
845 ROOT::Fit::Fitter fitter;
846 fitter.SetFCN(nPar, ac);
848 fitter.Config().ParSettings(0) = ROOT::Fit::ParameterSettings(
"a_{0}", data->GetMaximum() / z0->GetMaximum(), 0.0001, -10., 10.);
849 for(
int i = 1; i < nPar; ++i) {
851 fitter.Config().ParSettings(i) = ROOT::Fit::ParameterSettings(Form(
"a_{%d}", 2 * i), legendre->GetParameter(i), 0.0001, -10., 10.);
853 fitter.Config().MinimizerOptions().SetPrintLevel(0);
854 fitter.Config().SetMinimizer(
"Minuit2",
"Migrad");
855 if(!fitter.FitFCN()) {
856 std::cerr <<
"Fit failed" << std::endl;
861 auto fitResult = fitter.Result();
863 std::cout <<
"Reduced chi^2: " << fitResult.MinFcnValue() / (ac.
Np() - fitResult.NFreeParameters()) << std::endl;
864 std::vector<double> parameters(fitResult.GetParams(), fitResult.GetParams() + nPar);
865 const auto* errors = fitResult.GetErrors();
866 std::cout <<
"Parameters a_0: " << parameters[0] <<
" +- " << errors[0] <<
", a_2: " << parameters[1] <<
" +- " << errors[1] <<
", a_4: " << parameters[2] <<
" +- " << errors[2] << std::endl;
867 TMatrixD covariance(nPar, nPar);
868 fitResult.GetCovarianceMatrix(covariance);
871 auto* fit =
static_cast<TGraphErrors*
>(z0->Clone(
"fit"));
872 for(
int i = 0; i < fit->GetN(); ++i) {
873#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
874 fit->SetPointY(i, parameters[0] * ((1. - parameters[1] - parameters[2]) * z0->GetPointY(i) + parameters[1] * z2->GetPointY(i) + parameters[2] * z4->GetPointY(i)));
876 fit->SetPoint(i, fit->GetX()[i], parameters[0] * ((1. - parameters[1] - parameters[2]) * z0->GetY()[i] + parameters[1] * z2->GetY()[i] + parameters[2] * z4->GetY()[i]));
878 fit->SetPointError(i, 0., std::abs(parameters[0]) * TMath::Sqrt(TMath::Power((1. - parameters[1] - parameters[2]) * z0->GetErrorY(i), 2) + TMath::Power(parameters[1] * z2->GetErrorY(i), 2) + TMath::Power(parameters[2] * z4->GetErrorY(i), 2)));
881 auto* residual =
static_cast<TGraphErrors*
>(data->Clone(
"residual"));
882 for(
int i = 0; i < residual->GetN(); ++i) {
883#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
884 residual->SetPointY(i, data->GetPointY(i) - fit->GetPointY(i));
886 residual->SetPoint(i, residual->GetX()[i], data->GetY()[i] - fit->GetY()[i]);
888 residual->SetPointError(i, 0., TMath::Sqrt(TMath::Power(data->GetErrorY(i), 2) + TMath::Power(fit->GetErrorY(i), 2)));
893 auto* stats =
new TPaveText(0.35, 0.7, 0.75, 0.95,
"NDC");
894 stats->SetTextFont(133);
895 stats->SetTextSize(20);
896 stats->SetFillStyle(0);
897 stats->SetBorderSize(0);
898 stats->AddText(Form(
"a_{2} = %f #pm %f", parameters[1], errors[1]));
899 stats->AddText(Form(
"a_{4} = %f #pm %f", parameters[2], errors[2]));
900 stats->AddText(Form(
"#chi^{2}/NDF = %.2f", fitResult.MinFcnValue() / (ac.
Np() - fitResult.NFreeParameters())));
905 auto* canvas =
new TCanvas;
906 auto* resPad =
new TPad(
"resPad",
"resPad", 0., 0., 1., 0.3);
907 resPad->SetTopMargin(0.);
908 resPad->SetBottomMargin(0.22);
909 resPad->SetLeftMargin(0.1);
910 resPad->SetRightMargin(0.01);
912 auto* compPad =
new TPad(
"compPad",
"compPad", 0., 0.3, 1., 1.);
913 compPad->SetTopMargin(0.01);
914 compPad->SetBottomMargin(0.);
915 compPad->SetLeftMargin(0.1);
916 compPad->SetRightMargin(0.01);
922 auto* multiGraph =
new TMultiGraph;
923 fit->SetLineColor(kRed);
924 fit->SetFillColor(kRed);
925 fit->SetMarkerColor(kRed);
926 multiGraph->Add(fit,
"l3");
927 data->SetMarkerStyle(kFullDotLarge);
928 multiGraph->Add(data,
"p");
930 multiGraph->SetTitle(
";;Normalized Counts;");
932 multiGraph->GetXaxis()->SetRangeUser(0., 180.);
934 multiGraph->GetYaxis()->CenterTitle();
935 multiGraph->GetYaxis()->SetTitleSize(0.05);
936 multiGraph->GetYaxis()->SetTitleOffset(1.);
938 multiGraph->Draw(
"a");
944 residual->SetTitle(
";#vartheta [^{o}];Residual");
946 residual->GetXaxis()->CenterTitle();
947 residual->GetXaxis()->SetTitleSize(0.1);
948 residual->GetXaxis()->SetLabelSize(0.1);
949 residual->GetXaxis()->SetRangeUser(0., 180.);
951 residual->GetYaxis()->CenterTitle();
952 residual->GetYaxis()->SetTitleSize(0.1);
953 residual->GetYaxis()->SetTitleOffset(0.5);
954 residual->GetYaxis()->SetLabelSize(0.08);
956 residual->Draw(
"ap");
957 auto* zeroLine =
new TLine(0., 0., 180., 0.);
958 zeroLine->Draw(
"same");
960 fit->Write(
"A2a4Fit");
961 residual->Write(
"Residual");
962 multiGraph->Write(
"FitComparison");
963 canvas->Write(
"A2a4Canvas");