GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
AngularCorrelations.cxx
Go to the documentation of this file.
1#include <iostream>
2#include <iomanip>
3#include <vector>
4#include <string>
5#include <cassert>
6
7#include "TFile.h"
8#include "TH2.h"
9#include "TH1.h"
10#include "TGraphErrors.h"
11#include "TMultiGraph.h"
12#include "TCanvas.h"
13#include "TLine.h"
14#include "TLegend.h"
15#include "Fit/Fitter.h"
16#include "TMatrixD.h"
17#include "TPaveText.h"
18
19#include "ArgParser.h"
20#include "TUserSettings.h"
21#include "TGriffinAngles.h"
22#include "TPeakFitter.h"
23#include "TRWPeak.h"
24
25TGraph* MixingMethod(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4, int twoJhigh, int twoJmid, int twoJlow, std::vector<double>& bestParameters);
26std::vector<double> A2a4Method(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4);
27
28double GetYError(TGraphErrors* graph, const double& x)
29{
30 /// general function to get the error of a graph at point x (takes the maximum of the errors of the bracketing points)
31 for(int i = 0; i < graph->GetN(); ++i) {
32#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
33 // exact match: return the error at this point
34 if(graph->GetPointX(i) == x) { return graph->GetErrorY(i); }
35 // first point with larger x: take maximum of this point and previous point
36 // TGraphErrors::GetErrorY returns -1 if index is negative, so we don't need to check for this
37 if(graph->GetPointX(i) > x) { return std::max(graph->GetErrorY(i - 1), graph->GetErrorY(i)); }
38#else
39 double px, py;
40 graph->GetPoint(i, px, py);
41 // exact match: return the error at this point
42 if(px == x) { return graph->GetErrorY(i); }
43 // first point with larger x: take maximum of this point and previous point
44 // TGraphErrors::GetErrorY returns -1 if index is negative, so we don't need to check for this
45 if(px > x) { return std::max(graph->GetErrorY(i - 1), graph->GetErrorY(i)); }
46#endif
47 }
48 return 0.;
49}
50
51/// program to read in 2D matrices from AngularCorrelationHelper
52/// project and fit peaks to create angular correlation plots
53/// and to create chi-square plots
54int main(int argc, char** argv)
55{
56 // --------------------------------------------------------------------------------
57 // Reading and verifying inputs and settings.
58 // --------------------------------------------------------------------------------
59
60 // input parameters
61 bool help = false;
62 double projGateLow = -1.;
63 double projGateHigh = -1.;
64 double projBgLow = -1.;
65 double projBgHigh = -1.;
66 double timeRandomNorm = -1.;
67 double peakPos = -1.;
68 double peakLow = -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;
75
76 // set up the argument parser
77 ArgParser parser;
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\".");
92
93 parser.parse(argc, argv, true);
94
95 if(help) {
96 std::cout << parser << std::endl;
97 return 1;
98 }
99
100 // open the input root file and read any settings stored there, then add the ones potentially read in from command line
101 if(inputFile.empty()) {
102 std::cerr << "Need an input file!" << std::endl;
103 std::cout << parser << std::endl;
104 return 1;
105 }
106
107 TFile input(inputFile.c_str());
108
109 if(!input.IsOpen()) {
110 std::cerr << "Failed to open input file " << inputFile << std::endl;
111 std::cout << parser << std::endl;
112 return 1;
113 }
114
115 auto* settings = static_cast<TUserSettings*>(input.Get("UserSettings"));
116
117 // if we have a path for a settings file provided on command line, we either add it to the ones read
118 // from the root file or (if there weren't any) create a new instance from it
119 if(!settingsFile.empty()) {
120 if(settings == nullptr) {
121 settings = new TUserSettings(settingsFile);
122 } else {
123 settings->ReadSettings(settingsFile);
124 }
125 }
126
127 // check that we got settings either from the root file or a path provided on command line
128 if(settings == nullptr || settings->empty()) {
129 std::cerr << "Failed to get user settings from input file." << std::endl;
130 std::cout << parser << std::endl;
131 return 1;
132 }
133
134 // set all variables from settings if they haven't been set from command line
135 if(projGateLow == -1.) { projGateLow = settings->GetDouble("Projection.Low"); }
136 if(projGateHigh == -1.) { projGateHigh = settings->GetDouble("Projection.High"); }
137 if(timeRandomNorm == -1.) { timeRandomNorm = settings->GetDouble("TimeRandomNormalization"); }
138 if(peakPos == -1.) { peakPos = settings->GetDouble("Peak.Position"); }
139 if(peakLow == -1.) { peakLow = settings->GetDouble("Peak.Low"); }
140 if(peakHigh == -1.) { peakHigh = settings->GetDouble("Peak.High"); }
141 // for the background-peak positions and background gates we could have multiple, so we create vectors for them
142 std::vector<double> bgPeakPos;
143 std::vector<double> bgLow;
144 std::vector<double> bgHigh;
145 // we only fill the background gates from the settings if there were none provided on the command line
146 if(projBgLow == -1. || projBgHigh == -1. || projBgLow >= projBgHigh) {
147 for(int i = 1;; ++i) {
148 try {
149 auto low = settings->GetDouble(Form("Background.%d.Low", i), true);
150 auto high = settings->GetDouble(Form("Background.%d.High", i), true);
151 if(low >= high) {
152 std::cout << i << ". background gate, low edge not lower than the high edge: " << low << " >= " << high << std::endl;
153 break;
154 }
155 bgLow.push_back(low);
156 bgHigh.push_back(high);
157 } catch(std::out_of_range& e) {
158 break;
159 }
160 }
161 } else {
162 bgLow.push_back(projBgLow);
163 bgHigh.push_back(projBgHigh);
164 }
165 // background-peak positions can only be set via settings file
166 // we loop until we fail to find an entry
167 for(int i = 1;; ++i) {
168 try {
169 auto pos = settings->GetDouble(Form("BackgroundPeak.Position.%d", i), true);
170 if(pos <= peakLow || pos >= peakHigh) {
171 std::cout << i << ". background peak outside of fit range: " << pos << " <= " << peakLow << " or " << pos << " >= " << peakHigh << std::endl;
172 break;
173 }
174 bgPeakPos.push_back(pos);
175 } catch(std::out_of_range& e) {
176 break;
177 }
178 }
179
180 // parameter limits and fixed parameters for the peak and the background peaks
181 // if the limits are the same, the parameter is fixed to that value, if the high limit is lower than the low limit there is no limit
182 // right now we are hard coded to use TRWPeaks which have 6 parameters
183 std::vector<double> peakParameter(6, -2.);
184 std::vector<double> peakParameterLow(6, 0.);
185 std::vector<double> peakParameterHigh(6, -1.);
186 std::vector<double> bgPeakParameter(6, -2.);
187 std::vector<double> bgPeakParameterLow(6, 0.);
188 std::vector<double> bgPeakParameterHigh(6, -1.);
189 for(size_t i = 0; i < peakParameterLow.size(); ++i) {
190 peakParameter[i] = settings->GetDouble(Form("Peak.Parameter.%d", static_cast<int>(i)), -2.);
191 peakParameterLow[i] = settings->GetDouble(Form("Peak.Parameter.%d.Low", static_cast<int>(i)), 0.);
192 peakParameterHigh[i] = settings->GetDouble(Form("Peak.Parameter.%d.High", static_cast<int>(i)), -1.);
193 bgPeakParameter[i] = settings->GetDouble(Form("BackgroundPeak.Parameter.%d", static_cast<int>(i)), -2.);
194 bgPeakParameterLow[i] = settings->GetDouble(Form("BackgroundPeak.Parameter.%d.Low", static_cast<int>(i)), 0.);
195 bgPeakParameterHigh[i] = settings->GetDouble(Form("BackgroundPeak.Parameter.%d.High", static_cast<int>(i)), -1.);
196 // check that the result makes sense, i.e. if the limits are in the righ order that the parameter itself is within the limits
197 // only output a warning that the parameter is changed if it's not the default value
198 if(peakParameterLow[i] <= peakParameterHigh[i] && (peakParameter[i] < peakParameterLow[i] || peakParameterHigh[i] < peakParameter[i])) {
199 bool output = peakParameter[i] != -2.;
200 if(output) { std::cout << "Warning, " << i << ". peak parameter (" << peakParameter[i] << ") is out of range " << peakParameterLow[i] << " - " << peakParameterHigh[i] << ", resetting it to "; }
201 peakParameter[i] = (peakParameterHigh[i] + peakParameterLow[i]) / 2.;
202 if(output) { std::cout << peakParameter[i] << std::endl; }
203 }
204 if(bgPeakParameterLow[i] <= bgPeakParameterHigh[i] && (bgPeakParameter[i] < bgPeakParameterLow[i] || bgPeakParameterHigh[i] < bgPeakParameter[i])) {
205 bool output = bgPeakParameter[i] != -2.;
206 if(output) { std::cout << "Warning, " << i << ". background peak parameter (" << bgPeakParameter[i] << ") is out of range " << bgPeakParameterLow[i] << " - " << bgPeakParameterHigh[i] << ", resetting it to "; }
207 bgPeakParameter[i] = (bgPeakParameterHigh[i] + bgPeakParameterLow[i]) / 2.;
208 if(output) { std::cout << bgPeakParameter[i] << std::endl; }
209 }
210 }
211 // parameter limits for the background (A + B*(x-o) + C*(x-o)^2)
212 // same idea as above, lower limit = higher limit means fixed parameter, higher limit < lower limit means parameter not set
213 std::vector<double> backgroundParameter(4, -2.);
214 std::vector<double> backgroundParameterLow(4, 0.);
215 std::vector<double> backgroundParameterHigh(4, -1.);
216 backgroundParameter[0] = settings->GetDouble("Background.Offset", -2.);
217 backgroundParameterLow[0] = settings->GetDouble("Background.Offset.Low", 0.);
218 backgroundParameterHigh[0] = settings->GetDouble("Background.Offset.High", -1.);
219 backgroundParameter[1] = settings->GetDouble("Background.Linear", -2.);
220 backgroundParameterLow[1] = settings->GetDouble("Background.Linear.Low", 0.);
221 backgroundParameterHigh[1] = settings->GetDouble("Background.Linear.High", -1.);
222 backgroundParameter[2] = settings->GetDouble("Background.Quadratic", -2.);
223 backgroundParameterLow[2] = settings->GetDouble("Background.Quadratic.Low", 0.);
224 backgroundParameterHigh[2] = settings->GetDouble("Background.Quadratic.High", -1.);
225 backgroundParameter[3] = settings->GetDouble("Background.Xoffset", -2.);
226 backgroundParameterLow[3] = settings->GetDouble("Background.Xoffset.Low", 0.);
227 backgroundParameterHigh[3] = settings->GetDouble("Background.Xoffset.High", -1.);
228 for(size_t i = 0; i < backgroundParameter.size(); ++i) {
229 if(backgroundParameterLow[i] <= backgroundParameterHigh[i] && (backgroundParameter[i] < backgroundParameterLow[i] || backgroundParameterHigh[i] < backgroundParameter[i])) {
230 bool output = backgroundParameter[i] != -2.;
231 if(output) { std::cout << "Warning, " << i << ". background parameter (" << backgroundParameter[i] << ") is out of range " << backgroundParameterLow[i] << " - " << backgroundParameterHigh[i] << ", resetting it to "; }
232 backgroundParameter[i] = (backgroundParameterHigh[i] + backgroundParameterLow[i]) / 2.;
233 if(output) { std::cout << backgroundParameter[i] << std::endl; }
234 }
235 }
236
237 // the spins of the low, middle, and high levels, to be used for the mixing method
238 // two of these need to be vectors of length one (meaning the settings file should have an entry with "name: <value>,"), the third can have a length larger than 1
239 std::vector<int> twoJLow = settings->GetIntVector("TwoJ.Low");
240 std::vector<int> twoJMiddle = settings->GetIntVector("TwoJ.Middle");
241 std::vector<int> twoJHigh = settings->GetIntVector("TwoJ.High");
242
243 // confidence level (this value is used to draw a line at the confidence level for the mixing method)
244 double confidenceLevel = settings->GetDouble("ConfidenceLevel", 1.535); // 1.535 is 99% confidence level for 48 degrees of freedom
245
246 // check if all necessary settings have been provided
247 if(projGateLow >= projGateHigh) {
248 std::cerr << "Need a projection gate with a low edge that is smaller than the high edge, " << projGateLow << " >= " << projGateHigh << std::endl;
249 std::cout << parser << std::endl;
250 return 1;
251 }
252
253 if(peakLow >= peakHigh) {
254 std::cerr << "Need a fit range with a low edge that is smaller than the high edge, " << peakLow << " >= " << peakHigh << std::endl;
255 std::cout << parser << std::endl;
256 return 1;
257 }
258
259 if(peakPos >= peakHigh || peakPos <= peakLow) {
260 std::cerr << "Need a peak within the fit range, " << peakPos << " not within " << peakLow << " - " << peakHigh << std::endl;
261 std::cout << parser << std::endl;
262 return 1;
263 }
264
265 if(timeRandomNorm <= 0.) {
266 std::cerr << "Need a positive normalization factor for time random subtraction" << std::endl;
267 std::cout << parser << std::endl;
268 return 1;
269 }
270
271 // for the background gates we checked that the low edge is below the high edge before adding them to the vector
272 // so we only need to check that the vectors aren't empty (sizes should always be the same, but we check anyway)
273 if(bgLow.empty() || bgLow.size() != bgHigh.size()) {
274 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;
275 std::cout << parser << std::endl;
276 return 1;
277 }
278
279 // get the angles from the input file
280 auto* angles = static_cast<TGriffinAngles*>(input.Get("GriffinAngles"));
281
282 if(angles == nullptr) {
283 std::cerr << "Failed to find 'GriffinAngles' in '" << inputFile << "'" << std::endl;
284 std::cout << parser << std::endl;
285 return 1;
286 }
287
288 // --------------------------------------------------------------------------------
289 // Create the angular distribution.
290 // --------------------------------------------------------------------------------
291
292 // open output file and create graphs
293 TFile output(outputFile.c_str(), "recreate");
294
295 auto* rawAngularDistribution = new TGraphErrors(angles->NumberOfAngles());
296 rawAngularDistribution->SetName("RawAngularDistribution");
297 auto* angularDistribution = new TGraphErrors(angles->NumberOfAngles());
298 angularDistribution->SetName("AngularDistribution");
299 auto* mixedAngularDistribution = new TGraphErrors(angles->NumberOfAngles());
300 mixedAngularDistribution->SetName("MixedAngularDistribution");
301 auto* rawChiSquares = new TGraph(angles->NumberOfAngles());
302 rawChiSquares->SetName("RawChiSquares");
303 auto* mixedChiSquares = new TGraph(angles->NumberOfAngles());
304 mixedChiSquares->SetName("MixedChiSquares");
305
306 // write the user settings to the output file
307 settings->Write();
308
309#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
310 auto* fitDir = output.mkdir("fits", "Projections with fits", true);
311#else
312 auto* fitDir = output.mkdir("fits", "Projections with fits");
313#endif
314
315 // loop over all matrices
316 for(int i = 0; i < angles->NumberOfAngles(); ++i) {
317 // get the three histograms we need: prompt, time random, and event mixed
318 auto* prompt = static_cast<TH2*>(input.Get(Form("AngularCorrelation%d", i)));
319 if(prompt == nullptr) {
320 std::cerr << "Failed to find histogram '" << Form("AngularCorrelation%d", i) << "', should have " << angles->NumberOfAngles() << " angles in total!" << std::endl;
321 std::cout << parser << std::endl;
322 return 1;
323 }
324 auto* bg = static_cast<TH2*>(input.Get(Form("AngularCorrelationBG%d", i)));
325 if(bg == nullptr) {
326 std::cerr << "Failed to find histogram '" << Form("AngularCorrelationBG%d", i) << "', should have " << angles->NumberOfAngles() << " angles in total!" << std::endl;
327 std::cout << parser << std::endl;
328 return 1;
329 }
330 auto* mixed = static_cast<TH2*>(input.Get(Form("AngularCorrelationMixed%d", i)));
331 if(mixed == nullptr) {
332 std::cerr << "Failed to find histogram '" << Form("AngularCorrelationMixed%d", i) << "', should have " << angles->NumberOfAngles() << " angles in total!" << std::endl;
333 std::cout << parser << std::endl;
334 return 1;
335 }
336
337 // first subtract time random background from prompt data and enable proper error propagation for prompt and mixed data
338 prompt->Sumw2();
339 prompt->Add(bg, -timeRandomNorm);
340 mixed->Sumw2();
341
342 // project onto x-axis
343 auto* proj = prompt->ProjectionX(Form("proj%d", i), prompt->GetYaxis()->FindBin(projGateLow), prompt->GetYaxis()->FindBin(projGateHigh));
344 auto* projMixed = mixed->ProjectionX(Form("projMixed%d", i), mixed->GetYaxis()->FindBin(projGateLow), mixed->GetYaxis()->FindBin(projGateHigh));
345 // project background gate(s) onto x-axis
346 auto* projBg = prompt->ProjectionX(Form("projBg%d", i), prompt->GetYaxis()->FindBin(bgLow[0]), prompt->GetYaxis()->FindBin(bgHigh[0]));
347 auto* projMixedBg = mixed->ProjectionX(Form("projMixedBg%d", i), mixed->GetYaxis()->FindBin(bgLow[0]), mixed->GetYaxis()->FindBin(bgHigh[0]));
348 double bgGateWidth = prompt->GetYaxis()->FindBin(bgHigh[0]) - prompt->GetYaxis()->FindBin(bgLow[0]) + 1;
349 double mixedBgGateWidth = mixed->GetYaxis()->FindBin(bgHigh[0]) - mixed->GetYaxis()->FindBin(bgLow[0]) + 1;
350 for(size_t g = 1; g < bgLow.size(); ++g) {
351 projBg->Add(prompt->ProjectionX(Form("projBg%d", i), prompt->GetYaxis()->FindBin(bgLow[g]), prompt->GetYaxis()->FindBin(bgHigh[g])));
352 projMixedBg->Add(mixed->ProjectionX(Form("projMixedBg%d", i), mixed->GetYaxis()->FindBin(bgLow[g]), mixed->GetYaxis()->FindBin(bgHigh[g])));
353 bgGateWidth += prompt->GetYaxis()->FindBin(bgHigh[g]) - prompt->GetYaxis()->FindBin(bgLow[g]) + 1;
354 mixedBgGateWidth += mixed->GetYaxis()->FindBin(bgHigh[g]) - mixed->GetYaxis()->FindBin(bgLow[g]) + 1;
355 }
356
357 // subtract background gate (+1 because the projection includes first and last bin)
358 proj->Add(projBg, -(prompt->GetYaxis()->FindBin(projGateHigh) - prompt->GetYaxis()->FindBin(projGateLow) + 1) / bgGateWidth);
359 projMixed->Add(projMixedBg, -(mixed->GetYaxis()->FindBin(projGateHigh) - mixed->GetYaxis()->FindBin(projGateLow) + 1) / mixedBgGateWidth);
360
361 // fit the projections, we create separate peak fitters and peaks for the prompt and mixed histograms
362 TPeakFitter pf(peakLow, peakHigh);
363 TRWPeak peak(peakPos);
364 for(size_t p = 0; p < peakParameterLow.size(); ++p) {
365 if(peakParameterLow[p] == peakParameterHigh[p]) {
366 peak.GetFitFunction()->FixParameter(p, peakParameter[p]);
367 } else if(peakParameterLow[p] < peakParameterHigh[p]) {
368 peak.GetFitFunction()->SetParameter(p, peakParameter[p]);
369 peak.GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
370 }
371 }
372 pf.AddPeak(&peak);
373 for(auto bgPeak : bgPeakPos) {
374 auto* bgP = new TRWPeak(bgPeak);
375 for(size_t p = 0; p < bgPeakParameterLow.size(); ++p) {
376 if(bgPeakParameterLow[p] == bgPeakParameterHigh[p]) {
377 bgP->GetFitFunction()->FixParameter(p, bgPeakParameter[p]);
378 } else if(bgPeakParameterLow[p] < bgPeakParameterHigh[p]) {
379 bgP->GetFitFunction()->SetParameter(p, bgPeakParameter[p]);
380 bgP->GetFitFunction()->SetParLimits(p, bgPeakParameterLow[p], bgPeakParameterHigh[p]);
381 }
382 }
383 pf.AddPeak(bgP);
384 }
385 for(size_t p = 0; p < backgroundParameterLow.size(); ++p) {
386 if(backgroundParameterLow[p] == backgroundParameterHigh[p]) {
387 pf.GetBackground()->FixParameter(p, backgroundParameter[p]);
388 } else if(backgroundParameterLow[p] < backgroundParameterHigh[p]) {
389 pf.GetBackground()->SetParameter(p, backgroundParameter[p]);
390 pf.GetBackground()->SetParLimits(p, backgroundParameterLow[p], backgroundParameterHigh[p]);
391 }
392 }
393 pf.Fit(proj, "qretryfit");
394
395 TPeakFitter pfMixed(peakLow, peakHigh);
396 TRWPeak peakMixed(peakPos);
397 for(size_t p = 0; p < peakParameterLow.size(); ++p) {
398 if(peakParameterLow[p] == peakParameterHigh[p]) {
399 peakMixed.GetFitFunction()->FixParameter(p, peakParameter[p]);
400 } else if(peakParameterLow[p] < peakParameterHigh[p]) {
401 peakMixed.GetFitFunction()->SetParameter(p, peakParameter[p]);
402 peakMixed.GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
403 }
404 }
405 pfMixed.AddPeak(&peakMixed);
406 for(auto bgPeak : bgPeakPos) {
407 auto* bgP = new TRWPeak(bgPeak);
408 for(size_t p = 0; p < bgPeakParameterLow.size(); ++p) {
409 if(bgPeakParameterLow[p] == bgPeakParameterHigh[p]) {
410 bgP->GetFitFunction()->FixParameter(p, bgPeakParameter[p]);
411 } else if(bgPeakParameterLow[p] < bgPeakParameterHigh[p]) {
412 bgP->GetFitFunction()->SetParameter(p, bgPeakParameter[p]);
413 bgP->GetFitFunction()->SetParLimits(p, bgPeakParameterLow[p], bgPeakParameterHigh[p]);
414 }
415 }
416 pfMixed.AddPeak(bgP);
417 }
418 for(size_t p = 0; p < backgroundParameterLow.size(); ++p) {
419 if(backgroundParameterLow[p] == backgroundParameterHigh[p]) {
420 pfMixed.GetBackground()->FixParameter(p, backgroundParameter[p]);
421 } else if(backgroundParameterLow[p] < backgroundParameterHigh[p]) {
422 pfMixed.GetBackground()->SetParameter(p, backgroundParameter[p]);
423 pfMixed.GetBackground()->SetParLimits(p, backgroundParameterLow[p], backgroundParameterHigh[p]);
424 }
425 }
426 pfMixed.Fit(projMixed, "qretryfit");
427
428 // save the fitted histograms to the output file and the areas of the peaks
429 fitDir->cd();
430 proj->Write();
431 projMixed->Write();
432
433 // TODO: set an error for the angles?
434 rawAngularDistribution->SetPoint(i, angles->AverageAngle(i), peak.Area());
435 rawAngularDistribution->SetPointError(i, 0., peak.AreaErr());
436 mixedAngularDistribution->SetPoint(i, angles->AverageAngle(i), peakMixed.Area());
437 mixedAngularDistribution->SetPointError(i, 0., peakMixed.AreaErr());
438 angularDistribution->SetPoint(i, angles->AverageAngle(i), peak.Area() / peakMixed.Area());
439 angularDistribution->SetPointError(i, 0., peak.Area() / peakMixed.Area() * TMath::Sqrt(TMath::Power(peak.AreaErr() / peak.Area(), 2) + TMath::Power(peakMixed.AreaErr() / peakMixed.Area(), 2)));
440
441 rawChiSquares->SetPoint(i, angles->AverageAngle(i), peak.GetReducedChi2());
442 mixedChiSquares->SetPoint(i, angles->AverageAngle(i), peakMixed.GetReducedChi2());
443
444 std::cout << std::setw(3) << i << " of " << angles->NumberOfAngles() << " done\r" << std::flush;
445 }
446 std::cout << "Fitting of projections done." << std::endl;
447
448 // correct the raw and mixed graphs for the number of combinations that contribute to each angle
449 auto* rawAngularDistributionCorr = static_cast<TGraphErrors*>(rawAngularDistribution->Clone("RawAngularDistributionCorrected"));
450 for(int i = 0; i < rawAngularDistributionCorr->GetN(); ++i) {
451#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
452 if(angles->Count(rawAngularDistributionCorr->GetPointX(i)) != 0) {
453 rawAngularDistributionCorr->SetPointY(i, rawAngularDistributionCorr->GetPointY(i) / angles->Count(rawAngularDistributionCorr->GetPointX(i)));
454 rawAngularDistributionCorr->SetPointError(i, rawAngularDistributionCorr->GetErrorX(i), rawAngularDistributionCorr->GetErrorY(i) / angles->Count(rawAngularDistributionCorr->GetPointX(i)));
455 } else {
456 rawAngularDistributionCorr->SetPointY(i, 0);
457 }
458#else
459 double px, py;
460 rawAngularDistributionCorr->GetPoint(i, px, py);
461 if(angles->Count(px) != 0) {
462 rawAngularDistributionCorr->SetPoint(i, px, py / angles->Count(px));
463 rawAngularDistributionCorr->SetPointError(i, rawAngularDistributionCorr->GetErrorX(i), rawAngularDistributionCorr->GetErrorY(i) / angles->Count(px));
464 } else {
465 rawAngularDistributionCorr->SetPoint(i, px, 0);
466 }
467#endif
468 }
469 auto* mixedAngularDistributionCorr = static_cast<TGraphErrors*>(mixedAngularDistribution->Clone("MixedAngularDistributionCorrected"));
470 for(int i = 0; i < mixedAngularDistributionCorr->GetN(); ++i) {
471#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
472 if(angles->Count(mixedAngularDistributionCorr->GetPointX(i)) != 0) {
473 mixedAngularDistributionCorr->SetPointY(i, mixedAngularDistributionCorr->GetPointY(i) / angles->Count(mixedAngularDistributionCorr->GetPointX(i)));
474 mixedAngularDistributionCorr->SetPointError(i, mixedAngularDistributionCorr->GetErrorX(i), mixedAngularDistributionCorr->GetErrorY(i) / angles->Count(mixedAngularDistributionCorr->GetPointX(i)));
475 } else {
476 mixedAngularDistributionCorr->SetPointY(i, 0);
477 }
478#else
479 double px, py;
480 mixedAngularDistributionCorr->GetPoint(i, px, py);
481 if(angles->Count(px) != 0) {
482 mixedAngularDistributionCorr->SetPoint(i, px, py / angles->Count(px));
483 mixedAngularDistributionCorr->SetPointError(i, mixedAngularDistributionCorr->GetErrorX(i), mixedAngularDistributionCorr->GetErrorY(i) / angles->Count(px));
484 } else {
485 mixedAngularDistributionCorr->SetPoint(i, px, 0);
486 }
487#endif
488 }
489
490 // --------------------------------------------------------------------------------
491 // If a theory/simulation file has been provided, we use the a2/a4 and mixing
492 // methods to fit the angular distribution.
493 // --------------------------------------------------------------------------------
494
495 // check if we have a theory file
496 if(!theoryFile.empty()) {
497 TFile theory(theoryFile.c_str());
498
499 if(theory.IsOpen()) {
500 // read graphs from file
501 auto* z0 = static_cast<TGraphErrors*>(theory.Get("graph000"));
502 auto* z2 = static_cast<TGraphErrors*>(theory.Get("graph010"));
503 auto* z4 = static_cast<TGraphErrors*>(theory.Get("graph100"));
504
505 if(z0 != nullptr && z2 != nullptr && z4 != nullptr && z0->GetN() == z2->GetN() && z0->GetN() == z4->GetN()) {
506 // check if the sizes of the provided graphs match what we expect:
507 // if we have folding and grouping enabled we either need the graphs to match the angular distribution or have a size of either 51 (singles) or 49 (addback)
508 // otherwise the size needs to match the angular distribution
509 if(z0->GetN() == angularDistribution->GetN() || ((angles->Grouping() || angles->Folding()) && z0->GetN() == (angles->Addback() ? 49 : 51))) {
510 // if the theory needs to be folded and/or grouped, do so now
511 if(z0->GetN() != angularDistribution->GetN()) {
512 angles->FoldOrGroup(z0, z2, z4);
513 }
514 // calculate chi2 vs mixing graphs
515 std::vector<TGraph*> spin;
516 std::vector<double> spinLabel;
517 std::vector<std::vector<double>> parameters;
518 // first check which of the vectors we iterate over
519 if(twoJLow.size() > 1 && twoJMiddle.size() == 1 && twoJHigh.size() == 1) {
520 for(auto twoJ : twoJLow) {
521 parameters.emplace_back();
522 spin.push_back(MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJMiddle.at(0), twoJ, parameters.back()));
523 spinLabel.push_back(twoJ / 2.);
524 }
525 } else if(twoJLow.size() == 1 && twoJMiddle.size() > 1 && twoJHigh.size() == 1) {
526 for(auto twoJ : twoJMiddle) {
527 parameters.emplace_back();
528 spin.push_back(MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJ, twoJLow.at(0), parameters.back()));
529 spinLabel.push_back(twoJ / 2.);
530 }
531 } else if(twoJLow.size() == 1 && twoJMiddle.size() == 1 && twoJHigh.size() > 1) {
532 for(auto twoJ : twoJHigh) {
533 parameters.emplace_back();
534 spin.push_back(MixingMethod(angularDistribution, z0, z2, z4, twoJ, twoJMiddle.at(0), twoJLow.at(0), parameters.back()));
535 spinLabel.push_back(twoJ / 2.);
536 }
537 } else {
538 parameters.emplace_back();
539 spin.push_back(MixingMethod(angularDistribution, z0, z2, z4, twoJHigh.at(0), twoJMiddle.at(0), twoJLow.at(0), parameters.back()));
540 spinLabel.push_back(twoJHigh.at(0) / 2.);
541 }
542
543 // create canvas and plot graphs on it
544 auto* canvas = new TCanvas;
545
546 // determine minimum and maximum y-value
547 double min = TMath::MinElement(spin.at(0)->GetN(), spin.at(0)->GetY());
548 double max = TMath::MaxElement(spin.at(0)->GetN(), spin.at(0)->GetY());
549 for(size_t i = 1; i < spin.size(); ++i) {
550 min = TMath::Min(min, TMath::MinElement(spin.at(i)->GetN(), spin.at(i)->GetY()));
551 max = TMath::Max(max, TMath::MaxElement(spin.at(i)->GetN(), spin.at(i)->GetY()));
552 }
553 min = TMath::Min(min, confidenceLevel);
554
555 // find first graph with more than one data point
556 size_t first = 0;
557 for(first = 0; first < spin.size(); ++first) {
558 if(spin[first]->GetN() > 1) { break; }
559 }
560
561 spin[first]->SetTitle("");
562 spin[first]->SetMinimum(0.9 * min);
563 spin[first]->SetMaximum(1.1 * max);
564
565 for(size_t i = 0; i < spin.size(); ++i) {
566 spin[i]->SetLineColor(i + 1);
567 spin[i]->SetMarkerColor(i + 1);
568 spin[i]->SetLineWidth(2);
569 }
570
571 spin[first]->Draw("ac");
572 for(size_t i = 0; i < spin.size(); ++i) {
573 if(i == first) { continue; }
574 if(spin[i]->GetN() > 1) {
575 spin[i]->Draw("c");
576 } else {
577 spin[i]->Draw("*");
578 }
579 }
580
581 auto* confidenceLevelLine = new TLine(-1.5, confidenceLevel, 1.5, confidenceLevel);
582
583 confidenceLevelLine->Draw();
584
585#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
586 auto* legend = new TLegend(0.1, 0.3);
587#else
588 auto* legend = new TLegend(0.7, 0.6, 0.8, 0.9);
589#endif
590 for(size_t i = 0; i < spin.size(); ++i) {
591 if(spin[i]->GetN() == 1) {
592 legend->AddEntry(spin[i], Form("J = %.1f", spinLabel[i]), "p");
593 } else {
594 legend->AddEntry(spin[i], Form("J = %.1f", spinLabel[i]), "l");
595 }
596 }
597
598 legend->Draw();
599
600 canvas->SetLogy();
601
602 spin[first]->GetHistogram()->GetXaxis()->SetRangeUser(-1.5, 1.5);
603 spin[first]->GetHistogram()->GetXaxis()->SetTitle("atan(#delta) [rad]");
604 spin[first]->GetHistogram()->GetXaxis()->CenterTitle();
605 spin[first]->GetHistogram()->GetYaxis()->SetTitle("red. #chi^{2}");
606 spin[first]->GetHistogram()->GetYaxis()->CenterTitle();
607
608 // write graphs and canvas to output file
609 output.cd();
610 z0->Write("graph000");
611 z2->Write("graph010");
612 z4->Write("graph100");
613 for(size_t i = 0; i < spin.size(); ++i) {
614 spin[i]->Write(Form("spin%d", static_cast<int>(i)));
615 }
616 canvas->Write("MixingCanvas");
617
618 // create theory graphs with best fit for each spin, and write them to file
619 std::vector<TGraphErrors*> spinFit(spin.size(), new TGraphErrors(angularDistribution->GetN()));
620 auto* x = angularDistribution->GetX();
621 for(int p = 0; p < angularDistribution->GetN(); ++p) {
622 for(size_t i = 0; i < spinFit.size(); ++i) {
623 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])));
624 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))));
625 }
626 }
627 for(size_t i = 0; i < spin.size(); ++i) {
628 spinFit[i]->SetLineColor(i + 1);
629 spinFit[i]->SetMarkerColor(i + 1);
630 spinFit[i]->SetLineWidth(2);
631 spinFit[i]->Write(Form("SpinFit%d", static_cast<int>(i)));
632 output.WriteObject(&parameters[i], Form("Parameters%d", static_cast<int>(i)));
633 }
634
635 auto a2a4Parameters = A2a4Method(angularDistribution, z0, z2, z4);
636 output.WriteObject(&a2a4Parameters, "ParametersA2a4Fit");
637 } else { // if(z0->GetN() == angularDistribution->GetN() || ((angles->Grouping() || angles->Folding) && z0->GetN() == (angles->Addback() ? 49:51)))
638 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;
639 }
640 } else { // if(z0 != nullptr && z2 != nullptr && z4 != nullptr && z0->GetN() == z2->GetN() && z0->GetN() == z4->GetN())
641 std::cerr << "Failed to find z0 (" << z0 << ", \"graph000\"), z2 (" << z2 << ", \"graph010\"), or z4 (" << z4 << ", \"graph100\") in " << theoryFile << ", or they had mismatched sizes" << std::endl;
642 }
643 } else { // if(theory.IsOpen())
644 std::cerr << "Failed to open " << theoryFile << std::endl;
645 }
646 } else { // if(!theoryFile.empty())
647 std::cout << "No file with simulation results (--theory flag), so we won't produce chi2 vs mixing angle plot." << std::endl;
648 }
649 output.cd();
650
651 rawAngularDistribution->Write();
652 angularDistribution->Write();
653 mixedAngularDistribution->Write();
654 rawChiSquares->Write();
655 mixedChiSquares->Write();
656 rawAngularDistributionCorr->Write();
657 mixedAngularDistributionCorr->Write();
658
659 angles->Write();
660
661 output.Close();
662 input.Close();
663
664 return 0;
665}
666
667class Ac {
668public:
669 explicit Ac(TGraphErrors* data = nullptr, TGraphErrors* z0 = nullptr, TGraphErrors* z2 = nullptr, TGraphErrors* z4 = nullptr)
670 : fData(data), fZ0(z0), fZ2(z2), fZ4(z4) {}
671
672 void Data(TGraphErrors* data) { fData = data; }
673 void Z0(TGraphErrors* z0) { fZ0 = z0; }
674 void Z2(TGraphErrors* z2) { fZ2 = z2; }
675 void Z4(TGraphErrors* z4) { fZ4 = z4; }
676 void SetZ(TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4)
677 {
678 fZ0 = z0;
679 fZ2 = z2;
680 fZ4 = z4;
681 }
682
683 int Np() { return fData->GetN(); }
684
685 double operator()(const double* p)
686 {
687 double chi2 = 0;
688 for(int point = 0; point < fData->GetN(); ++point) {
689 // get data values
690 double x = 0.;
691 double y = 0.;
692 fData->GetPoint(point, x, y);
693 double yError = fData->GetErrorY(point);
694
695 // get simulation values
696 // for the y-values we can use Eval (uses linear interpolation)
697 double functionValue = p[0] * ((1. - p[1] - p[2]) * fZ0->Eval(x) + p[1] * fZ2->Eval(x) + p[2] * fZ4->Eval(x));
698 // for the y uncertainties we need to find the index (for each graph)
699 double errorSquare = TMath::Power(yError, 2) + TMath::Power(p[0], 2) * (TMath::Power((1. - p[1] - p[2]) * GetYError(fZ0, x), 2) + TMath::Power(p[1] * GetYError(fZ2, x), 2) + TMath::Power(p[2] * GetYError(fZ4, x), 2));
700
701 // calculate chi^2
702 chi2 += TMath::Power((y - functionValue), 2) / errorSquare;
703 }
704 return chi2;
705 }
706
707private:
708 TGraphErrors* fData{nullptr};
709 TGraphErrors* fZ0{nullptr};
710 TGraphErrors* fZ2{nullptr};
711 TGraphErrors* fZ4{nullptr};
712};
713
714TGraph* MixingMethod(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4, int twoJhigh, int twoJmid, int twoJlow, std::vector<double>& bestParameters)
715{
716 TGraph* result = nullptr;
717 Ac ac(data, z0, z2, z4);
718 ROOT::Fit::Fitter fitter;
719 int nPar = 3;
720 fitter.SetFCN(nPar, ac);
721 for(int i = 0; i < nPar; ++i) {
722 // parameter settings arguments are parameter name, initial value, step size, minimum, and maximum
723 fitter.Config().ParSettings(i) = ROOT::Fit::ParameterSettings(Form("a_{%d}", 2 * i), 0.5, 0.0001, -10., 10.);
724 }
725 fitter.Config().MinimizerOptions().SetPrintLevel(0);
726 fitter.Config().SetMinimizer("Minuit2", "Migrad"); // or simplex?
727
728 // j1 is the spin of the highest level
729 // j2 is the spin of the middle level
730 // j3 is the spin of the bottom level
731 double j1 = 0.5 * twoJhigh;
732 double j2 = 0.5 * twoJmid;
733 double j3 = 0.5 * twoJlow;
734 // l1 is the transition between j1 and j2
735 // a is the lowest allowed spin
736 // b is the mixing spin
737 int l1a = TMath::Abs(twoJhigh - twoJmid) / 2;
738 if(l1a == 0) { l1a = 1; }
739 int l1b = l1a + 1;
740 // l2 is the transition between j2 and j3
741 // a is the lowest allowed spin
742 // b is the mixing spin
743 int l2a = TMath::Abs(twoJmid - twoJlow) / 2;
744 if(l2a == 0) { l2a = 1; }
745 int l2b = l2a + 1;
746
747 // run some quick checks on the mixing ratios
748 if((twoJhigh == 0 && twoJmid == 0) || (twoJmid == 0 && twoJlow == 0)) {
749 std::cout << "!!!!!!!!!!!!!!! ERROR !!!!!!!!!!!!!!!" << std::endl
750 << "Can't have gamma transition between J=0 states (high " << twoJhigh << ", mid " << twoJmid << ", low " << twoJlow << ")." << std::endl
751 << "Aborting..." << std::endl;
752 return result;
753 }
754 if(l1a == TMath::Abs(twoJhigh + twoJmid) / 2) {
755 //std::cout<<"!!!!!!!!!!!!!!! ALERT !!!!!!!!!!!!!!!"<<std::endl
756 // <<"Only one angular momentum allowed for high->middle transition (l1a "<<l1a<<" == "<<TMath::Abs(twoJhigh+twoJmid)/2<<")."<<std::endl
757 // <<"That mixing ratio (delta1) will be fixed at zero."<<std::endl
758 // <<"!!!!!!!!!!!!! END ALERT !!!!!!!!!!!!!"<<std::endl;
759 l1b = l1a;
760 }
761 if(l2a == TMath::Abs(twoJmid + twoJlow) / 2) {
762 //std::cout<<"!!!!!!!!!!!!!!! ALERT !!!!!!!!!!!!!!!"<<std::endl
763 // <<"Only one angular momentum allowed for middle->low transition (l2a "<<l2a<<" == "<<TMath::Abs(twoJmid+twoJlow)/2<<")."<<std::endl
764 // <<"That mixing ratio (delta2) will be fixed at zero."<<std::endl
765 // <<"!!!!!!!!!!!!! END ALERT !!!!!!!!!!!!!"<<std::endl;
766 l2b = l2a;
767 }
768
769 // -------------------------------------------------------------------//
770 // Constrained fitting
771 // -------------------------------------------------------------------//
772 // The basic idea is to select a particular set of physical quantities,
773 // calculate a2/a4, fix the a2/a4 parameters, and fit the scaling
774 // factor a0. Then output the specifications for that set of physical
775 // quantities and the chi^2 for further analysis.
776 // -------------------------------------------------------------------//
777
778 // delta runs from -infinity to infinity (unless constrained by known physics)
779 // in this case, it then makes more sense to sample evenly from tan^{-1}(delta)
780 // The next few lines are where you may want to include limits to significantly speed up calculations
781 // mixing for the high-middle transition
782 double mixingAngle1Minimum = -TMath::Pi() / 2;
783 double mixingAngle1Maximum = TMath::Pi() / 2;
784 int steps1 = 100;
785 double stepSize1 = (mixingAngle1Maximum - mixingAngle1Minimum) / steps1;
786 // mixing for the middle-low transition
787 double mixingAngle2Minimum = -TMath::Pi() / 2;
788 double mixingAngle2Maximum = TMath::Pi() / 2;
789 int steps2 = 100;
790 double stepSize2 = (mixingAngle2Maximum - mixingAngle2Minimum) / steps2;
791
792 // if appropriate, constrain the delta values
793 if(l1a == l1b) {
794 mixingAngle1Minimum = 0;
795 steps1 = 1;
796 }
797 if(l2a == l2b) {
798 mixingAngle2Minimum = 0;
799 steps2 = 1;
800 }
801
802 result = new TGraph(steps1);
803 double minChi2 = 1e6;
804 for(int i = 0; i < steps1; i++) {
805 double mixangle1 = mixingAngle1Minimum + i * stepSize1;
806 double delta1 = TMath::Tan(mixangle1);
807 for(int j = 0; j < steps2; j++) {
808 double mixangle2 = mixingAngle2Minimum + j * stepSize2;
809 double delta2 = TMath::Tan(mixangle2);
810 // calculate a2
811 double a2 = TGRSIFunctions::CalculateA2(j1, j2, j3, l1a, l1b, l2a, l2b, delta1, delta2);
812 // fix a2
813 fitter.Config().ParSettings(1).Set("a_{2}", a2);
814 fitter.Config().ParSettings(1).Fix();
815 // calculate a4
816 double a4 = TGRSIFunctions::CalculateA4(j1, j2, j3, l1a, l1b, l2a, l2b, delta1, delta2);
817 // fix a4
818 fitter.Config().ParSettings(2).Set("a_{4}", a4);
819 fitter.Config().ParSettings(2).Fix();
820 if(!fitter.FitFCN()) {
821 std::cerr << i << ", " << j << ": Fit failed using a2 " << a2 << ", a4 " << a4 << std::endl;
822 continue;
823 }
824 auto fitResult = fitter.Result();
825 // MinFcnValue() is the minimum chi2, ac.Np gives the number of data points
826 result->SetPoint(i, mixangle1, fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters()));
827 if(fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters()) < minChi2) {
828 minChi2 = fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters());
829 bestParameters.assign(fitResult.GetParams(), fitResult.GetParams() + nPar);
830 }
831 }
832 }
833 return result;
834}
835
836std::vector<double> A2a4Method(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4)
837{
838 assert(data->GetN() == z0->GetN());
839 assert(data->GetN() == z2->GetN());
840 assert(data->GetN() == z4->GetN());
841 // create a copy of the data with cos(theta) as x-axis and fit it with a legenre polynomial
842 // this is to get initial conditions for our fit
843 auto* cosTheta = new TGraphErrors(*data);
844 auto* x = cosTheta->GetX();
845 for(int i = 0; i < cosTheta->GetN(); ++i) {
846 x[i] = TMath::Cos(x[i] / 180. * TMath::Pi());
847 }
848
849 int nPar = 3;
850 auto* legendre = new TF1("legendre", TGRSIFunctions::LegendrePolynomial, -1., 1., nPar);
851 legendre->SetParNames("a_{0}", "a_{2}", "a_{4}");
852 legendre->SetParameters(1., 0.5, 0.5);
853 cosTheta->Fit(legendre, "QN0");
854
855 // the actual fitting
856 Ac ac(data, z0, z2, z4);
857 ROOT::Fit::Fitter fitter;
858 fitter.SetFCN(nPar, ac);
859 // this is a good guess for the initial scale
860 fitter.Config().ParSettings(0) = ROOT::Fit::ParameterSettings("a_{0}", data->GetMaximum() / z0->GetMaximum(), 0.0001, -10., 10.);
861 for(int i = 1; i < nPar; ++i) {
862 // parameter settings arguments are parameter name, initial value, step size, minimum, and maximum
863 fitter.Config().ParSettings(i) = ROOT::Fit::ParameterSettings(Form("a_{%d}", 2 * i), legendre->GetParameter(i), 0.0001, -10., 10.);
864 }
865 fitter.Config().MinimizerOptions().SetPrintLevel(0);
866 fitter.Config().SetMinimizer("Minuit2", "Migrad"); // or simplex?
867 if(!fitter.FitFCN()) {
868 std::cerr << "Fit failed" << std::endl;
869 return {};
870 }
871
872 // get the fit result and print chi^2 and parameters
873 auto fitResult = fitter.Result();
874 // MinFcnValue() is the minimum chi2, ac.Np gives the number of data points
875 std::cout << "Reduced chi^2: " << fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters()) << std::endl;
876 std::vector<double> parameters(fitResult.GetParams(), fitResult.GetParams() + nPar);
877 const auto* errors = fitResult.GetErrors();
878 std::cout << "Parameters a_0: " << parameters[0] << " +- " << errors[0] << ", a_2: " << parameters[1] << " +- " << errors[1] << ", a_4: " << parameters[2] << " +- " << errors[2] << std::endl;
879 TMatrixD covariance(nPar, nPar);
880 fitResult.GetCovarianceMatrix(covariance);
881
882 // create fit and residual graphs
883 auto* fit = static_cast<TGraphErrors*>(z0->Clone("fit"));
884 for(int i = 0; i < fit->GetN(); ++i) {
885#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
886 fit->SetPointY(i, parameters[0] * ((1. - parameters[1] - parameters[2]) * z0->GetPointY(i) + parameters[1] * z2->GetPointY(i) + parameters[2] * z4->GetPointY(i)));
887#else
888 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]));
889#endif
890 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)));
891 }
892
893 auto* residual = static_cast<TGraphErrors*>(data->Clone("residual"));
894 for(int i = 0; i < residual->GetN(); ++i) {
895#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
896 residual->SetPointY(i, data->GetPointY(i) - fit->GetPointY(i));
897#else
898 residual->SetPoint(i, residual->GetX()[i], data->GetY()[i] - fit->GetY()[i]);
899#endif
900 residual->SetPointError(i, 0., TMath::Sqrt(TMath::Power(data->GetErrorY(i), 2) + TMath::Power(fit->GetErrorY(i), 2)));
901 }
902
903 // This text box will display the fit statistics
904 // with the margins of the pads, the center in x is at 0.55 (0.545 to be exact_, so we center around that point
905 auto* stats = new TPaveText(0.35, 0.7, 0.75, 0.95, "NDC");
906 stats->SetTextFont(133);
907 stats->SetTextSize(20);
908 stats->SetFillStyle(0);
909 stats->SetBorderSize(0);
910 stats->AddText(Form("a_{2} = %f #pm %f", parameters[1], errors[1]));
911 stats->AddText(Form("a_{4} = %f #pm %f", parameters[2], errors[2]));
912 stats->AddText(Form("#chi^{2}/NDF = %.2f", fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters())));
913
914 // create canvas and two pads (big one for comparison and small one for residuals)
915 // wider left margin for y-axis labels and title
916 // same for the bottom margin of the residuals
917 auto* canvas = new TCanvas;
918 auto* resPad = new TPad("resPad", "resPad", 0., 0., 1., 0.3);
919 resPad->SetTopMargin(0.);
920 resPad->SetBottomMargin(0.22);
921 resPad->SetLeftMargin(0.1);
922 resPad->SetRightMargin(0.01);
923 resPad->Draw();
924 auto* compPad = new TPad("compPad", "compPad", 0., 0.3, 1., 1.);
925 compPad->SetTopMargin(0.01);
926 compPad->SetBottomMargin(0.);
927 compPad->SetLeftMargin(0.1);
928 compPad->SetRightMargin(0.01);
929 compPad->Draw();
930
931 // plot comparison of fit and data
932 compPad->cd();
933
934 auto* multiGraph = new TMultiGraph;
935 fit->SetLineColor(kRed);
936 fit->SetFillColor(kRed);
937 fit->SetMarkerColor(kRed);
938 multiGraph->Add(fit, "l3"); // 3 = filled contour between upper and lower error bars
939 data->SetMarkerStyle(kFullDotLarge);
940 multiGraph->Add(data, "p");
941
942 multiGraph->SetTitle(";;Normalized Counts;");
943
944 multiGraph->GetXaxis()->SetRangeUser(0., 180.);
945
946 multiGraph->GetYaxis()->CenterTitle();
947 multiGraph->GetYaxis()->SetTitleSize(0.05);
948 multiGraph->GetYaxis()->SetTitleOffset(1.);
949
950 multiGraph->Draw("a");
951
952 stats->Draw();
953
954 // plot residuals
955 resPad->cd();
956 residual->SetTitle(";#vartheta [^{o}];Residual");
957
958 residual->GetXaxis()->CenterTitle();
959 residual->GetXaxis()->SetTitleSize(0.1);
960 residual->GetXaxis()->SetLabelSize(0.1);
961 residual->GetXaxis()->SetRangeUser(0., 180.);
962
963 residual->GetYaxis()->CenterTitle();
964 residual->GetYaxis()->SetTitleSize(0.1);
965 residual->GetYaxis()->SetTitleOffset(0.5);
966 residual->GetYaxis()->SetLabelSize(0.08);
967
968 residual->Draw("ap");
969 auto* zeroLine = new TLine(0., 0., 180., 0.);
970 zeroLine->Draw("same");
971
972 fit->Write("A2a4Fit");
973 residual->Write("Residual");
974 multiGraph->Write("FitComparison");
975 canvas->Write("A2a4Canvas");
976
977 return parameters;
978}
int main(int argc, char **argv)
double GetYError(TGraphErrors *graph, const double &x)
TGraph * MixingMethod(TGraphErrors *data, TGraphErrors *z0, TGraphErrors *z2, TGraphErrors *z4, int twoJhigh, int twoJmid, int twoJlow, std::vector< double > &bestParameters)
std::vector< double > A2a4Method(TGraphErrors *data, TGraphErrors *z0, TGraphErrors *z2, TGraphErrors *z4)
double bgLow
double bgHigh
void Z2(TGraphErrors *z2)
TGraphErrors * fData
void SetZ(TGraphErrors *z0, TGraphErrors *z2, TGraphErrors *z4)
void Z4(TGraphErrors *z4)
double operator()(const double *p)
void Z0(TGraphErrors *z0)
TGraphErrors * fZ4
TGraphErrors * fZ0
void Data(TGraphErrors *data)
Ac(TGraphErrors *data=nullptr, TGraphErrors *z0=nullptr, TGraphErrors *z2=nullptr, TGraphErrors *z4=nullptr)
TGraphErrors * fZ2
void parse(int argc, char **argv, bool firstPass)
Definition ArgParser.h:333
ArgParseConfigT< T > & option(const std::string flag, T *output_location, bool firstPass)
Definition ArgParser.h:412
TF1 * GetBackground() const
Definition TPeakFitter.h:76
void AddPeak(TSinglePeak *peak)
Definition TPeakFitter.h:41
TFitResultPtr Fit(TH1 *fit_hist, Option_t *opt="")
TF1 * GetFitFunction() const
Definition TSinglePeak.h:78
Double_t GetReducedChi2() const
Definition TSinglePeak.h:93
Double_t AreaErr() const
Definition TSinglePeak.h:58
Double_t Area() const
Definition TSinglePeak.h:57
double CalculateA2(double j1, double j2, double j3, double l1a, double l1b, double l2a, double l2b, double delta1, double delta2)
Double_t LegendrePolynomial(Double_t *x, Double_t *p)
double CalculateA4(double j1, double j2, double j3, double l1a, double l1b, double l2a, double l2b, double delta1, double delta2)