MaCh3  2.4.2
Reference Guide
OscProcessor.cpp
Go to the documentation of this file.
1 #include "OscProcessor.h"
2 
4 // ROOT includes
5 #include "TArrow.h"
6 #include "TEllipse.h"
7 #include "TLatex.h"
8 #include "TVector2.h"
10 
11 #pragma GCC diagnostic ignored "-Wfloat-conversion"
12 
13 // ****************************
14 OscProcessor::OscProcessor(const std::string &InputFile) : MCMCProcessor(InputFile) {
15 // ****************************
16  //KS: WARNING this only work when you project from Chain, will nor work when you try SetBranchAddress etc. Turn it on only if you know how to use it
17  PlotJarlskog = false;
18 
25 }
26 
27 // ****************************
28 // The destructor
30 // ****************************
31 }
32 
33 // ***************
34 // Read the Osc cov file and get the input central values and errors
36 // ***************
37  // KS: Check if OscParams were enabled, in future we will also get
38  for(size_t i = 0; i < ParameterGroup.size(); i++) {
39  if(ParameterGroup[i] == "Osc"){
40  OscEnabled = true;
41  break;
42  }
43  }
44 
45  if(OscEnabled)
46  {
47  for (int i = 0; i < nDraw; ++i)
48  {
49  //Those keep which parameter type we run currently and relative number
50  const int ParamEnum = ParamType[i];
51  const int ParamNo = i - ParamTypeStartPos[ParameterEnum(ParamEnum)];
52  const std::string CurrentName = ParamNames[ParamEnum][ParamNo].Data();
53 
55  if (CurrentName == "sin2th_13") {
56  Sin2Theta13Index = i;
57  Sin2Theta13Name = CurrentName;
58  } else if (CurrentName == "sin2th_12") {
59  Sin2Theta12Index = i;
60  Sin2Theta12Name = CurrentName;
61  } else if (CurrentName == "sin2th_23") {
62  Sin2Theta23Index = i;
63  Sin2Theta23Name = CurrentName;
64  } else if (CurrentName == "delta_cp") {
65  DeltaCPIndex = i;
66  DeltaCPName = CurrentName;
67  } else if (CurrentName == "delm2_23") {
68  DeltaM2_23Index = i;
69  DeltaM2_23Name = CurrentName;
70  }
71  }
72  } else{
73  MACH3LOG_WARN("Didn't find oscillation parameters");
74  }
75 
77  {
78  Chain->SetAlias("J_cp", "TMath::Sqrt(sin2th_13)*TMath::Sqrt(1.-sin2th_13)*TMath::Sqrt(1.-sin2th_13)*TMath::Sqrt(sin2th_12)*TMath::Sqrt(1.-sin2th_12)*TMath::Sqrt(sin2th_23)*TMath::Sqrt(1.-sin2th_23)*TMath::Sin(delta_cp)");
79  BranchNames.push_back("J_cp");
80  ParamType.push_back(kXSecPar);
81  nParam[kXSecPar]++;
82  nDraw++;
83 
85  ParamCentral[kXSecPar].push_back( 0. );
86  ParamErrors[kXSecPar].push_back( 1. );
87  // Push back the name
88  ParamNames[kXSecPar].push_back("J_cp");
89  ParamFlat[kXSecPar].push_back( false );
90  } else if(PlotJarlskog && !OscEnabled) {
91  MACH3LOG_ERROR("Trying to enable Jarlskog without oscillations");
92  throw MaCh3Exception(__FILE__,__LINE__);
93  }
94 }
95 
96 // ***************
97 // Calculate Jarlskog Invariant using oscillation parameters
98 double OscProcessor::CalcJarlskog(const double s2th13, const double s2th23, const double s2th12, const double dcp) const {
99 // ***************
100  const double s13 = std::sqrt(s2th13);
101  const double s23 = std::sqrt(s2th23);
102  const double s12 = std::sqrt(s2th12);
103  const double sdcp = std::sin(dcp);
104  const double c13 = std::sqrt(1.-s2th13);
105  const double c12 = std::sqrt(1.-s2th12);
106  const double c23 = std::sqrt(1.-s2th23);
107 
108  const double j = s13*c13*c13*s12*c12*s23*c23*sdcp;
109 
110  return j;
111 }
112 
113 // ***************
114 double OscProcessor::SamplePriorForParam(const int paramIndex, const std::unique_ptr<TRandom3>& randGen, const std::vector<double>& FlatBounds) const {
115 // ***************
116  TString Title = "";
117  double Prior = 1.0, PriorError = 1.0;
118  bool FlatPrior = false;
119 
120  // Get info for this parameter
121  GetNthParameter(paramIndex, Prior, PriorError, Title);
122 
123  ParameterEnum ParType = ParamType[paramIndex];
124  int ParamTemp = paramIndex - ParamTypeStartPos[ParType];
125  FlatPrior = ParamFlat[ParType][ParamTemp];
126 
127  if (FlatPrior) {
128  return randGen->Uniform(FlatBounds[0], FlatBounds[1]);
129  } else {
130  // Gaussian prior centered at Prior with width PriorError
131  return randGen->Gaus(Prior, PriorError);
132  }
133 }
134 
135 // ***************
136 // Perform Several Jarlskog Plotting
138 // ***************
139  if(!OscEnabled ||
145  {
146  MACH3LOG_WARN("Will not {}, as oscillation parameters are missing", __func__);
147  return;
148  }
149  MACH3LOG_INFO("Starting {}", __func__);
150 
151  bool DoReweight = false;
152 
153  double s2th13, s2th23, s2th12, dcp, dm2 = M3::_BAD_DOUBLE_;
154  double weight = 1.0;
155  std::pair<double, double> Sin13_NewPrior;
156 
157  // Now read the MCMC file
158  TFile *TempFile = M3::Open((MCMCFile + ".root"), "open", __FILE__, __LINE__);
159 
160  // Get the settings for the MCMC
161  TMacro *Config = TempFile->Get<TMacro>("Reweight_Config");
162 
163  if (Config != nullptr) {
164  MACH3LOG_INFO("Found Reweight_Config in chain");
165 
166  // Print the reweight configuration for user info
167  YAML::Node Settings = TMacroToYAML(*Config);
168  // Simple check: only enable DoReweight if it's a 1D sin2th_13 Gaussian reweight since Savage Dickey process later on generates values from the Gaussian
169  if(CheckNodeExists(Settings, "ReweightMCMC")) {
170  YAML::Node firstReweight = Settings["ReweightMCMC"].begin()->second;
171  int dimension = GetFromManager<int>(firstReweight["ReweightDim"], 1);
172  std::string reweightType = GetFromManager<std::string>(firstReweight["ReweightType"], "");
173  auto paramNames = GetFromManager<std::vector<std::string>>(firstReweight["ReweightVar"], {});
174  if (dimension == 1 && reweightType == "Gaussian" && paramNames.size() == 1){
175  Sin13_NewPrior = Get<std::pair<double, double>>(firstReweight["ReweightPrior"],__FILE__,__LINE__);
176  DoReweight = true;
177  } else {
178  MACH3LOG_INFO("No valid reweighting configuration (1D Gaussian on sin2th_13 only) found for Jarlskog analysis");
179  }
180  } else {
181  MACH3LOG_INFO("No reweighting configuration found for Jarlskog analysis");
182  }
183 }
184 
185  TempFile->Close();
186  delete TempFile;
187 
188  TDirectory *JarlskogDir = OutputFile->mkdir("Jarlskog");
189  JarlskogDir->cd();
190 
191  unsigned int step = 0;
192  Chain->SetBranchStatus("*", false);
193 
194  Chain->SetBranchStatus(Sin2Theta13Name.c_str(), true);
195  Chain->SetBranchAddress(Sin2Theta13Name.c_str(), &s2th13);
196 
197  Chain->SetBranchStatus(Sin2Theta23Name.c_str(), true);
198  Chain->SetBranchAddress(Sin2Theta23Name.c_str(), &s2th23);
199 
200  Chain->SetBranchStatus(Sin2Theta12Name.c_str(), true);
201  Chain->SetBranchAddress(Sin2Theta12Name.c_str(), &s2th12);
202 
203  Chain->SetBranchStatus(DeltaCPName.c_str(), true);
204  Chain->SetBranchAddress(DeltaCPName.c_str(), &dcp);
205 
206  Chain->SetBranchStatus(DeltaM2_23Name.c_str(), true);
207  Chain->SetBranchAddress(DeltaM2_23Name.c_str(), &dm2);
208 
209  Chain->SetBranchStatus("step", true);
210  Chain->SetBranchAddress("step", &step);
211 
212  if(DoReweight) {
213  Chain->SetBranchStatus("Weight", true);
214  Chain->SetBranchAddress("Weight", &weight);
215  } else {
216  MACH3LOG_WARN("Not applying reweighting weight");
217  weight = 1.0;
218  }
219 
220  // Original histograms
221  auto jarl = std::make_unique<TH1D>("jarl", "jarl", 1000, -0.05, 0.05);
222  jarl->SetDirectory(nullptr);
223  auto jarl_th23 = std::make_unique<TH2D>("jarl_th23", "jarl_th23", 500, -0.05, 0.05, 500, 0.3, 0.7);
224  jarl_th23->SetDirectory(nullptr);
225  auto jarl_dcp = std::make_unique<TH2D>("jarl_dcp", "jarl_dcp", 500, -0.05, 0.05, 500, -1. * TMath::Pi(), TMath::Pi());
226  jarl_dcp->SetDirectory(nullptr);
227 
228  jarl->SetTitle("Jarlskog Invariant;J #equiv s_{13}c_{13}^{2}s_{12}c_{12}s_{23}c_{23}sin#delta;Posterior probability");
229  jarl_th23->SetTitle("Jarlskog Invariant;J #equiv s_{13}c_{13}^{2}s_{12}c_{12}s_{23}c_{23}sin#delta;Posterior probability");
230 
231  // Clones
232  auto jarl_IH = M3::Clone(jarl.get(), "jarl_IH");
233  auto jarl_NH = M3::Clone(jarl.get(), "jarl_NH");
234 
235  auto jarl_th23_IH = M3::Clone(jarl_th23.get(), "jarl_th23_IH");
236  auto jarl_th23_NH = M3::Clone(jarl_th23.get(), "jarl_th23_NH");
237 
238  auto jarl_dcp_IH = M3::Clone(jarl_dcp.get(), "jarl_dcp_IH");
239  auto jarl_dcp_NH = M3::Clone(jarl_dcp.get(), "jarl_dcp_NH");
240 
241  auto jarl_flatsindcp = M3::Clone(jarl.get(), "jarl_flatsindcp");
242  auto jarl_IH_flatsindcp = M3::Clone(jarl.get(), "jarl_IH_flatsindcp");
243  auto jarl_NH_flatsindcp = M3::Clone(jarl.get(), "jarl_NH_flatsindcp");
244 
245  auto jarl_th23_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_flatsindcp");
246  auto jarl_th23_IH_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_IH_flatsindcp");
247  auto jarl_th23_NH_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_NH_flatsindcp");
248 
249  auto jarl_prior = M3::Clone(jarl.get(), "jarl_prior");
250  auto jarl_prior_flatsindcp = M3::Clone(jarl.get(), "jarl_prior_flatsindcp");
251  std::unique_ptr<TH1D> jarl_wRC_prior, jarl_wRC_prior_flatsindcp, jarl_wRC_prior_t2kth23;
252  // Only use this if chain has reweigh weight [mostly coming from Reactor Constrains]
253  if(DoReweight){
254  jarl_wRC_prior = M3::Clone(jarl.get(), "jarl_wRC_prior");
255  jarl_wRC_prior_flatsindcp = M3::Clone(jarl.get(), "jarl_wRC_prior_flatsindcp");
256  jarl_wRC_prior_t2kth23 = M3::Clone(jarl.get(), "jarl_wRC_prior_flatsindcp");
257  }
258 
259  // to apply a prior that is flat in sin(dcp) intead of dcp
260  auto prior3 = std::make_unique<TF1>("prior3", "TMath::Abs(TMath::Cos(x))");
261 
262  // T2K prior is flat (and uncorrelated) in dcp, sin^2(th13), sin^2(th23)
263  auto randGen = std::make_unique<TRandom3>(0);
264  const Long64_t countwidth = nEntries/5;
265 
266  for(int i = 0; i < nEntries; ++i) {
267  if (i % countwidth == 0) {
270  } else {
271  Chain->GetEntry(i);
272  }
273 
274  if(step < BurnInCut) continue; // burn-in cut
275 
276  const double j = CalcJarlskog(s2th13, s2th23, s2th12, dcp);
277  const double prior_weight = prior3->Eval(dcp);
278 
279  jarl->Fill(j, weight);
280  jarl_th23->Fill(j, s2th23, weight);
281  jarl_dcp->Fill(j, dcp, weight);
282 
283  jarl_flatsindcp->Fill(j, prior_weight*weight);
284  jarl_th23_flatsindcp->Fill(j, s2th23, prior_weight*weight);
285 
286  const double prior_s2th13 = SamplePriorForParam(Sin2Theta13Index, randGen, {0.,1.});
287  const double prior_s2th23 = SamplePriorForParam(Sin2Theta23Index, randGen, {0.,1.});
288  const double prior_s2th12 = SamplePriorForParam(Sin2Theta12Index, randGen, {0.,1.});
289  const double prior_dcp = SamplePriorForParam(DeltaCPIndex, randGen, {-1.*TMath::Pi(),TMath::Pi()});
290  // KS: This is hardcoded but we always assume flat in delta CP so probably fine
291  const double prior_sindcp = randGen->Uniform(-1., 1.);
292 
293  const double prior_s13 = std::sqrt(prior_s2th13);
294  const double prior_s23 = std::sqrt(prior_s2th23);
295  const double prior_s12 = std::sqrt(prior_s2th12);
296  const double prior_sdcp = std::sin(prior_dcp);
297  const double prior_c13 = std::sqrt(1.-prior_s2th13);
298  const double prior_c12 = std::sqrt(1.-prior_s2th12);
299  const double prior_c23 = std::sqrt(1.-prior_s2th23);
300  const double prior_j = prior_s13*prior_c13*prior_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sdcp;
301  const double prior_flatsindcp_j = prior_s13*prior_c13*prior_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sindcp;
302 
303  jarl_prior->Fill(prior_j);
304  jarl_prior_flatsindcp->Fill(prior_flatsindcp_j);
305 
306  if(DoReweight) {
307  const double prior_wRC_s2th13 = randGen->Gaus(Sin13_NewPrior.first, Sin13_NewPrior.second);
308  const double prior_wRC_s13 = std::sqrt(prior_wRC_s2th13);
309  const double prior_wRC_c13 = std::sqrt(1.-prior_wRC_s2th13);
310  const double prior_wRC_j = prior_wRC_s13*prior_wRC_c13*prior_wRC_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sdcp;
311  const double prior_wRC_flatsindcp_j = prior_wRC_s13*prior_wRC_c13*prior_wRC_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sindcp;
312  const double s23 = std::sqrt(s2th23);
313  const double c23 = std::sqrt(1.-s2th23);
314 
315  jarl_wRC_prior->Fill(prior_wRC_j);
316  jarl_wRC_prior_flatsindcp->Fill(prior_wRC_flatsindcp_j);
317  jarl_wRC_prior_t2kth23->Fill(prior_wRC_s13*prior_wRC_c13*prior_wRC_c13*prior_s12*prior_c12*s23*c23*prior_sdcp);
318  }
319 
320  if(dm2 > 0.) {
321  jarl_NH->Fill(j, weight);
322  jarl_th23_NH->Fill(j, s2th23, weight);
323  jarl_dcp_NH->Fill(j, dcp, weight);
324  jarl_NH_flatsindcp->Fill(j, prior_weight*weight);
325  jarl_th23_NH_flatsindcp->Fill(j, s2th23, prior_weight*weight);
326  }
327  else if(dm2 < 0.) {
328  jarl_IH->Fill(j, weight);
329  jarl_th23_IH->Fill(j, s2th23, weight);
330  jarl_dcp_IH->Fill(j, dcp, weight);
331  jarl_IH_flatsindcp->Fill(j, prior_weight*weight);
332  jarl_th23_IH_flatsindcp->Fill(j, s2th23, prior_weight*weight);
333  }
334  }
335 
336  jarl->Write("jarlskog_both");
337  jarl_NH->Write("jarlskog_NH");
338  jarl_IH->Write("jarlskog_IH");
339  jarl_th23->Write("jarlskog_th23_both");
340  jarl_th23_NH->Write("jarlskog_th23_NH");
341  jarl_th23_IH->Write("jarlskog_th23_IH");
342 
343  jarl_dcp->Write("jarlskog_dcp_both");
344  jarl_dcp_NH->Write("jarlskog_dcp_NH");
345  jarl_dcp_IH->Write("jarlskog_dcp_IH");
346 
347 
348  jarl_flatsindcp->Write("jarlskog_both_flatsindcp");
349  jarl_NH_flatsindcp->Write("jarlskog_NH_flatsindcp");
350  jarl_IH_flatsindcp->Write("jarlskog_IH_flatsindcp");
351  jarl_th23_flatsindcp->Write("jarlskog_th23_both_flatsindcp");
352  jarl_th23_NH_flatsindcp->Write("jarlskog_th23_NH_flatsindcp");
353  jarl_th23_IH_flatsindcp->Write("jarlskog_th23_IH_flatsindcp");
354 
355  jarl_prior->Write("jarl_prior");
356  jarl_prior_flatsindcp->Write("jarl_prior_flatsindcp");
357  if(DoReweight) {
358  jarl_wRC_prior->Write("jarl_wRC_prior");
359  jarl_wRC_prior_flatsindcp->Write("jarl_wRC_prior_flatsindcp");
360  jarl_wRC_prior_t2kth23->Write("jarl_wRC_prior_t2kth23");
361  }
362 
363  MakeJarlskogPlot(jarl, jarl_flatsindcp,
364  jarl_NH, jarl_NH_flatsindcp,
365  jarl_IH, jarl_IH_flatsindcp);
366 
367  // Perform Savage Dickey analysis
368  if(DoReweight) {
369  SavageDickeyPlot(jarl, jarl_wRC_prior, "Jarlskog flat #delta_{CP}", 0);
370  SavageDickeyPlot(jarl_flatsindcp, jarl_wRC_prior_flatsindcp, "Jarlskog flat sin#delta_{CP}", 0);
371  } else {
372  SavageDickeyPlot(jarl, jarl_prior, "Jarlskog flat #delta_{CP}", 0);
373  SavageDickeyPlot(jarl_flatsindcp, jarl_prior_flatsindcp, "Jarlskog flat sin#delta_{CP}", 0);
374  }
375 
376  JarlskogDir->Close();
377  delete JarlskogDir;
378 
379  Chain->SetBranchStatus("*", true);
380  OutputFile->cd();
381 }
382 
383 
384 // ***************
385 void OscProcessor::MakeJarlskogPlot(const std::unique_ptr<TH1D>& jarl,
386  const std::unique_ptr<TH1D>& jarl_flatsindcp,
387  const std::unique_ptr<TH1D>& jarl_NH,
388  const std::unique_ptr<TH1D>& jarl_NH_flatsindcp,
389  const std::unique_ptr<TH1D>& jarl_IH,
390  const std::unique_ptr<TH1D>& jarl_IH_flatsindcp) {
391 // ***************
392  MACH3LOG_INFO("Starting {}", __func__);
393  int originalErrorLevel = gErrorIgnoreLevel;
394  gErrorIgnoreLevel = kFatal;
395 
396  // 1-->NH, 0-->both, -1-->IH
397  for(int hierarchy = -1; hierarchy <= 1; hierarchy++)
398  {
399  std::unique_ptr<TH1D> j_hist;
400  std::unique_ptr<TH1D> j_hist_sdcp;
401  if(hierarchy == 1) {
402  j_hist = M3::Clone(jarl_NH.get(), "");
403  j_hist_sdcp = M3::Clone(jarl_NH_flatsindcp.get(), "");
404  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
405  } else if(hierarchy == 0) {
406  j_hist = M3::Clone(jarl.get(), "");
407  j_hist_sdcp = M3::Clone(jarl_flatsindcp.get(), "");
408  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
409  } else if(hierarchy == -1) {
410  j_hist = M3::Clone(jarl_IH.get(), "");
411  j_hist_sdcp = M3::Clone(jarl_IH_flatsindcp.get(), "");
412  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
413  } else {
414  MACH3LOG_ERROR("Invalid hierarchy option. 1 for NH, 0 for both, -1 for IH");
415  throw MaCh3Exception(__FILE__ , __LINE__ );
416  }
417 
418  j_hist->Rebin(7);
419  j_hist_sdcp->Rebin(7);
420 
421  j_hist->SetLineColor(kAzure-2);
422  j_hist_sdcp->SetLineColor(kOrange+1);
423  j_hist->SetLineWidth(2);
424  j_hist_sdcp->SetLineWidth(2);
425 
426  auto StyleAxis = [](TH1* h) {
427  auto xAxis = h->GetXaxis();
428  auto yAxis = h->GetYaxis();
429 
430  xAxis->SetLabelSize(0.04);
431  xAxis->SetLabelFont(132);
432  xAxis->SetTitleSize(0.04);
433  xAxis->SetTitleOffset(0.80);
434  xAxis->SetTitleFont(132);
435  xAxis->SetNdivisions(505);
436  xAxis->SetTickSize(0.04);
437 
438  yAxis->SetLabelSize(0.04);
439  yAxis->SetLabelFont(132);
440  yAxis->SetTitleSize(0.04);
441  yAxis->SetTitleOffset(1.2);
442  yAxis->SetTitleFont(132);
443  yAxis->SetNdivisions(505);
444  yAxis->SetTickSize(0.04);
445  };
446 
447  StyleAxis(j_hist.get());
448 
449  j_hist->GetXaxis()->SetRangeUser(-0.04,0.04);
450  j_hist->Scale(1./j_hist->Integral());
451  j_hist_sdcp->Scale(1./j_hist_sdcp->Integral());
452 
453  std::unique_ptr<TH1D> j_hist_copy = M3::Clone(j_hist.get(), "j_hist_copy");
454  std::unique_ptr<TH1D> j_hist_1sig = M3::Clone(j_hist.get(), "j_hist_1sig");
455  std::unique_ptr<TH1D> j_hist_2sig = M3::Clone(j_hist.get(), "j_hist_2sig");
456  std::unique_ptr<TH1D> j_hist_3sig = M3::Clone(j_hist.get(), "j_hist_3sig");
457 
458  //upper and lower edges
459  double j_bf = j_hist_copy->GetXaxis()->GetBinCenter(j_hist_copy->GetMaximumBin());
460  double j_1sig_low = 9999999.;
461  double j_1sig_up = -9999999.;
462  double j_2sig_low = 9999999.;;
463  double j_2sig_up = -9999999.;
464  double j_3sig_low = 9999999.;;
465  double j_3sig_up = -9999999.;
466 
467 
468  std::unique_ptr<TH1D> j_hist_sdcp_copy = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_copy");
469  std::unique_ptr<TH1D> j_hist_sdcp_1sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_1sig");
470  std::unique_ptr<TH1D> j_hist_sdcp_2sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_2sig");
471  std::unique_ptr<TH1D> j_hist_sdcp_3sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_3sig");
472 
473  //upper and lower edges
474  double j_sdcp_1sig_low = 9999999.;
475  double j_sdcp_1sig_up = -9999999.;
476  double j_sdcp_2sig_low = 9999999.;;
477  double j_sdcp_2sig_up = -9999999.;
478  double j_sdcp_3sig_low = 9999999.;;
479  double j_sdcp_3sig_up = -9999999.;
480 
481  double contlevel1 = 0.68;
482  double contlevel2 = 0.90;
483  double contlevel4 = 0.99;
484  double contlevel5 = 0.9973;
485  double integral, tsum = 0.;
486 
487  integral = j_hist_copy->Integral();
488 
489  while((tsum/integral)<contlevel5) {
490  double tmax = j_hist_copy->GetMaximum();
491  int bin = j_hist_copy->GetMaximumBin();
492  double xval = j_hist_copy->GetXaxis()->GetBinCenter(bin);
493  double xwidth = j_hist_copy->GetXaxis()->GetBinWidth(bin);
494  if((tsum/integral)<contlevel1) {
495  j_hist_copy->SetBinContent(bin,-1.0);
496  j_hist_1sig->SetBinContent(bin,0.);
497  j_hist_2sig->SetBinContent(bin,0.);
498  j_hist_3sig->SetBinContent(bin,0.);
499  if(xval<j_1sig_low && xval<j_bf) j_1sig_low = xval - xwidth/2.;
500  if(xval>j_1sig_up && xval>j_bf) j_1sig_up = xval + xwidth/2.;
501  }
502  if((tsum/integral)<contlevel2 && (tsum / integral > contlevel1) ) {
503  j_hist_copy->SetBinContent(bin,-5.0);
504  j_hist_2sig->SetBinContent(bin,0.);
505  j_hist_3sig->SetBinContent(bin,0.);
506  if(xval<j_2sig_low && xval<j_bf) j_2sig_low = xval - xwidth/2.;
507  if(xval>j_2sig_up && xval>j_bf) j_2sig_up = xval + xwidth/2.;
508  }
509  if((tsum/integral)<contlevel4 && (tsum / integral > contlevel1) ) {
510  j_hist_copy->SetBinContent(bin,-9.0);
511  j_hist_3sig->SetBinContent(bin,0.);
512  if(xval < j_3sig_low && xval <j_bf) j_3sig_low = xval - xwidth/2.;
513  if(xval > j_3sig_up && xval > j_bf) j_3sig_up = xval + xwidth/2.;
514  }
515  tsum+=tmax;
516  }
517 
518  integral = j_hist_sdcp_copy->Integral();
519  tsum = 0.;
520 
521  while((tsum/integral)<contlevel5) {
522  double tmax = j_hist_sdcp_copy->GetMaximum();
523  int bin = j_hist_sdcp_copy->GetMaximumBin();
524  double xval = j_hist_sdcp_copy->GetXaxis()->GetBinCenter(bin);
525  double xwidth = j_hist_sdcp_copy->GetXaxis()->GetBinWidth(bin);
526  if((tsum/integral)<contlevel1) {
527  j_hist_sdcp_copy->SetBinContent(bin,-1.0);
528  j_hist_sdcp_1sig->SetBinContent(bin,0.);
529  j_hist_sdcp_2sig->SetBinContent(bin,0.);
530  j_hist_sdcp_3sig->SetBinContent(bin,0.);
531  if(xval<j_sdcp_1sig_low && xval<j_bf) j_sdcp_1sig_low = xval - xwidth/2.;
532  if(xval>j_sdcp_1sig_up && xval>j_bf) j_sdcp_1sig_up = xval + xwidth/2.;
533  }
534  if((tsum/integral)<contlevel2 && (tsum / integral > contlevel1) ) {
535  j_hist_sdcp_copy->SetBinContent(bin,-5.0);
536  j_hist_sdcp_2sig->SetBinContent(bin,0.);
537  j_hist_sdcp_3sig->SetBinContent(bin,0.);
538  if(xval<j_sdcp_2sig_low && xval<j_bf) j_sdcp_2sig_low = xval - xwidth/2.;
539  if(xval>j_sdcp_2sig_up && xval>j_bf) j_sdcp_2sig_up = xval + xwidth/2.;
540  }
541  if((tsum/integral)<contlevel4 && (tsum / integral > contlevel1) ) {
542  j_hist_sdcp_copy->SetBinContent(bin,-9.0);
543  j_hist_sdcp_3sig->SetBinContent(bin,0.);
544  if(xval<j_sdcp_3sig_low && xval<j_bf) j_sdcp_3sig_low = xval - xwidth/2.;
545  if(xval>j_sdcp_3sig_up && xval>j_bf) j_sdcp_3sig_up = xval + xwidth/2.;
546  }
547  tsum+=tmax;
548  }
549 
550  j_hist_1sig->SetLineStyle(9);
551  j_hist_sdcp_1sig->SetLineStyle(9);
552  j_hist_2sig->SetLineStyle(7);
553  j_hist_sdcp_2sig->SetLineStyle(7);
554  j_hist_3sig->SetLineStyle(2);
555  j_hist_sdcp_3sig->SetLineStyle(2);
556 
557  auto ldash = std::make_unique<TH1D>("ldash", "solid", 10, -0.04, 0.04);
558  auto sdash = std::make_unique<TH1D>("sdash", "dashed", 10, -0.04, 0.04);
559  auto fdash = std::make_unique<TH1D>("fdash", "fdashed",10, -0.04, 0.04);
560  ldash->SetLineColor(kBlack);
561  sdash->SetLineColor(kBlack);
562  fdash->SetLineColor(kBlack);
563  ldash->SetLineWidth(2);
564  sdash->SetLineWidth(2);
565  fdash->SetLineWidth(2);
566  ldash->SetLineStyle(9);
567  sdash->SetLineStyle(7);
568  fdash->SetLineStyle(2);
569 
570  double vertUp = 0.5 * j_hist->GetMaximum();
571  auto jline_1sig_low = std::make_unique<TLine>(j_1sig_low, 0., j_1sig_low, vertUp);
572  auto jline_2sig_low = std::make_unique<TLine>(j_2sig_low, 0., j_2sig_low, vertUp);
573  auto jline_3sig_low = std::make_unique<TLine>(j_3sig_low, 0., j_3sig_low, vertUp);
574 
575  auto jline_1sig_up = std::make_unique<TLine>(j_1sig_up, 0., j_1sig_up,vertUp);
576  auto jline_2sig_up = std::make_unique<TLine>(j_2sig_up, 0., j_2sig_up,vertUp);
577  auto jline_3sig_up = std::make_unique<TLine>(j_3sig_up, 0., j_3sig_up,vertUp);
578 
579  auto jline_sdcp_1sig_low = std::make_unique<TLine>(j_sdcp_1sig_low, 0., j_sdcp_1sig_low, vertUp);
580  auto jline_sdcp_2sig_low = std::make_unique<TLine>(j_sdcp_2sig_low, 0., j_sdcp_2sig_low, vertUp);
581  auto jline_sdcp_3sig_low = std::make_unique<TLine>(j_sdcp_3sig_low, 0., j_sdcp_3sig_low, vertUp);
582 
583  auto jline_sdcp_1sig_up = std::make_unique<TLine>(j_sdcp_1sig_up, 0., j_sdcp_1sig_up, vertUp);
584  auto jline_sdcp_2sig_up = std::make_unique<TLine>(j_sdcp_2sig_up, 0., j_sdcp_2sig_up, vertUp);
585  auto jline_sdcp_3sig_up = std::make_unique<TLine>(j_sdcp_3sig_up, 0., j_sdcp_3sig_up, vertUp);
586 
587  double arrowLength = 0.003;
588  double arrowHeight = vertUp;
589 
590  auto MakeArrow = [&](double x, Color_t color, Width_t width) -> std::unique_ptr<TArrow> {
591  auto arrow = std::make_unique<TArrow>(x, arrowHeight, x - arrowLength, arrowHeight, 0.02, ">");
592  arrow->SetLineColor(color);
593  arrow->SetLineWidth(width);
594  return arrow;
595  };
596 
597  auto j_arrow_1sig_up = MakeArrow(j_1sig_up, j_hist_1sig->GetLineColor(), j_hist_1sig->GetLineWidth());
598  auto j_arrow_2sig_up = MakeArrow(j_2sig_up, j_hist_2sig->GetLineColor(), j_hist_2sig->GetLineWidth());
599  auto j_arrow_3sig_up = MakeArrow(j_3sig_up, j_hist_3sig->GetLineColor(), j_hist_3sig->GetLineWidth());
600 
601  auto j_sdcp_arrow_1sig_up = MakeArrow(j_sdcp_1sig_up, j_hist_sdcp_1sig->GetLineColor(), j_hist_sdcp_1sig->GetLineWidth());
602  auto j_sdcp_arrow_2sig_up = MakeArrow(j_sdcp_2sig_up, j_hist_sdcp_2sig->GetLineColor(), j_hist_sdcp_2sig->GetLineWidth());
603  auto j_sdcp_arrow_3sig_up = MakeArrow(j_sdcp_3sig_up, j_hist_sdcp_3sig->GetLineColor(), j_hist_sdcp_3sig->GetLineWidth());
604 
605  MACH3LOG_DEBUG("j_1sig_low = {:.4f}, j_2sig_low = {:.4f}, j_3sig_low = {:.4f}", j_1sig_low, j_2sig_low, j_3sig_low);
606  MACH3LOG_DEBUG("j_1sig_up = {:.4f}, j_2sig_up = {:.4f}, j_3sig_up = {:.4f}", j_1sig_up, j_2sig_up, j_3sig_up);
607 
608  auto CopyLineStyle = [](const TH1D* src, TLine* dst) {
609  dst->SetLineColor(src->GetLineColor());
610  dst->SetLineStyle(src->GetLineStyle());
611  dst->SetLineWidth(src->GetLineWidth());
612  };
613 
614  CopyLineStyle(j_hist_1sig.get(), jline_1sig_low.get());
615  CopyLineStyle(j_hist_1sig.get(), jline_1sig_up.get());
616  CopyLineStyle(j_hist_2sig.get(), jline_2sig_low.get());
617  CopyLineStyle(j_hist_2sig.get(), jline_2sig_up.get());
618  CopyLineStyle(j_hist_3sig.get(), jline_3sig_low.get());
619  CopyLineStyle(j_hist_3sig.get(), jline_3sig_up.get());
620 
621  CopyLineStyle(j_hist_sdcp_1sig.get(), jline_sdcp_1sig_low.get());
622  CopyLineStyle(j_hist_sdcp_1sig.get(), jline_sdcp_1sig_up.get());
623  CopyLineStyle(j_hist_sdcp_2sig.get(), jline_sdcp_2sig_low.get());
624  CopyLineStyle(j_hist_sdcp_2sig.get(), jline_sdcp_2sig_up.get());
625  CopyLineStyle(j_hist_sdcp_3sig.get(), jline_sdcp_3sig_low.get());
626  CopyLineStyle(j_hist_sdcp_3sig.get(), jline_sdcp_3sig_up.get());
627 
628  auto leg = std::make_unique<TLegend>(0.45, 0.60, 0.75, 0.90);
629  leg->SetTextSize(0.05);
630  leg->SetFillStyle(0);
631  leg->SetNColumns(1);
632  leg->SetTextFont(132);
633  leg->SetBorderSize(0);
634 
635  leg->AddEntry(j_hist.get(), "Prior flat in #delta_{CP}", "l");
636  leg->AddEntry(j_hist_sdcp.get(), "Prior flat in sin#delta_{CP}", "l");
637  leg->AddEntry(ldash.get(), "68% CI", "l");
638  leg->AddEntry(sdash.get(), "90% CI", "l");
639  leg->AddEntry(fdash.get(), "99% CI", "l");
640 
641  j_hist->GetYaxis()->SetRangeUser(0., j_hist->GetMaximum()*1.15);
642  j_hist->Draw("h");
643  j_hist_sdcp->Draw("same h");
644 
645  jline_sdcp_1sig_up->Draw("same");
646  jline_sdcp_2sig_up->Draw("same");
647  jline_sdcp_3sig_up->Draw("same");
648  jline_1sig_up->Draw("same");
649  jline_2sig_up->Draw("same");
650  jline_3sig_up->Draw("same");
651 
652  j_arrow_1sig_up->Draw();
653  j_arrow_2sig_up->Draw();
654  j_arrow_3sig_up->Draw();
655  j_sdcp_arrow_1sig_up->Draw();
656  j_sdcp_arrow_2sig_up->Draw();
657  j_sdcp_arrow_3sig_up->Draw();
658  leg->Draw("same");
659 
660  auto ttext = std::make_unique<TText>();
661  ttext->SetNDC(); // Use normalized device coordinates
662  ttext->SetTextSize(0.03); // Adjust size as needed
663  ttext->SetTextAlign(13); // Align left-top
664 
665  if (hierarchy == 1) ttext->DrawText(0.15, 0.85, "Normal Ordering");
666  else if (hierarchy == 0) ttext->DrawText(0.15, 0.85, "Both Orderings");
667  else if (hierarchy == -1) ttext->DrawText(0.15, 0.85, "Inverted Ordering");
668 
669  gPad->RedrawAxis();
670  Posterior->Update();
671  gPad->Update();
672 
673  Posterior->Print(CanvasName);
674 
675  if(hierarchy == 1) Posterior->Write("jarl1D_NH_comp");
676  else if(hierarchy == 0) Posterior->Write("jarl1D_both_comp");
677  else if(hierarchy == -1) Posterior->Write("jarl1D_IH_comp");
678  }
679 
680  gErrorIgnoreLevel = originalErrorLevel;
681 }
682 
683 // ***************
685 // ***************
687  {
688  MACH3LOG_WARN("Will not {}, as oscillation parameters are missing", __func__);
689  return;
690  }
691  MACH3LOG_INFO("Starting {}", __func__);
692 
693  // get best fit for delta CP
694  const double best_fit = (*Means_HPD)(DeltaCPIndex);
695 
696  const double sigma_p = (*Errors_HPD_Positive)(DeltaCPIndex);
697  const double sigma_n = (*Errors_HPD_Negative)(DeltaCPIndex);
698  // make sure result is between -pi and pi
699  auto wrap_pi = [](double x) {
700  while (x > TMath::Pi()) x -= 2*TMath::Pi();
701  while (x < -TMath::Pi()) x += 2*TMath::Pi();
702  return x;
703  };
704 
705  std::array<double, 6> bds;
706  bds[0] = wrap_pi(best_fit - 3.0 * sigma_n); // -3σ
707  bds[1] = wrap_pi(best_fit - 2.0 * sigma_n); // -2σ
708  bds[2] = wrap_pi(best_fit - 1.0 * sigma_n); // -1σ
709  bds[3] = wrap_pi(best_fit + 1.0 * sigma_p); // +1σ
710  bds[4] = wrap_pi(best_fit + 2.0 * sigma_p); // +2σ
711  bds[5] = wrap_pi(best_fit + 3.0 * sigma_p); // +3σ
712 
713  constexpr double radius = 0.4;
714  constexpr double rad_to_deg = 180.0 / TMath::Pi();
715 
716  // ROOT expects TEllipse angles in degrees, counterclockwise from the x-axis.
717  // If phimax < phimin, ROOT draws counterclockwise across the full circle, causing overlaps.
718  // This ensures threesigA slice stays within the intended range.
719  auto normalize_angle = [](double rad) {
720  // If rad is negative, add 2*pi to wrap into [0, 2*pi)
721  if (rad < 0) rad += 2.0 * TMath::Pi();
722  return rad;
723  };
724 
725  TEllipse onesig (0.5, 0.5, radius, radius, bds[2] * rad_to_deg, bds[4] * rad_to_deg);
726  TEllipse twosigA (0.5, 0.5, radius, radius, bds[1] * rad_to_deg, bds[2] * rad_to_deg);
727  TEllipse twosigB (0.5, 0.5, radius, radius, bds[3] * rad_to_deg, bds[4] * rad_to_deg);
728 
729  // three sigma slices
730  TEllipse threesigA(0.5, 0.5, radius, radius, bds[0] * rad_to_deg, normalize_angle(bds[1]) * rad_to_deg);
731  TEllipse threesigB(0.5, 0.5, radius, radius, bds[4] * rad_to_deg, bds[5] * rad_to_deg);
732 
733  // Remaining slices
734  TEllipse rest(0.5, 0.5, radius, radius, bds[5]*rad_to_deg, bds[0]*rad_to_deg);
735  TEllipse restA(0.5, 0.5, radius, radius, bds[5]*rad_to_deg, 180.0);
736  TEllipse restB(0.5, 0.5, radius, radius, -180.0, bds[0]*rad_to_deg);
737 
738  onesig.SetFillColor(13);
739  twosigA.SetFillColor(12);
740  twosigB.SetFillColor(12);
741  threesigA.SetFillColor(11);
742  threesigB.SetFillColor(11);
743  TLine line1(0.5 - radius, 0.5, 0.5 + radius, 0.5);
744  line1.SetLineWidth(3);
745 
746  TLine line2(0.5, 0.5 - radius, 0.5, 0.5 + radius);
747  line2.SetLineWidth(3);
748 
749  TArrow bf(0.5, 0.5, 0.5 + radius * cos(best_fit),0.5 + radius * sin(best_fit),0.04, "|>");
750  bf.SetLineWidth(3);
751  bf.SetLineColor(kRed);
752  bf.SetFillColor(kRed);
753 
754  TCanvas canvas("canvas", "canvas", 0, 0, 1000, 1000);
755  onesig.Draw();
756  twosigA.Draw();
757  twosigB.Draw();
758  threesigA.Draw();
759  threesigB.Draw();
760 
761  // Check if the rest wraps around the circle
762  if (bds[5] > 0) {
763  // Single rest slice
764  rest.Draw();
765  } else {
766  // Split rest into two slices
767  restA.Draw();
768  restB.Draw();
769  }
770 
771  line1.Draw();
772  line2.Draw();
773  bf.Draw();
774 
775  TLegend leg(0.0, 0.8, 0.23, 0.95);
776  leg.AddEntry(&bf, "Best Fit", "L");
777  leg.AddEntry(&onesig, "1#sigma", "F");
778  leg.AddEntry(&twosigA, "2#sigma", "F");
779  leg.AddEntry(&threesigA, "3#sigma", "F");
780  leg.Draw();
781 
782  // KS: Simple lambda to avoid copy-pasting
783  auto draw_text = [](auto& txt, Color_t color = kBlack) {
784  txt.SetTextAlign(22);
785  txt.SetTextColor(color);
786  txt.SetTextFont(43);
787  txt.SetTextSize(40);
788  txt.SetTextAngle(0);
789  txt.Draw();
790  };
791 
792  //KS: If best fit point is somehow very close text we simply not plot it
793  // Define a threshold for "too close"
794  constexpr double too_close_threshold = 0.1;
795 
796  // Position of tbf
797  const double tbf_x = 0.5 + (radius + 0.02) * cos(best_fit);
798  const double tbf_y = 0.5 + (radius + 0.02) * sin(best_fit);
799 
800  // Function to calculate distance between two points
801  auto distance = [](double x1, double y1, double x2, double y2) {
802  return std::sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2));
803  };
804 
805  // Check and draw right 0
806  constexpr double t0_x = 0.5 + radius + 0.02;
807  constexpr double t0_y = 0.5;
808  TText t0(t0_x, t0_y, "0");
809  if (distance(tbf_x, tbf_y, t0_x, t0_y) > too_close_threshold) {
810  draw_text(t0);
811  }
812 
813  // Check and draw left pi
814  constexpr double tp_x = 0.5 - radius - 0.02;
815  constexpr double tp_y = 0.5;
816  TLatex tp(tp_x, tp_y, "#pi");
817  if (distance(tbf_x, tbf_y, tp_x, tp_y) > too_close_threshold) {
818  draw_text(tp);
819  }
820 
821  // Check and draw top pi/2
822  constexpr double tp2_x = 0.5;
823  constexpr double tp2_y = 0.5 + radius + 0.04;
824  TLatex tp2(tp2_x, tp2_y, "#frac{#pi}{2}");
825  if (distance(tbf_x, tbf_y, tp2_x, tp2_y) > too_close_threshold) {
826  draw_text(tp2);
827  }
828 
829  // Check and draw bottom -pi/2
830  constexpr double tmp2_x = 0.5;
831  constexpr double tmp2_y = 0.5 - radius - 0.04;
832  TLatex tmp2(tmp2_x, tmp2_y, "-#frac{#pi}{2}");
833  if (distance(tbf_x, tbf_y, tmp2_x, tmp2_y) > too_close_threshold) {
834  draw_text(tmp2);
835  }
836 
837  TLatex tbf(0.5 + (radius + 0.02) * cos(best_fit),
838  0.5 + (radius + 0.02) * sin(best_fit),
839  fmt::format("{:.2f}", best_fit).c_str());
840  draw_text(tbf, kRed);
841 
842  canvas.Print(CanvasName);
843 }
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:140
ParameterEnum
Definition: MCMCProcessor.h:45
@ kXSecPar
Definition: MCMCProcessor.h:46
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
YAML::Node TMacroToYAML(const TMacro &macro)
KS: Convert a ROOT TMacro object to a YAML node.
Definition: YamlHelper.h:152
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Definition: YamlHelper.h:60
Class responsible for processing MCMC chains, performing diagnostics, generating plots,...
Definition: MCMCProcessor.h:58
void GetNthParameter(const int param, double &Prior, double &PriorError, TString &Title) const
Get properties of parameter by passing it number.
std::vector< std::vector< double > > ParamCentral
Parameters central values which we are going to analyse.
std::vector< std::vector< double > > ParamErrors
Uncertainty on a single parameter.
std::vector< int > nParam
Number of parameters per type.
std::vector< std::vector< bool > > ParamFlat
Whether Param has flat prior or not.
std::unique_ptr< TCanvas > Posterior
Fancy canvas used for our beautiful plots.
TFile * OutputFile
The output file.
TChain * Chain
Main chain storing all steps etc.
std::string MCMCFile
Name of MCMC file.
int nDraw
Number of all parameters used in the analysis.
std::vector< ParameterEnum > ParamType
Make an enum for which class this parameter belongs to so we don't have to keep string comparing.
TString CanvasName
Name of canvas which help to save to the sample pdf.
std::vector< std::string > ParameterGroup
std::vector< TString > BranchNames
std::vector< std::vector< TString > > ParamNames
Name of parameters which we are going to analyse.
std::vector< int > ParamTypeStartPos
int nEntries
KS: For merged chains number of entries will be different from nSteps.
void SavageDickeyPlot(std::unique_ptr< TH1D > &PriorHist, std::unique_ptr< TH1D > &PosteriorHist, const std::string &Title, const double EvaluationPoint) const
Produce Savage Dickey plot.
unsigned int BurnInCut
Value of burn in cut.
Custom exception class used throughout MaCh3.
double SamplePriorForParam(const int paramIndex, const std::unique_ptr< TRandom3 > &randGen, const std::vector< double > &FlatBounds) const
Draw Prior value.
int DeltaCPIndex
Index of in the parameter list.
Definition: OscProcessor.h:80
double CalcJarlskog(const double s2th13, const double s2th23, const double s2th12, const double dcp) const
Calculate Jarlskog Invariant using oscillation parameters.
std::string Sin2Theta13Name
Name of the parameter representing .
Definition: OscProcessor.h:63
std::string Sin2Theta12Name
Name of the parameter representing .
Definition: OscProcessor.h:65
int Sin2Theta12Index
Index of in the parameter list.
Definition: OscProcessor.h:76
OscProcessor(const std::string &InputFile)
Constructs an OscProcessor object with the specified input file and options.
int DeltaM2_23Index
Index of in the parameter list.
Definition: OscProcessor.h:82
std::string DeltaCPName
Name of the parameter representing (the CP-violating phase).
Definition: OscProcessor.h:69
bool PlotJarlskog
Will plot Jarlskog Invariant using information in the chain.
Definition: OscProcessor.h:57
int Sin2Theta13Index
Index of in the parameter list.
Definition: OscProcessor.h:74
std::string Sin2Theta23Name
Name of the parameter representing .
Definition: OscProcessor.h:67
void LoadAdditionalInfo() override
Read the Osc cov file and get the input central values and errors Here we allow Jarlskog Shenanigans.
void MakeJarlskogPlot(const std::unique_ptr< TH1D > &jarl, const std::unique_ptr< TH1D > &jarl_flatsindcp, const std::unique_ptr< TH1D > &jarl_NH, const std::unique_ptr< TH1D > &jarl_NH_flatsindcp, const std::unique_ptr< TH1D > &jarl_IH, const std::unique_ptr< TH1D > &jarl_IH_flatsindcp)
Perform Jarlskog Plotting.
int Sin2Theta23Index
Index of in the parameter list.
Definition: OscProcessor.h:78
virtual ~OscProcessor()
Destroys the OscProcessor object.
std::string DeltaM2_23Name
Name of the parameter representing (mass-squared difference).
Definition: OscProcessor.h:71
bool OscEnabled
Will plot Jarlskog Invariant using information in the chain.
Definition: OscProcessor.h:60
void MakePiePlot()
Make fancy Pie plot for delta CP.
void PerformJarlskogAnalysis()
Perform Several Jarlskog Plotting.
std::unique_ptr< ObjectType > Clone(const ObjectType *obj, const std::string &name="")
KS: Creates a copy of a ROOT-like object and wraps it in a smart pointer.
constexpr static const double _BAD_DOUBLE_
Default value used for double initialisation.
Definition: Core.h:53
TFile * Open(const std::string &Name, const std::string &Type, const std::string &File, const int Line)
Opens a ROOT file with the given name and mode.
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55
void PrintProgressBar(const Long64_t Done, const Long64_t All)
KS: Simply print progress bar.
Definition: Monitor.cpp:228
void EstimateDataTransferRate(TChain *chain, const Long64_t entry)
KS: Check what CPU you are using.
Definition: Monitor.cpp:211