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#include <fstream>
7
8#include "TFile.h"
9#include "TH2.h"
10#include "TH1.h"
11#include "TGraphErrors.h"
12#include "TMultiGraph.h"
13#include "TCanvas.h"
14#include "TLine.h"
15#include "TLegend.h"
16#include "Fit/Fitter.h"
17#include "TMatrixD.h"
18#include "TPaveText.h"
19
20#include "ArgParser.h"
21#include "TUserSettings.h"
22#include "TGriffinAngles.h"
23#include "TPeakFitter.h"
24#include "TRWPeak.h"
25#include "TRedirect.h"
26#include "TGRSIFunctions.h"
27
28TGraph* MixingMethod(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4, int twoJhigh, int twoJmid, int twoJlow, std::vector<double>& bestParameters, std::ofstream& logFile);
29std::vector<double> A2a4Method(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4);
30
31double GetYError(TGraphErrors* graph, const double& x)
32{
33 /// general function to get the error of a graph at point x (takes the maximum of the errors of the bracketing points)
34 for(int i = 0; i < graph->GetN(); ++i) {
35#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
36 // exact match: return the error at this point
37 if(graph->GetPointX(i) == x) { return graph->GetErrorY(i); }
38 // first point with larger x: take maximum of this point and previous point
39 // TGraphErrors::GetErrorY returns -1 if index is negative, so we don't need to check for this
40 if(graph->GetPointX(i) > x) { return std::max(graph->GetErrorY(i - 1), graph->GetErrorY(i)); }
41#else
42 double px, py;
43 graph->GetPoint(i, px, py);
44 // exact match: return the error at this point
45 if(px == x) { return graph->GetErrorY(i); }
46 // first point with larger x: take maximum of this point and previous point
47 // TGraphErrors::GetErrorY returns -1 if index is negative, so we don't need to check for this
48 if(px > x) { return std::max(graph->GetErrorY(i - 1), graph->GetErrorY(i)); }
49#endif
50 }
51 return 0.;
52}
53
54/// program to read in 2D matrices from AngularCorrelationHelper
55/// project and fit peaks to create angular correlation plots
56/// and to create chi-square plots
57int main(int argc, char** argv)
58{
59 // --------------------------------------------------------------------------------
60 // Reading and verifying inputs and settings.
61 // --------------------------------------------------------------------------------
62
63 // input parameters
64 bool help = false;
65 double projGateLow = -1.;
66 double projGateHigh = -1.;
67 double projBgLow = -1.;
68 double projBgHigh = -1.;
69 double timeRandomNorm = -1.;
70 double peakPos = -1.;
71 double peakLow = -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;
80
81 // set up the argument parser
82 ArgParser parser;
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\".");
99
100 parser.parse(argc, argv, true);
101
102 if(help) {
103 std::cout << parser << std::endl;
104 return 1;
105 }
106
107 // open the input root file and read any settings stored there, then add the ones potentially read in from command line
108 if(inputFile.empty()) {
109 std::cerr << "Need an input file!" << std::endl;
110 std::cout << parser << std::endl;
111 return 1;
112 }
113
114 TFile input(inputFile.c_str());
115
116 if(!input.IsOpen()) {
117 std::cerr << "Failed to open input file " << inputFile << std::endl;
118 std::cout << parser << std::endl;
119 return 1;
120 }
121
122 auto* settings = static_cast<TUserSettings*>(input.Get("UserSettings"));
123
124 // if we have a path for a settings file provided on command line, we either add it to the ones read
125 // from the root file or (if there weren't any) create a new instance from it
126 if(!settingsFile.empty()) {
127 if(settings == nullptr) {
128 settings = new TUserSettings(settingsFile);
129 } else {
130 settings->ReadSettings(settingsFile);
131 }
132 }
133
134 // check that we got settings either from the root file or a path provided on command line
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;
138 return 1;
139 }
140
141 // set all variables from settings if they haven't been set from command line
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()) {
152 try {
153 theoryFile = settings->GetString("Theory", true);
154 } catch(std::out_of_range&) {}
155 }
156 if(outputFile == "AngularCorrelations.root") { outputFile = settings->GetString("Output", "AngularCorrelations.root"); }
157
158 // for the background-peak positions and background gates we could have multiple, so we create vectors for them
159 std::vector<std::tuple<double, double, double>> bgPeakPos; // position, low, high
160 std::vector<double> bgLow;
161 std::vector<double> bgHigh;
162 // we only fill the background gates from the settings if there were none provided on the command line
163 if(projBgLow == -1. || projBgHigh == -1. || projBgLow >= projBgHigh) {
164 for(int i = 1;; ++i) {
165 try {
166 auto low = settings->GetDouble(Form("Background.%d.Low", i), true);
167 auto high = settings->GetDouble(Form("Background.%d.High", i), true);
168 if(low >= high) {
169 std::cout << i << ". background gate, low edge not lower than the high edge: " << low << " >= " << high << std::endl;
170 break;
171 }
172 bgLow.push_back(low);
173 bgHigh.push_back(high);
174 } catch(std::out_of_range&) {
175 break;
176 }
177 }
178 } else {
179 bgLow.push_back(projBgLow);
180 bgHigh.push_back(projBgHigh);
181 }
182 // background-peak positions can only be set via settings file
183 // we loop until we fail to find an entry
184 for(int i = 1;; ++i) {
185 try {
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;
189 break;
190 }
191 // read low and high limit for this background peak, defaults to -1 if not set in the settings file
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&) {
196 break;
197 }
198 }
199
200 // parameter limits and fixed parameters for the peak and the background peaks
201 // 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
202 // right now we are hard coded to use TRWPeaks which have 6 parameters
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.);
216 // check that the result makes sense, i.e. if the limits are in the righ order that the parameter itself is within the limits
217 // only output a warning that the parameter is changed if it's not the default value
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; }
223 }
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; }
229 }
230 // if we don't have limits for the peak position, fix it
231 if(i == 1 && peakParameterLow[i] == 0. && peakParameterHigh[i] == -1.) {
232 peakParameter[i] = peakPos;
233 peakParameterLow[i] = peakPos;
234 peakParameterHigh[i] = peakPos;
235 }
236 }
237 // parameter limits for the background (A + B*(x-o) + C*(x-o)^2)
238 // same idea as above, lower limit = higher limit means fixed parameter, higher limit < lower limit means parameter not set
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; }
260 }
261 }
262
263 // the spins of the low, middle, and high levels, to be used for the mixing method
264 // 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
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");
268
269 // confidence level (this value is used to draw a line at the confidence level for the mixing method)
270 double confidenceLevel = settings->GetDouble("ConfidenceLevel", 1.535); // 1.535 is 99% confidence level for 48 degrees of freedom
271
272 // check if all necessary settings have been provided
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;
276 return 1;
277 }
278
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;
282 return 1;
283 }
284
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;
288 return 1;
289 }
290
291 if(timeRandomNorm <= 0.) {
292 std::cerr << "Need a positive normalization factor for time random subtraction" << std::endl;
293 std::cout << parser << std::endl;
294 return 1;
295 }
296
297 // for the background gates we checked that the low edge is below the high edge before adding them to the vector
298 // so we only need to check that the vectors aren't empty (sizes should always be the same, but we check anyway)
299 if(bgLow.empty() || bgLow.size() != bgHigh.size()) {
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;
302 return 1;
303 }
304
305 // get the angles from the input file
306 auto* angles = static_cast<TGriffinAngles*>(input.Get("GriffinAngles"));
307
308 if(angles == nullptr) {
309 std::cerr << "Failed to find 'GriffinAngles' in '" << inputFile << "'" << std::endl;
310 std::cout << parser << std::endl;
311 return 1;
312 }
313
314 // get the log file name from the output file name
315 auto logFileName = outputFile.substr(0, outputFile.find_last_of('.')) + ".log";
316 std::ofstream logFile(logFileName.c_str());
317
318 // --------------------------------------------------------------------------------
319 // Create the angular distribution.
320 // --------------------------------------------------------------------------------
321
322 // open output file and create graphs
323 TFile output(outputFile.c_str(), "recreate");
324
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");
335
336 // write the user settings to the output file
337 settings->Write();
338
339#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
340 auto* fitDir = output.mkdir("fits", "Projections with fits", true);
341#else
342 auto* fitDir = output.mkdir("fits", "Projections with fits");
343#endif
344
345 logFile << "# columns are:" << std::endl
346 << "#ID p/m centroid +- uncertainty area +- uncertainty FWHM +- uncertainty red. chi^2" << std::endl;
347
348 // loop over all matrices
349 for(int i = 0; i < angles->NumberOfAngles(); ++i) {
350 // get the three histograms we need: prompt, time random, and event mixed
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;
355 return 1;
356 }
357 auto* bg = static_cast<TH2*>(input.Get(Form("%s%d", bgName.c_str(), i)));
358 if(bg == nullptr) {
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;
361 return 1;
362 }
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;
367 return 1;
368 }
369
370 // first subtract time random background from prompt data and enable proper error propagation for prompt and mixed data
371 prompt->Sumw2();
372 prompt->Add(bg, -timeRandomNorm);
373 mixed->Sumw2();
374
375 // project onto x-axis
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));
378 // project background gate(s) onto x-axis
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;
388 }
389
390 // subtract background gate (+1 because the projection includes first and last bin)
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);
393
394 // fit the projections, we create separate peak fitters and peaks for the prompt and mixed histograms
395
396 TPeakFitter pf(peakLow, peakHigh);
397 TRWPeak peak(peakPos);
398 for(size_t p = 0; p < peakParameterLow.size(); ++p) {
399 if(peakParameterLow[p] == peakParameterHigh[p]) {
400 peak.GetFitFunction()->FixParameter(p, peakParameter[p]);
401 } else if(peakParameterLow[p] < peakParameterHigh[p]) {
402 peak.GetFitFunction()->SetParameter(p, peakParameter[p]);
403 peak.GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
404 }
405 }
406 pf.AddPeak(&peak);
407 for(auto bgPeak : bgPeakPos) {
408 auto* bgP = new TRWPeak(std::get<0>(bgPeak));
409 // if we have limits for the position of this peak, apply them
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));
412 }
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]);
419 }
420 }
421 pf.AddPeak(bgP);
422 }
423 for(size_t p = 0; p < backgroundParameterLow.size(); ++p) {
424 if(backgroundParameterLow[p] == backgroundParameterHigh[p]) {
425 pf.GetBackground()->FixParameter(p, backgroundParameter[p]);
426 } else if(backgroundParameterLow[p] < backgroundParameterHigh[p]) {
427 pf.GetBackground()->SetParameter(p, backgroundParameter[p]);
428 pf.GetBackground()->SetParLimits(p, backgroundParameterLow[p], backgroundParameterHigh[p]);
429 }
430 }
431 {
432 TRedirect redirect("/dev/null");
433 pf.Fit(proj, "qretryfit");
434 }
435
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() << " "
440 << std::setw(10) << peak.GetReducedChi2() << std::endl;
441
442 TPeakFitter pfMixed(peakLow, peakHigh);
443 TRWPeak peakMixed(peakPos);
444 for(size_t p = 0; p < peakParameterLow.size(); ++p) {
445 if(peakParameterLow[p] == peakParameterHigh[p]) {
446 peakMixed.GetFitFunction()->FixParameter(p, peakParameter[p]);
447 } else if(peakParameterLow[p] < peakParameterHigh[p]) {
448 peakMixed.GetFitFunction()->SetParameter(p, peakParameter[p]);
449 peakMixed.GetFitFunction()->SetParLimits(p, peakParameterLow[p], peakParameterHigh[p]);
450 }
451 }
452 pfMixed.AddPeak(&peakMixed);
453 for(auto bgPeak : bgPeakPos) {
454 auto* bgP = new TRWPeak(std::get<0>(bgPeak));
455 // if we have limits for the position of this peak, apply them
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));
458 }
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]);
465 }
466 }
467 pfMixed.AddPeak(bgP);
468 }
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]);
475 }
476 }
477 {
478 TRedirect redirect("/dev/null");
479 pfMixed.Fit(projMixed, "qretryfit");
480 }
481
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() << " "
486 << std::setw(8) << peakMixed.GetReducedChi2() << std::endl;
487
488 // save the fitted histograms to the output file and the areas of the peaks
489 fitDir->cd();
490 proj->Write();
491 projMixed->Write();
492
493 // TODO: set an error for the angles?
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)));
500
501 rawChiSquares->SetPoint(i, angles->AverageAngle(i), peak.GetReducedChi2());
502 mixedChiSquares->SetPoint(i, angles->AverageAngle(i), peakMixed.GetReducedChi2());
503
504 std::cout << "Angle " << std::setw(3) << i << " of " << angles->NumberOfAngles() << " done\r" << std::flush;
505 }
506 std::cout << "Fitting of projections done." << std::endl;
507
508 // correct the raw and mixed graphs for the number of combinations that contribute to each angle
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)));
515 } else {
516 rawAngularDistributionCorr->SetPointY(i, 0);
517 }
518#else
519 double px, py;
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));
524 } else {
525 rawAngularDistributionCorr->SetPoint(i, px, 0);
526 }
527#endif
528 }
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)));
535 } else {
536 mixedAngularDistributionCorr->SetPointY(i, 0);
537 }
538#else
539 double px, py;
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));
544 } else {
545 mixedAngularDistributionCorr->SetPoint(i, px, 0);
546 }
547#endif
548 }
549
550 // change to output file
551 output.cd();
552
553 // --------------------------------------------------------------------------------
554 // If a theory/simulation file has been provided, we use the a2/a4 and mixing
555 // methods to fit the angular distribution.
556 // --------------------------------------------------------------------------------
557
558 // check if we have a theory file
559 if(!theoryFile.empty()) {
560 TFile theory(theoryFile.c_str());
561
562 if(theory.IsOpen()) {
563 // read graphs from file
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"));
567
568 if(z0 != nullptr && z2 != nullptr && z4 != nullptr && z0->GetN() == z2->GetN() && z0->GetN() == z4->GetN()) {
569 // check if the sizes of the provided graphs match what we expect:
570 // 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)
571 // otherwise the size needs to match the angular distribution
572 if(z0->GetN() == angularDistribution->GetN() || ((angles->Grouping() || angles->Folding()) && z0->GetN() == (angles->Addback() ? 49 : 51))) {
573 // if the theory needs to be folded and/or grouped, do so now
574 if(z0->GetN() != angularDistribution->GetN()) {
575 angles->FoldOrGroup(z0, z2, z4);
576 }
577 // calculate chi2 vs mixing graphs
578 std::vector<TGraph*> spin;
579 std::vector<double> spinLabel;
580 std::vector<std::vector<double>> parameters;
581 logFile << std::endl;
582 // first check which of the vectors we iterate over
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();
587 output.cd();
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.);
590 }
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();
595 output.cd();
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.);
598 }
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();
603 output.cd();
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.);
606 }
607 } else {
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();
610 output.cd();
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.);
613 }
614
615 // create canvas and plot graphs on it
616 auto* canvas = new TCanvas;
617
618 // determine minimum and maximum y-value
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()));
624 }
625 min = TMath::Min(min, confidenceLevel);
626
627 // find first graph with more than one data point
628 size_t first = 0;
629 for(first = 0; first < spin.size(); ++first) {
630 if(spin[first]->GetN() > 1) { break; }
631 }
632
633 spin[first]->SetTitle("");
634 spin[first]->SetMinimum(0.9 * min);
635 spin[first]->SetMaximum(1.1 * max);
636
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);
641 }
642
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) {
647 spin[i]->Draw("c");
648 } else {
649 spin[i]->Draw("*");
650 }
651 }
652
653 auto* confidenceLevelLine = new TLine(-1.5, confidenceLevel, 1.5, confidenceLevel);
654
655 confidenceLevelLine->Draw();
656
657#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
658 auto* legend = new TLegend(0.1, 0.3);
659#else
660 auto* legend = new TLegend(0.7, 0.6, 0.8, 0.9);
661#endif
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");
665 } else {
666 legend->AddEntry(spin[i], Form("J = %.1f", spinLabel[i]), "l");
667 }
668 }
669
670 legend->Draw();
671
672 canvas->SetLogy();
673
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();
679
680 // write graphs and canvas to output file
681 output.cd();
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)));
687 }
688 canvas->Write("MixingCanvas");
689
690 // create theory graphs with best fit for each spin, and write them to file
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))));
697 }
698 }
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(&parameters[i], Form("Parameters%d", static_cast<int>(i)));
705 }
706
707 output.cd();
708 auto a2a4Parameters = A2a4Method(angularDistribution, z0, z2, z4);
709 output.WriteObject(&a2a4Parameters, "ParametersA2a4Fit");
710 } else { // if(z0->GetN() == angularDistribution->GetN() || ((angles->Grouping() || angles->Folding) && z0->GetN() == (angles->Addback() ? 49:51)))
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;
712 }
713 } else { // if(z0 != nullptr && z2 != nullptr && z4 != nullptr && z0->GetN() == z2->GetN() && z0->GetN() == z4->GetN())
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;
715 }
716 } else { // if(theory.IsOpen())
717 std::cerr << "Failed to open " << theoryFile << std::endl;
718 }
719 } else { // if(!theoryFile.empty())
720 std::cout << "No file with simulation results (--theory flag), so we won't produce chi2 vs mixing angle plot." << std::endl;
721 }
722
723 rawAngularDistribution->Write();
724 angularDistribution->Write();
725 mixedAngularDistribution->Write();
726 rawChiSquares->Write();
727 mixedChiSquares->Write();
728 rawAngularDistributionCorr->Write();
729 mixedAngularDistributionCorr->Write();
730
731 angles->Write();
732
733 output.Close();
734 input.Close();
735 logFile.close();
736
737 return 0;
738}
739
740class Ac {
741public:
742 explicit Ac(TGraphErrors* data = nullptr, TGraphErrors* z0 = nullptr, TGraphErrors* z2 = nullptr, TGraphErrors* z4 = nullptr)
743 : fData(data), fZ0(z0), fZ2(z2), fZ4(z4) {}
744
745 void Data(TGraphErrors* data) { fData = data; }
746 void Z0(TGraphErrors* z0) { fZ0 = z0; }
747 void Z2(TGraphErrors* z2) { fZ2 = z2; }
748 void Z4(TGraphErrors* z4) { fZ4 = z4; }
749 void SetZ(TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4)
750 {
751 fZ0 = z0;
752 fZ2 = z2;
753 fZ4 = z4;
754 }
755
756 int Np() { return fData->GetN(); }
757
758 double operator()(const double* p)
759 {
760 double chi2 = 0;
761 for(int point = 0; point < fData->GetN(); ++point) {
762 // get data values
763 double x = 0.;
764 double y = 0.;
765 fData->GetPoint(point, x, y);
766 double yError = fData->GetErrorY(point);
767
768 // get simulation values
769 // for the y-values we can use Eval (uses linear interpolation)
770 double functionValue = p[0] * ((1. - p[1] - p[2]) * fZ0->Eval(x) + p[1] * fZ2->Eval(x) + p[2] * fZ4->Eval(x));
771 // for the y uncertainties we need to find the index (for each graph)
772 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));
773
774 // calculate chi^2
775 chi2 += TMath::Power((y - functionValue), 2) / errorSquare;
776 }
777 return chi2;
778 }
779
780private:
781 TGraphErrors* fData{nullptr};
782 TGraphErrors* fZ0{nullptr};
783 TGraphErrors* fZ2{nullptr};
784 TGraphErrors* fZ4{nullptr};
785};
786
787TMultiGraph* PlotCanvas(TGraphErrors* data, TGraphErrors* fit, TGraphErrors* residual, const std::vector<double>& parameters, const std::vector<double>& errors, const double& redChiSquare, const char* extraText = nullptr)
788{
789 /// This function plots the data with the fit and residuals on the currently active canvas.
790
791 // suppress any error messages
792 TRedirect redirect("/dev/null");
793
794 // This text box will display the fit statistics
795 // 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
796 auto* stats = new TPaveText(0.35, 0.7, 0.75, 0.95, "NDC");
797 stats->SetTextFont(133);
798 stats->SetTextSize(20);
799 stats->SetFillStyle(0);
800 stats->SetBorderSize(0);
801 for(int p = 0; p < 3; ++p) {
802 if(errors[p] != 0.) {
803 stats->AddText(Form("a_{%d} = %f #pm %f", 2 * p, parameters[p], errors[p]));
804 } else {
805 stats->AddText(Form("a_{%d} = %f", 2 * p, parameters[p]));
806 }
807 }
808 stats->AddText(Form("#chi^{2}/NDF = %.2f", redChiSquare));
809 if(extraText != nullptr) {
810 stats->AddText(extraText);
811 }
812
813 // create canvas and two pads (big one for comparison and small one for residuals)
814 // wider left margin for y-axis labels and title
815 // same for the bottom margin of the residuals
816 auto* resPad = new TPad("resPad", "resPad", 0., 0., 1., 0.3);
817 resPad->SetTopMargin(0.);
818 resPad->SetBottomMargin(0.22);
819 resPad->SetLeftMargin(0.1);
820 resPad->SetRightMargin(0.01);
821 resPad->Draw();
822 auto* compPad = new TPad("compPad", "compPad", 0., 0.3, 1., 1.);
823 compPad->SetTopMargin(0.01);
824 compPad->SetBottomMargin(0.);
825 compPad->SetLeftMargin(0.1);
826 compPad->SetRightMargin(0.01);
827 compPad->Draw();
828
829 // plot comparison of fit and data
830 compPad->cd();
831
832 auto* multiGraph = new TMultiGraph;
833 fit->SetLineColor(kRed);
834 fit->SetFillColor(kRed);
835 fit->SetMarkerColor(kRed);
836 multiGraph->Add(fit, "l3"); // 3 = filled contour between upper and lower error bars
837 data->SetMarkerStyle(kFullDotLarge);
838 multiGraph->Add(data, "p");
839
840 multiGraph->SetTitle(";;Normalized Counts;");
841
842 multiGraph->GetXaxis()->SetRangeUser(0., 180.);
843
844 multiGraph->GetYaxis()->CenterTitle();
845 multiGraph->GetYaxis()->SetTitleSize(0.05);
846 multiGraph->GetYaxis()->SetTitleOffset(1.);
847
848 multiGraph->Draw("a");
849
850 stats->Draw();
851
852 // plot residuals
853 resPad->cd();
854 residual->SetTitle(";#vartheta [^{o}];Residual");
855
856 residual->GetXaxis()->CenterTitle();
857 residual->GetXaxis()->SetTitleSize(0.1);
858 residual->GetXaxis()->SetLabelSize(0.1);
859 residual->GetXaxis()->SetRangeUser(0., 180.);
860
861 residual->GetYaxis()->CenterTitle();
862 residual->GetYaxis()->SetTitleSize(0.1);
863 residual->GetYaxis()->SetTitleOffset(0.5);
864 residual->GetYaxis()->SetLabelSize(0.08);
865
866 residual->Draw("ap");
867 auto* zeroLine = new TLine(residual->GetXaxis()->GetXmin(), 0., residual->GetXaxis()->GetXmax(), 0.);
868 zeroLine->Draw("same");
869
870 return multiGraph;
871}
872
873TGraph* MixingMethod(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4, int twoJhigh, int twoJmid, int twoJlow, std::vector<double>& bestParameters, std::ofstream& logFile)
874{
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;
880 int nPar = 3;
881 fitter.SetFCN(nPar, ac);
882 for(int i = 0; i < nPar; ++i) {
883 // parameter settings arguments are parameter name, initial value, step size, minimum, and maximum
884 fitter.Config().ParSettings(i) = ROOT::Fit::ParameterSettings(Form("a_{%d}", 2 * i), 0.5, 0.0001, -10., 10.);
885 }
886 fitter.Config().MinimizerOptions().SetPrintLevel(0);
887 fitter.Config().SetMinimizer("Minuit2", "Migrad"); // or simplex?
888
889 // j1 is the spin of the highest level
890 // j2 is the spin of the middle level
891 // j3 is the spin of the bottom level
892 double j1 = 0.5 * twoJhigh;
893 double j2 = 0.5 * twoJmid;
894 double j3 = 0.5 * twoJlow;
895 // l1 is the transition between j1 and j2
896 // a is the lowest allowed spin
897 // b is the mixing spin
898 int l1a = TMath::Abs(twoJhigh - twoJmid) / 2;
899 if(l1a == 0) { l1a = 1; }
900 int l1b = l1a + 1;
901 // l2 is the transition between j2 and j3
902 // a is the lowest allowed spin
903 // b is the mixing spin
904 int l2a = TMath::Abs(twoJmid - twoJlow) / 2;
905 if(l2a == 0) { l2a = 1; }
906 int l2b = l2a + 1;
907
908 // run some quick checks on the mixing ratios
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;
913 return result;
914 }
915 if(l1a == TMath::Abs(twoJhigh + twoJmid) / 2) {
916 //std::cout<<"!!!!!!!!!!!!!!! ALERT !!!!!!!!!!!!!!!"<<std::endl
917 // <<"Only one angular momentum allowed for high->middle transition (l1a "<<l1a<<" == "<<TMath::Abs(twoJhigh+twoJmid)/2<<")."<<std::endl
918 // <<"That mixing ratio (delta1) will be fixed at zero."<<std::endl
919 // <<"!!!!!!!!!!!!! END ALERT !!!!!!!!!!!!!"<<std::endl;
920 l1b = l1a;
921 }
922 if(l2a == TMath::Abs(twoJmid + twoJlow) / 2) {
923 //std::cout<<"!!!!!!!!!!!!!!! ALERT !!!!!!!!!!!!!!!"<<std::endl
924 // <<"Only one angular momentum allowed for middle->low transition (l2a "<<l2a<<" == "<<TMath::Abs(twoJmid+twoJlow)/2<<")."<<std::endl
925 // <<"That mixing ratio (delta2) will be fixed at zero."<<std::endl
926 // <<"!!!!!!!!!!!!! END ALERT !!!!!!!!!!!!!"<<std::endl;
927 l2b = l2a;
928 }
929
930 // -------------------------------------------------------------------//
931 // Constrained fitting
932 // -------------------------------------------------------------------//
933 // The basic idea is to select a particular set of physical quantities,
934 // calculate a2/a4, fix the a2/a4 parameters, and fit the scaling
935 // factor a0. Then output the specifications for that set of physical
936 // quantities and the chi^2 for further analysis.
937 // -------------------------------------------------------------------//
938
939 // delta runs from -infinity to infinity (unless constrained by known physics)
940 // in this case, it then makes more sense to sample evenly from tan^{-1}(delta)
941
942 // mixing for the high-middle transition
943 double mixingAngle1Minimum = -TMath::Pi() / 2;
944 double mixingAngle1Maximum = TMath::Pi() / 2;
945 int steps1 = 100;
946 double stepSize1 = (mixingAngle1Maximum - mixingAngle1Minimum) / steps1;
947 // mixing for the middle-low transition
948 double mixingAngle2Minimum = -TMath::Pi() / 2;
949 double mixingAngle2Maximum = TMath::Pi() / 2;
950 int steps2 = 100;
951 double stepSize2 = (mixingAngle2Maximum - mixingAngle2Minimum) / steps2;
952
953 // if appropriate, constrain the delta values
954 if(l1a == l1b) {
955 mixingAngle1Minimum = 0;
956 steps1 = 1;
957 }
958 if(l2a == l2b) {
959 mixingAngle2Minimum = 0;
960 steps2 = 1;
961 }
962
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);
974 // calculate a2
975 double a2 = TGRSIFunctions::CalculateA2(j1, j2, j3, l1a, l1b, l2a, l2b, delta1, delta2);
976 // fix a2
977 fitter.Config().ParSettings(1).Set("a_{2}", a2);
978 fitter.Config().ParSettings(1).Fix();
979 // calculate a4
980 double a4 = TGRSIFunctions::CalculateA4(j1, j2, j3, l1a, l1b, l2a, l2b, delta1, delta2);
981 // fix a4
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;
986 continue;
987 }
988 auto fitResult = fitter.Result();
989 // MinFcnValue() is the minimum chi2, ac.Np gives the number of data points
990 double chi2 = fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters());
991 // is it correct to always plot vs mixangle1? Or should we use mixangle2 if mixangle1 had only 1 step?
992 // and what if both angles have multiple steps?
993 result->SetPoint(i, mixangle1, chi2);
994 if(chi2 < minChi2) {
995 minChi2 = chi2;
996 bestParameters = fitResult.Parameters();
997 bestErrors = fitResult.Errors();
998 bestMixingAngle1 = mixangle1;
999 bestMixingAngle2 = mixangle2;
1000 }
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;
1002 }
1003 }
1004
1005 // find mixing ratio with minimum chi2 and its uncertainty
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);
1017 break;
1018 }
1019 }
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.);
1023 }
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);
1031 break;
1032 }
1033 }
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.);
1037 }
1038
1039 // plot the best fit, if there are multiple minima this is most likely the first one found
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)));
1044#else
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]));
1046#endif
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)));
1048 }
1049
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));
1054#else
1055 residual->SetPoint(i, residual->GetX()[i], data->GetY()[i] - fit->GetY()[i]);
1056#endif
1057 residual->SetPointError(i, 0., TMath::Sqrt(TMath::Power(data->GetErrorY(i), 2) + TMath::Power(fit->GetErrorY(i), 2)));
1058 }
1059
1060 auto* canvas = new TCanvas;
1061 TMultiGraph* multiGraph = nullptr;
1062
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);
1066 } else {
1067 if(steps1 == 1) {
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));
1073 } else {
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));
1076 }
1077 }
1078
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));
1083
1084 return result;
1085}
1086
1087std::vector<double> A2a4Method(TGraphErrors* data, TGraphErrors* z0, TGraphErrors* z2, TGraphErrors* z4)
1088{
1089 /// This method does a free fit of a_0, a_2, and a_4 to get the best possible result.
1090 /// The resulting parameters do not necessarily correspond to a meaningful physical result.
1091
1092 assert(data->GetN() == z0->GetN());
1093 assert(data->GetN() == z2->GetN());
1094 assert(data->GetN() == z4->GetN());
1095 // create a copy of the data with cos(theta) as x-axis and fit it with a legenre polynomial
1096 // this is to get initial conditions for our fit
1097 auto* cosTheta = new TGraphErrors(*data);
1098 auto* x = cosTheta->GetX();
1099 for(int i = 0; i < cosTheta->GetN(); ++i) {
1100 x[i] = TMath::Cos(x[i] / 180. * TMath::Pi());
1101 }
1102
1103 int nPar = 3;
1104 auto* legendre = new TF1("legendre", TGRSIFunctions::LegendrePolynomial, -1., 1., nPar);
1105 legendre->SetParNames("a_{0}", "a_{2}", "a_{4}");
1106 legendre->SetParameters(1., 0.5, 0.5);
1107 cosTheta->Fit(legendre, "QN0");
1108
1109 // the actual fitting
1110 Ac ac(data, z0, z2, z4);
1111 ROOT::Fit::Fitter fitter;
1112 fitter.SetFCN(nPar, ac);
1113 // this is a good guess for the initial scale
1114 fitter.Config().ParSettings(0) = ROOT::Fit::ParameterSettings("a_{0}", data->GetMaximum() / z0->GetMaximum(), 0.0001, -10., 10.);
1115 for(int i = 1; i < nPar; ++i) {
1116 // parameter settings arguments are parameter name, initial value, step size, minimum, and maximum
1117 fitter.Config().ParSettings(i) = ROOT::Fit::ParameterSettings(Form("a_{%d}", 2 * i), legendre->GetParameter(i), 0.0001, -10., 10.);
1118 }
1119 fitter.Config().MinimizerOptions().SetPrintLevel(0);
1120 fitter.Config().SetMinimizer("Minuit2", "Migrad"); // or simplex?
1121 if(!fitter.FitFCN()) {
1122 std::cerr << "Fit failed" << std::endl;
1123 return {};
1124 }
1125
1126 // get the fit result and print chi^2 and parameters
1127 auto fitResult = fitter.Result();
1128 // MinFcnValue() is the minimum chi2, ac.Np gives the number of data points
1129 std::cout << "Best reduced chi^2 from free fit at a2 and a4: " << fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters()) << std::endl;
1130 auto parameters = fitResult.Parameters();
1131 auto errors = fitResult.Errors();
1132 std::cout << "Parameters a_0: " << parameters[0] << " +- " << errors[0] << ", a_2: " << parameters[1] << " +- " << errors[1] << ", a_4: " << parameters[2] << " +- " << errors[2] << std::endl;
1133 TMatrixD covariance(nPar, nPar);
1134 fitResult.GetCovarianceMatrix(covariance);
1135
1136 // create fit and residual graphs
1137 auto* fit = static_cast<TGraphErrors*>(z0->Clone("fit"));
1138 for(int i = 0; i < fit->GetN(); ++i) {
1139#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
1140 fit->SetPointY(i, parameters[0] * ((1. - parameters[1] - parameters[2]) * z0->GetPointY(i) + parameters[1] * z2->GetPointY(i) + parameters[2] * z4->GetPointY(i)));
1141#else
1142 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]));
1143#endif
1144 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)));
1145 }
1146
1147 auto* residual = static_cast<TGraphErrors*>(data->Clone("residual"));
1148 for(int i = 0; i < residual->GetN(); ++i) {
1149#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 20, 0)
1150 residual->SetPointY(i, data->GetPointY(i) - fit->GetPointY(i));
1151#else
1152 residual->SetPoint(i, residual->GetX()[i], data->GetY()[i] - fit->GetY()[i]);
1153#endif
1154 residual->SetPointError(i, 0., TMath::Sqrt(TMath::Power(data->GetErrorY(i), 2) + TMath::Power(fit->GetErrorY(i), 2)));
1155 }
1156
1157 double redChiSquare = fitResult.MinFcnValue() / (ac.Np() - fitResult.NFreeParameters());
1158
1159 auto* canvas = new TCanvas;
1160 auto* multiGraph = PlotCanvas(data, fit, residual, parameters, errors, redChiSquare);
1161
1162 fit->Write("A2a4Fit");
1163 residual->Write("Residual");
1164 multiGraph->Write("FitComparison");
1165 canvas->Write("A2a4Canvas");
1166
1167 return parameters;
1168}
double bgHigh
TMultiGraph * PlotCanvas(TGraphErrors *data, TGraphErrors *fit, TGraphErrors *residual, const std::vector< double > &parameters, const std::vector< double > &errors, const double &redChiSquare, const char *extraText=nullptr)
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::ofstream &logFile)
std::vector< double > A2a4Method(TGraphErrors *data, TGraphErrors *z0, TGraphErrors *z2, TGraphErrors *z4)
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:72
void AddPeak(TSinglePeak *peak)
Definition TPeakFitter.h:37
TFitResultPtr Fit(TH1 *fit_hist, Option_t *opt="")
void Centroid(const Double_t &centroid) override
Definition TRWPeak.cxx:5
Double_t CentroidErr() const override
Definition TRWPeak.cxx:67
TF1 * GetFitFunction() const
Definition TSinglePeak.h:72
virtual Double_t FWHMErr()
Double_t GetReducedChi2() const
Definition TSinglePeak.h:87
Double_t AreaErr() const
Definition TSinglePeak.h:52
virtual Double_t FWHM()
Double_t Area() const
Definition TSinglePeak.h:51
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)