MaCh3  2.5.1
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"
9 #include "TComplex.h"
11 
12 #pragma GCC diagnostic ignored "-Wfloat-conversion"
13 
14 // ****************************
15 OscProcessor::OscProcessor(const std::string &InputFile) : MCMCProcessor(InputFile) {
16 // ****************************
17  //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
18  PlotJarlskog = false;
19 
26 }
27 
28 // ****************************
29 // The destructor
31 // ****************************
32 }
33 
34 // ***************
35 // Read the Osc cov file and get the input central values and errors
37 // ***************
38  // KS: Check if OscParams were enabled, in future we will also get
39  for(size_t i = 0; i < ParameterGroup.size(); i++) {
40  if(ParameterGroup[i] == "Osc"){
41  OscEnabled = true;
42  break;
43  }
44  }
45 
46  if(OscEnabled)
47  {
48  for (int i = 0; i < nDraw; ++i)
49  {
50  //Those keep which parameter type we run currently and relative number
51  const int ParamEnum = ParamType[i];
52  const int ParamNo = i - ParamTypeStartPos[ParameterEnum(ParamEnum)];
53  const std::string CurrentName = ParamNames[ParamEnum][ParamNo].Data();
54 
56  if (CurrentName == "sin2th_13") {
57  Sin2Theta13Index = i;
58  Sin2Theta13Name = CurrentName;
59  } else if (CurrentName == "sin2th_12") {
60  Sin2Theta12Index = i;
61  Sin2Theta12Name = CurrentName;
62  } else if (CurrentName == "sin2th_23") {
63  Sin2Theta23Index = i;
64  Sin2Theta23Name = CurrentName;
65  } else if (CurrentName == "delta_cp") {
66  DeltaCPIndex = i;
67  DeltaCPName = CurrentName;
68  } else if (CurrentName == "delm2_23") {
69  DeltaM2_23Index = i;
70  DeltaM2_23Name = CurrentName;
71  }
72  }
73  } else{
74  MACH3LOG_WARN("Didn't find oscillation parameters");
75  }
76 
78  {
79  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)");
80  BranchNames.push_back("J_cp");
81  ParamType.push_back(kXSecPar);
82  nParam[kXSecPar]++;
83  nDraw++;
84 
86  ParamCentral[kXSecPar].push_back( 0. );
87  ParamErrors[kXSecPar].push_back( 1. );
88  // Push back the name
89  ParamNames[kXSecPar].push_back("J_cp");
90  ParamFlat[kXSecPar].push_back( false );
91  } else if(PlotJarlskog && !OscEnabled) {
92  MACH3LOG_ERROR("Trying to enable Jarlskog without oscillations");
93  throw MaCh3Exception(__FILE__,__LINE__);
94  }
95 }
96 
97 // ***************
98 // Calculate Jarlskog Invariant using oscillation parameters
99 double OscProcessor::CalcJarlskog(const double s2th13, const double s2th23, const double s2th12, const double dcp) const {
100 // ***************
101  const double s13 = std::sqrt(s2th13);
102  const double s23 = std::sqrt(s2th23);
103  const double s12 = std::sqrt(s2th12);
104  const double sdcp = std::sin(dcp);
105  const double c13 = std::sqrt(1.-s2th13);
106  const double c12 = std::sqrt(1.-s2th12);
107  const double c23 = std::sqrt(1.-s2th23);
108 
109  const double j = s13*c13*c13*s12*c12*s23*c23*sdcp;
110 
111  return j;
112 }
113 
114 // ***************
115 double OscProcessor::SamplePriorForParam(const int paramIndex, const std::unique_ptr<TRandom3>& randGen, const std::vector<double>& FlatBounds) const {
116 // ***************
117  TString Title = "";
118  double Prior = 1.0, PriorError = 1.0;
119  bool FlatPrior = false;
120 
121  // Get info for this parameter
122  GetNthParameter(paramIndex, Prior, PriorError, Title);
123 
124  ParameterEnum ParType = ParamType[paramIndex];
125  int ParamTemp = paramIndex - ParamTypeStartPos[ParType];
126  FlatPrior = ParamFlat[ParType][ParamTemp];
127 
128  if (FlatPrior) {
129  return randGen->Uniform(FlatBounds[0], FlatBounds[1]);
130  } else {
131  // Gaussian prior centered at Prior with width PriorError
132  return randGen->Gaus(Prior, PriorError);
133  }
134 }
135 
136 // ***************
137 void OscProcessor::Get1DReactorConstraintInfo(std::pair<double, double>& Sin13_NewPrior, bool& DoReweight) const {
138 // ***************
139  // Now read the MCMC file
140  TFile *TempFile = M3::Open((MCMCFile + ".root"), "open", __FILE__, __LINE__);
141 
142  // Get the settings for the MCMC
143  TMacro *Config = TempFile->Get<TMacro>("Reweight_Config");
144 
145  if (Config != nullptr) {
146  MACH3LOG_INFO("Found Reweight_Config in chain");
147 
148  // Print the reweight configuration for user info
149  YAML::Node Settings = TMacroToYAML(*Config);
150  // 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
151  if(CheckNodeExists(Settings, "ReweightMCMC")) {
152  YAML::Node firstReweight = Settings["ReweightMCMC"].begin()->second;
153  int dimension = GetFromManager<int>(firstReweight["ReweightDim"], 1);
154  std::string reweightType = GetFromManager<std::string>(firstReweight["ReweightType"], "");
155  auto paramNames = GetFromManager<std::vector<std::string>>(firstReweight["ReweightVar"], {});
156  if (dimension == 1 && reweightType == "Gaussian" && paramNames.size() == 1){
157  Sin13_NewPrior = Get<std::pair<double, double>>(firstReweight["ReweightPrior"],__FILE__,__LINE__);
158  DoReweight = true;
159  } else {
160  MACH3LOG_INFO("No valid reweighting configuration (1D Gaussian on sin2th_13 only) found for Jarlskog analysis");
161  }
162  } else {
163  MACH3LOG_INFO("No reweighting configuration found for Jarlskog analysis");
164  }
165  }
166 
167  TempFile->Close();
168  delete TempFile;
169 }
170 
171 // ***************
172 // Perform Several Jarlskog Plotting
174 // ***************
175  if(!OscEnabled ||
181  {
182  MACH3LOG_WARN("Will not {}, as oscillation parameters are missing", __func__);
183  return;
184  }
185  MACH3LOG_INFO("Starting {}", __func__);
186 
187  double s2th13, s2th23, s2th12, dcp, dm2 = M3::_BAD_DOUBLE_;
188  double weight = 1.0;
189 
190  bool DoReweight = false;
191  std::pair<double, double> Sin13_NewPrior;
192  Get1DReactorConstraintInfo(Sin13_NewPrior, DoReweight);
193 
194  TDirectory *JarlskogDir = OutputFile->mkdir("Jarlskog");
195  JarlskogDir->cd();
196 
197  unsigned int step = 0;
198  Chain->SetBranchStatus("*", false);
199 
200  Chain->SetBranchStatus(Sin2Theta13Name.c_str(), true);
201  Chain->SetBranchAddress(Sin2Theta13Name.c_str(), &s2th13);
202 
203  Chain->SetBranchStatus(Sin2Theta23Name.c_str(), true);
204  Chain->SetBranchAddress(Sin2Theta23Name.c_str(), &s2th23);
205 
206  Chain->SetBranchStatus(Sin2Theta12Name.c_str(), true);
207  Chain->SetBranchAddress(Sin2Theta12Name.c_str(), &s2th12);
208 
209  Chain->SetBranchStatus(DeltaCPName.c_str(), true);
210  Chain->SetBranchAddress(DeltaCPName.c_str(), &dcp);
211 
212  Chain->SetBranchStatus(DeltaM2_23Name.c_str(), true);
213  Chain->SetBranchAddress(DeltaM2_23Name.c_str(), &dm2);
214 
215  Chain->SetBranchStatus("step", true);
216  Chain->SetBranchAddress("step", &step);
217 
218  if(DoReweight) {
219  Chain->SetBranchStatus("Weight", true);
220  Chain->SetBranchAddress("Weight", &weight);
221  } else {
222  MACH3LOG_WARN("Not applying reweighting weight");
223  weight = 1.0;
224  }
225 
226  // Original histograms
227  auto jarl = std::make_unique<TH1D>("jarl", "jarl", 1000, -0.05, 0.05);
228  jarl->SetDirectory(nullptr);
229  auto jarl_th23 = std::make_unique<TH2D>("jarl_th23", "jarl_th23", 500, -0.05, 0.05, 500, 0.3, 0.7);
230  jarl_th23->SetDirectory(nullptr);
231  auto jarl_dcp = std::make_unique<TH2D>("jarl_dcp", "jarl_dcp", 500, -0.05, 0.05, 500, -1. * TMath::Pi(), TMath::Pi());
232  jarl_dcp->SetDirectory(nullptr);
233 
234  jarl->SetTitle("Jarlskog Invariant;J #equiv s_{13}c_{13}^{2}s_{12}c_{12}s_{23}c_{23}sin#delta;Posterior probability");
235  jarl_th23->SetTitle("Jarlskog Invariant;J #equiv s_{13}c_{13}^{2}s_{12}c_{12}s_{23}c_{23}sin#delta;Posterior probability");
236 
237  // Clones
238  auto jarl_IH = M3::Clone(jarl.get(), "jarl_IH");
239  auto jarl_NH = M3::Clone(jarl.get(), "jarl_NH");
240 
241  auto jarl_th23_IH = M3::Clone(jarl_th23.get(), "jarl_th23_IH");
242  auto jarl_th23_NH = M3::Clone(jarl_th23.get(), "jarl_th23_NH");
243 
244  auto jarl_dcp_IH = M3::Clone(jarl_dcp.get(), "jarl_dcp_IH");
245  auto jarl_dcp_NH = M3::Clone(jarl_dcp.get(), "jarl_dcp_NH");
246 
247  auto jarl_flatsindcp = M3::Clone(jarl.get(), "jarl_flatsindcp");
248  auto jarl_IH_flatsindcp = M3::Clone(jarl.get(), "jarl_IH_flatsindcp");
249  auto jarl_NH_flatsindcp = M3::Clone(jarl.get(), "jarl_NH_flatsindcp");
250 
251  auto jarl_th23_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_flatsindcp");
252  auto jarl_th23_IH_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_IH_flatsindcp");
253  auto jarl_th23_NH_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_NH_flatsindcp");
254 
255  auto jarl_prior = M3::Clone(jarl.get(), "jarl_prior");
256  auto jarl_prior_flatsindcp = M3::Clone(jarl.get(), "jarl_prior_flatsindcp");
257  std::unique_ptr<TH1D> jarl_wRC_prior, jarl_wRC_prior_flatsindcp, jarl_wRC_prior_t2kth23;
258  // Only use this if chain has reweigh weight [mostly coming from Reactor Constrains]
259  if(DoReweight){
260  jarl_wRC_prior = M3::Clone(jarl.get(), "jarl_wRC_prior");
261  jarl_wRC_prior_flatsindcp = M3::Clone(jarl.get(), "jarl_wRC_prior_flatsindcp");
262  jarl_wRC_prior_t2kth23 = M3::Clone(jarl.get(), "jarl_wRC_prior_flatsindcp");
263  }
264 
265  // to apply a prior that is flat in sin(dcp) intead of dcp
266  auto prior3 = std::make_unique<TF1>("prior3", "TMath::Abs(TMath::Cos(x))");
267 
268  // T2K prior is flat (and uncorrelated) in dcp, sin^2(th13), sin^2(th23)
269  auto randGen = std::make_unique<TRandom3>(0);
270  const Long64_t countwidth = nEntries/5;
271 
272  for(int i = 0; i < nEntries; ++i) {
273  if (i % countwidth == 0) {
276  } else {
277  Chain->GetEntry(i);
278  }
279 
280  if(step < BurnInCut) continue; // burn-in cut
281 
282  const double j = CalcJarlskog(s2th13, s2th23, s2th12, dcp);
283  const double prior_weight = prior3->Eval(dcp);
284 
285  jarl->Fill(j, weight);
286  jarl_th23->Fill(j, s2th23, weight);
287  jarl_dcp->Fill(j, dcp, weight);
288 
289  jarl_flatsindcp->Fill(j, prior_weight*weight);
290  jarl_th23_flatsindcp->Fill(j, s2th23, prior_weight*weight);
291 
292  const double prior_s2th13 = SamplePriorForParam(Sin2Theta13Index, randGen, {0.,1.});
293  const double prior_s2th23 = SamplePriorForParam(Sin2Theta23Index, randGen, {0.,1.});
294  const double prior_s2th12 = SamplePriorForParam(Sin2Theta12Index, randGen, {0.,1.});
295  const double prior_dcp = SamplePriorForParam(DeltaCPIndex, randGen, {-1.*TMath::Pi(),TMath::Pi()});
296  // KS: This is hardcoded but we always assume flat in delta CP so probably fine
297  const double prior_sindcp = randGen->Uniform(-1., 1.);
298 
299  const double prior_s13 = std::sqrt(prior_s2th13);
300  const double prior_s23 = std::sqrt(prior_s2th23);
301  const double prior_s12 = std::sqrt(prior_s2th12);
302  const double prior_sdcp = std::sin(prior_dcp);
303  const double prior_c13 = std::sqrt(1.-prior_s2th13);
304  const double prior_c12 = std::sqrt(1.-prior_s2th12);
305  const double prior_c23 = std::sqrt(1.-prior_s2th23);
306  const double prior_j = prior_s13*prior_c13*prior_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sdcp;
307  const double prior_flatsindcp_j = prior_s13*prior_c13*prior_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sindcp;
308 
309  jarl_prior->Fill(prior_j);
310  jarl_prior_flatsindcp->Fill(prior_flatsindcp_j);
311 
312  if(DoReweight) {
313  const double prior_wRC_s2th13 = randGen->Gaus(Sin13_NewPrior.first, Sin13_NewPrior.second);
314  const double prior_wRC_s13 = std::sqrt(prior_wRC_s2th13);
315  const double prior_wRC_c13 = std::sqrt(1.-prior_wRC_s2th13);
316  const double prior_wRC_j = prior_wRC_s13*prior_wRC_c13*prior_wRC_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sdcp;
317  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;
318  const double s23 = std::sqrt(s2th23);
319  const double c23 = std::sqrt(1.-s2th23);
320 
321  jarl_wRC_prior->Fill(prior_wRC_j);
322  jarl_wRC_prior_flatsindcp->Fill(prior_wRC_flatsindcp_j);
323  jarl_wRC_prior_t2kth23->Fill(prior_wRC_s13*prior_wRC_c13*prior_wRC_c13*prior_s12*prior_c12*s23*c23*prior_sdcp);
324  }
325 
326  if(dm2 > 0.) {
327  jarl_NH->Fill(j, weight);
328  jarl_th23_NH->Fill(j, s2th23, weight);
329  jarl_dcp_NH->Fill(j, dcp, weight);
330  jarl_NH_flatsindcp->Fill(j, prior_weight*weight);
331  jarl_th23_NH_flatsindcp->Fill(j, s2th23, prior_weight*weight);
332  }
333  else if(dm2 < 0.) {
334  jarl_IH->Fill(j, weight);
335  jarl_th23_IH->Fill(j, s2th23, weight);
336  jarl_dcp_IH->Fill(j, dcp, weight);
337  jarl_IH_flatsindcp->Fill(j, prior_weight*weight);
338  jarl_th23_IH_flatsindcp->Fill(j, s2th23, prior_weight*weight);
339  }
340  }
341 
342  jarl->Write("jarlskog_both");
343  jarl_NH->Write("jarlskog_NH");
344  jarl_IH->Write("jarlskog_IH");
345  jarl_th23->Write("jarlskog_th23_both");
346  jarl_th23_NH->Write("jarlskog_th23_NH");
347  jarl_th23_IH->Write("jarlskog_th23_IH");
348 
349  jarl_dcp->Write("jarlskog_dcp_both");
350  jarl_dcp_NH->Write("jarlskog_dcp_NH");
351  jarl_dcp_IH->Write("jarlskog_dcp_IH");
352 
353 
354  jarl_flatsindcp->Write("jarlskog_both_flatsindcp");
355  jarl_NH_flatsindcp->Write("jarlskog_NH_flatsindcp");
356  jarl_IH_flatsindcp->Write("jarlskog_IH_flatsindcp");
357  jarl_th23_flatsindcp->Write("jarlskog_th23_both_flatsindcp");
358  jarl_th23_NH_flatsindcp->Write("jarlskog_th23_NH_flatsindcp");
359  jarl_th23_IH_flatsindcp->Write("jarlskog_th23_IH_flatsindcp");
360 
361  jarl_prior->Write("jarl_prior");
362  jarl_prior_flatsindcp->Write("jarl_prior_flatsindcp");
363  if(DoReweight) {
364  jarl_wRC_prior->Write("jarl_wRC_prior");
365  jarl_wRC_prior_flatsindcp->Write("jarl_wRC_prior_flatsindcp");
366  jarl_wRC_prior_t2kth23->Write("jarl_wRC_prior_t2kth23");
367  }
368 
369  MakeJarlskogPlot(jarl, jarl_flatsindcp,
370  jarl_NH, jarl_NH_flatsindcp,
371  jarl_IH, jarl_IH_flatsindcp);
372 
373  // Perform Savage Dickey analysis
374  if(DoReweight) {
375  SavageDickeyPlot(jarl, jarl_wRC_prior, "Jarlskog flat #delta_{CP}", 0);
376  SavageDickeyPlot(jarl_flatsindcp, jarl_wRC_prior_flatsindcp, "Jarlskog flat sin#delta_{CP}", 0);
377  } else {
378  SavageDickeyPlot(jarl, jarl_prior, "Jarlskog flat #delta_{CP}", 0);
379  SavageDickeyPlot(jarl_flatsindcp, jarl_prior_flatsindcp, "Jarlskog flat sin#delta_{CP}", 0);
380  }
381 
382  JarlskogDir->Close();
383  delete JarlskogDir;
384 
385  Chain->SetBranchStatus("*", true);
386  OutputFile->cd();
387 }
388 
389 
390 // ***************
391 void OscProcessor::MakeJarlskogPlot(const std::unique_ptr<TH1D>& jarl,
392  const std::unique_ptr<TH1D>& jarl_flatsindcp,
393  const std::unique_ptr<TH1D>& jarl_NH,
394  const std::unique_ptr<TH1D>& jarl_NH_flatsindcp,
395  const std::unique_ptr<TH1D>& jarl_IH,
396  const std::unique_ptr<TH1D>& jarl_IH_flatsindcp) {
397 // ***************
398  MACH3LOG_INFO("Starting {}", __func__);
399  int originalErrorLevel = gErrorIgnoreLevel;
400  gErrorIgnoreLevel = kFatal;
401 
402  // 1-->NH, 0-->both, -1-->IH
403  for(int hierarchy = -1; hierarchy <= 1; hierarchy++)
404  {
405  std::unique_ptr<TH1D> j_hist;
406  std::unique_ptr<TH1D> j_hist_sdcp;
407  if(hierarchy == 1) {
408  j_hist = M3::Clone(jarl_NH.get(), "");
409  j_hist_sdcp = M3::Clone(jarl_NH_flatsindcp.get(), "");
410  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
411  } else if(hierarchy == 0) {
412  j_hist = M3::Clone(jarl.get(), "");
413  j_hist_sdcp = M3::Clone(jarl_flatsindcp.get(), "");
414  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
415  } else if(hierarchy == -1) {
416  j_hist = M3::Clone(jarl_IH.get(), "");
417  j_hist_sdcp = M3::Clone(jarl_IH_flatsindcp.get(), "");
418  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
419  } else {
420  MACH3LOG_ERROR("Invalid hierarchy option. 1 for NH, 0 for both, -1 for IH");
421  throw MaCh3Exception(__FILE__ , __LINE__ );
422  }
423 
424  j_hist->Rebin(7);
425  j_hist_sdcp->Rebin(7);
426 
427  j_hist->SetLineColor(kAzure-2);
428  j_hist_sdcp->SetLineColor(kOrange+1);
429  j_hist->SetLineWidth(2);
430  j_hist_sdcp->SetLineWidth(2);
431 
432  auto StyleAxis = [](TH1* h) {
433  auto xAxis = h->GetXaxis();
434  auto yAxis = h->GetYaxis();
435 
436  xAxis->SetLabelSize(0.04);
437  xAxis->SetLabelFont(132);
438  xAxis->SetTitleSize(0.04);
439  xAxis->SetTitleOffset(0.80);
440  xAxis->SetTitleFont(132);
441  xAxis->SetNdivisions(505);
442  xAxis->SetTickSize(0.04);
443 
444  yAxis->SetLabelSize(0.04);
445  yAxis->SetLabelFont(132);
446  yAxis->SetTitleSize(0.04);
447  yAxis->SetTitleOffset(1.2);
448  yAxis->SetTitleFont(132);
449  yAxis->SetNdivisions(505);
450  yAxis->SetTickSize(0.04);
451  };
452 
453  StyleAxis(j_hist.get());
454 
455  j_hist->GetXaxis()->SetRangeUser(-0.04,0.04);
456  j_hist->Scale(1./j_hist->Integral());
457  j_hist_sdcp->Scale(1./j_hist_sdcp->Integral());
458 
459  std::unique_ptr<TH1D> j_hist_copy = M3::Clone(j_hist.get(), "j_hist_copy");
460  std::unique_ptr<TH1D> j_hist_1sig = M3::Clone(j_hist.get(), "j_hist_1sig");
461  std::unique_ptr<TH1D> j_hist_2sig = M3::Clone(j_hist.get(), "j_hist_2sig");
462  std::unique_ptr<TH1D> j_hist_3sig = M3::Clone(j_hist.get(), "j_hist_3sig");
463 
464  //upper and lower edges
465  double j_bf = j_hist_copy->GetXaxis()->GetBinCenter(j_hist_copy->GetMaximumBin());
466  double j_1sig_low = 9999999.;
467  double j_1sig_up = -9999999.;
468  double j_2sig_low = 9999999.;;
469  double j_2sig_up = -9999999.;
470  double j_3sig_low = 9999999.;;
471  double j_3sig_up = -9999999.;
472 
473 
474  std::unique_ptr<TH1D> j_hist_sdcp_copy = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_copy");
475  std::unique_ptr<TH1D> j_hist_sdcp_1sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_1sig");
476  std::unique_ptr<TH1D> j_hist_sdcp_2sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_2sig");
477  std::unique_ptr<TH1D> j_hist_sdcp_3sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_3sig");
478 
479  //upper and lower edges
480  double j_sdcp_1sig_low = 9999999.;
481  double j_sdcp_1sig_up = -9999999.;
482  double j_sdcp_2sig_low = 9999999.;;
483  double j_sdcp_2sig_up = -9999999.;
484  double j_sdcp_3sig_low = 9999999.;;
485  double j_sdcp_3sig_up = -9999999.;
486 
487  double contlevel1 = 0.68;
488  double contlevel2 = 0.90;
489  double contlevel4 = 0.99;
490  double contlevel5 = 0.9973;
491  double integral, tsum = 0.;
492 
493  integral = j_hist_copy->Integral();
494 
495  while((tsum/integral)<contlevel5) {
496  double tmax = j_hist_copy->GetMaximum();
497  int bin = j_hist_copy->GetMaximumBin();
498  double xval = j_hist_copy->GetXaxis()->GetBinCenter(bin);
499  double xwidth = j_hist_copy->GetXaxis()->GetBinWidth(bin);
500  if((tsum/integral)<contlevel1) {
501  j_hist_copy->SetBinContent(bin,-1.0);
502  j_hist_1sig->SetBinContent(bin,0.);
503  j_hist_2sig->SetBinContent(bin,0.);
504  j_hist_3sig->SetBinContent(bin,0.);
505  if(xval<j_1sig_low && xval<j_bf) j_1sig_low = xval - xwidth/2.;
506  if(xval>j_1sig_up && xval>j_bf) j_1sig_up = xval + xwidth/2.;
507  }
508  if((tsum/integral)<contlevel2 && (tsum / integral > contlevel1) ) {
509  j_hist_copy->SetBinContent(bin,-5.0);
510  j_hist_2sig->SetBinContent(bin,0.);
511  j_hist_3sig->SetBinContent(bin,0.);
512  if(xval<j_2sig_low && xval<j_bf) j_2sig_low = xval - xwidth/2.;
513  if(xval>j_2sig_up && xval>j_bf) j_2sig_up = xval + xwidth/2.;
514  }
515  if((tsum/integral)<contlevel4 && (tsum / integral > contlevel1) ) {
516  j_hist_copy->SetBinContent(bin,-9.0);
517  j_hist_3sig->SetBinContent(bin,0.);
518  if(xval < j_3sig_low && xval <j_bf) j_3sig_low = xval - xwidth/2.;
519  if(xval > j_3sig_up && xval > j_bf) j_3sig_up = xval + xwidth/2.;
520  }
521  tsum+=tmax;
522  }
523 
524  integral = j_hist_sdcp_copy->Integral();
525  tsum = 0.;
526 
527  while((tsum/integral)<contlevel5) {
528  double tmax = j_hist_sdcp_copy->GetMaximum();
529  int bin = j_hist_sdcp_copy->GetMaximumBin();
530  double xval = j_hist_sdcp_copy->GetXaxis()->GetBinCenter(bin);
531  double xwidth = j_hist_sdcp_copy->GetXaxis()->GetBinWidth(bin);
532  if((tsum/integral)<contlevel1) {
533  j_hist_sdcp_copy->SetBinContent(bin,-1.0);
534  j_hist_sdcp_1sig->SetBinContent(bin,0.);
535  j_hist_sdcp_2sig->SetBinContent(bin,0.);
536  j_hist_sdcp_3sig->SetBinContent(bin,0.);
537  if(xval<j_sdcp_1sig_low && xval<j_bf) j_sdcp_1sig_low = xval - xwidth/2.;
538  if(xval>j_sdcp_1sig_up && xval>j_bf) j_sdcp_1sig_up = xval + xwidth/2.;
539  }
540  if((tsum/integral)<contlevel2 && (tsum / integral > contlevel1) ) {
541  j_hist_sdcp_copy->SetBinContent(bin,-5.0);
542  j_hist_sdcp_2sig->SetBinContent(bin,0.);
543  j_hist_sdcp_3sig->SetBinContent(bin,0.);
544  if(xval<j_sdcp_2sig_low && xval<j_bf) j_sdcp_2sig_low = xval - xwidth/2.;
545  if(xval>j_sdcp_2sig_up && xval>j_bf) j_sdcp_2sig_up = xval + xwidth/2.;
546  }
547  if((tsum/integral)<contlevel4 && (tsum / integral > contlevel1) ) {
548  j_hist_sdcp_copy->SetBinContent(bin,-9.0);
549  j_hist_sdcp_3sig->SetBinContent(bin,0.);
550  if(xval<j_sdcp_3sig_low && xval<j_bf) j_sdcp_3sig_low = xval - xwidth/2.;
551  if(xval>j_sdcp_3sig_up && xval>j_bf) j_sdcp_3sig_up = xval + xwidth/2.;
552  }
553  tsum+=tmax;
554  }
555 
556  j_hist_1sig->SetLineStyle(9);
557  j_hist_sdcp_1sig->SetLineStyle(9);
558  j_hist_2sig->SetLineStyle(7);
559  j_hist_sdcp_2sig->SetLineStyle(7);
560  j_hist_3sig->SetLineStyle(2);
561  j_hist_sdcp_3sig->SetLineStyle(2);
562 
563  auto ldash = std::make_unique<TH1D>("ldash", "solid", 10, -0.04, 0.04);
564  auto sdash = std::make_unique<TH1D>("sdash", "dashed", 10, -0.04, 0.04);
565  auto fdash = std::make_unique<TH1D>("fdash", "fdashed",10, -0.04, 0.04);
566  ldash->SetLineColor(kBlack);
567  sdash->SetLineColor(kBlack);
568  fdash->SetLineColor(kBlack);
569  ldash->SetLineWidth(2);
570  sdash->SetLineWidth(2);
571  fdash->SetLineWidth(2);
572  ldash->SetLineStyle(9);
573  sdash->SetLineStyle(7);
574  fdash->SetLineStyle(2);
575 
576  double vertUp = 0.5 * j_hist->GetMaximum();
577  auto jline_1sig_low = std::make_unique<TLine>(j_1sig_low, 0., j_1sig_low, vertUp);
578  auto jline_2sig_low = std::make_unique<TLine>(j_2sig_low, 0., j_2sig_low, vertUp);
579  auto jline_3sig_low = std::make_unique<TLine>(j_3sig_low, 0., j_3sig_low, vertUp);
580 
581  auto jline_1sig_up = std::make_unique<TLine>(j_1sig_up, 0., j_1sig_up,vertUp);
582  auto jline_2sig_up = std::make_unique<TLine>(j_2sig_up, 0., j_2sig_up,vertUp);
583  auto jline_3sig_up = std::make_unique<TLine>(j_3sig_up, 0., j_3sig_up,vertUp);
584 
585  auto jline_sdcp_1sig_low = std::make_unique<TLine>(j_sdcp_1sig_low, 0., j_sdcp_1sig_low, vertUp);
586  auto jline_sdcp_2sig_low = std::make_unique<TLine>(j_sdcp_2sig_low, 0., j_sdcp_2sig_low, vertUp);
587  auto jline_sdcp_3sig_low = std::make_unique<TLine>(j_sdcp_3sig_low, 0., j_sdcp_3sig_low, vertUp);
588 
589  auto jline_sdcp_1sig_up = std::make_unique<TLine>(j_sdcp_1sig_up, 0., j_sdcp_1sig_up, vertUp);
590  auto jline_sdcp_2sig_up = std::make_unique<TLine>(j_sdcp_2sig_up, 0., j_sdcp_2sig_up, vertUp);
591  auto jline_sdcp_3sig_up = std::make_unique<TLine>(j_sdcp_3sig_up, 0., j_sdcp_3sig_up, vertUp);
592 
593  double arrowLength = 0.003;
594  double arrowHeight = vertUp;
595 
596  auto MakeArrow = [&](double x, Color_t color, Width_t width) -> std::unique_ptr<TArrow> {
597  auto arrow = std::make_unique<TArrow>(x, arrowHeight, x - arrowLength, arrowHeight, 0.02, ">");
598  arrow->SetLineColor(color);
599  arrow->SetLineWidth(width);
600  return arrow;
601  };
602 
603  auto j_arrow_1sig_up = MakeArrow(j_1sig_up, j_hist_1sig->GetLineColor(), j_hist_1sig->GetLineWidth());
604  auto j_arrow_2sig_up = MakeArrow(j_2sig_up, j_hist_2sig->GetLineColor(), j_hist_2sig->GetLineWidth());
605  auto j_arrow_3sig_up = MakeArrow(j_3sig_up, j_hist_3sig->GetLineColor(), j_hist_3sig->GetLineWidth());
606 
607  auto j_sdcp_arrow_1sig_up = MakeArrow(j_sdcp_1sig_up, j_hist_sdcp_1sig->GetLineColor(), j_hist_sdcp_1sig->GetLineWidth());
608  auto j_sdcp_arrow_2sig_up = MakeArrow(j_sdcp_2sig_up, j_hist_sdcp_2sig->GetLineColor(), j_hist_sdcp_2sig->GetLineWidth());
609  auto j_sdcp_arrow_3sig_up = MakeArrow(j_sdcp_3sig_up, j_hist_sdcp_3sig->GetLineColor(), j_hist_sdcp_3sig->GetLineWidth());
610 
611  MACH3LOG_DEBUG("j_1sig_low = {:.4f}, j_2sig_low = {:.4f}, j_3sig_low = {:.4f}", j_1sig_low, j_2sig_low, j_3sig_low);
612  MACH3LOG_DEBUG("j_1sig_up = {:.4f}, j_2sig_up = {:.4f}, j_3sig_up = {:.4f}", j_1sig_up, j_2sig_up, j_3sig_up);
613 
614  auto CopyLineStyle = [](const TH1D* src, TLine* dst) {
615  dst->SetLineColor(src->GetLineColor());
616  dst->SetLineStyle(src->GetLineStyle());
617  dst->SetLineWidth(src->GetLineWidth());
618  };
619 
620  CopyLineStyle(j_hist_1sig.get(), jline_1sig_low.get());
621  CopyLineStyle(j_hist_1sig.get(), jline_1sig_up.get());
622  CopyLineStyle(j_hist_2sig.get(), jline_2sig_low.get());
623  CopyLineStyle(j_hist_2sig.get(), jline_2sig_up.get());
624  CopyLineStyle(j_hist_3sig.get(), jline_3sig_low.get());
625  CopyLineStyle(j_hist_3sig.get(), jline_3sig_up.get());
626 
627  CopyLineStyle(j_hist_sdcp_1sig.get(), jline_sdcp_1sig_low.get());
628  CopyLineStyle(j_hist_sdcp_1sig.get(), jline_sdcp_1sig_up.get());
629  CopyLineStyle(j_hist_sdcp_2sig.get(), jline_sdcp_2sig_low.get());
630  CopyLineStyle(j_hist_sdcp_2sig.get(), jline_sdcp_2sig_up.get());
631  CopyLineStyle(j_hist_sdcp_3sig.get(), jline_sdcp_3sig_low.get());
632  CopyLineStyle(j_hist_sdcp_3sig.get(), jline_sdcp_3sig_up.get());
633 
634  auto leg = std::make_unique<TLegend>(0.45, 0.60, 0.75, 0.90);
635  leg->SetTextSize(0.05);
636  leg->SetFillStyle(0);
637  leg->SetNColumns(1);
638  leg->SetTextFont(132);
639  leg->SetBorderSize(0);
640 
641  leg->AddEntry(j_hist.get(), "Prior flat in #delta_{CP}", "l");
642  leg->AddEntry(j_hist_sdcp.get(), "Prior flat in sin#delta_{CP}", "l");
643  leg->AddEntry(ldash.get(), "68% CI", "l");
644  leg->AddEntry(sdash.get(), "90% CI", "l");
645  leg->AddEntry(fdash.get(), "99% CI", "l");
646 
647  j_hist->GetYaxis()->SetRangeUser(0., j_hist->GetMaximum()*1.15);
648  j_hist->Draw("h");
649  j_hist_sdcp->Draw("same h");
650 
651  jline_sdcp_1sig_up->Draw("same");
652  jline_sdcp_2sig_up->Draw("same");
653  jline_sdcp_3sig_up->Draw("same");
654  jline_1sig_up->Draw("same");
655  jline_2sig_up->Draw("same");
656  jline_3sig_up->Draw("same");
657 
658  j_arrow_1sig_up->Draw();
659  j_arrow_2sig_up->Draw();
660  j_arrow_3sig_up->Draw();
661  j_sdcp_arrow_1sig_up->Draw();
662  j_sdcp_arrow_2sig_up->Draw();
663  j_sdcp_arrow_3sig_up->Draw();
664  leg->Draw("same");
665 
666  auto ttext = std::make_unique<TText>();
667  ttext->SetNDC(); // Use normalized device coordinates
668  ttext->SetTextSize(0.03); // Adjust size as needed
669  ttext->SetTextAlign(13); // Align left-top
670 
671  if (hierarchy == 1) ttext->DrawText(0.15, 0.85, "Normal Ordering");
672  else if (hierarchy == 0) ttext->DrawText(0.15, 0.85, "Both Orderings");
673  else if (hierarchy == -1) ttext->DrawText(0.15, 0.85, "Inverted Ordering");
674 
675  gPad->RedrawAxis();
676  Posterior->Update();
677  gPad->Update();
678 
679  Posterior->Print(CanvasName);
680 
681  if(hierarchy == 1) Posterior->Write("jarl1D_NH_comp");
682  else if(hierarchy == 0) Posterior->Write("jarl1D_both_comp");
683  else if(hierarchy == -1) Posterior->Write("jarl1D_IH_comp");
684  }
685 
686  gErrorIgnoreLevel = originalErrorLevel;
687 }
688 
689 // ***************
691 // ***************
693  {
694  MACH3LOG_WARN("Will not {}, as oscillation parameters are missing", __func__);
695  return;
696  }
697  MACH3LOG_INFO("Starting {}", __func__);
698 
699  // get best fit for delta CP
700  const double best_fit = (*Means_HPD)(DeltaCPIndex);
701 
702  const double sigma_p = (*Errors_HPD_Positive)(DeltaCPIndex);
703  const double sigma_n = (*Errors_HPD_Negative)(DeltaCPIndex);
704  // make sure result is between -pi and pi
705  auto wrap_pi = [](double x) {
706  while (x > TMath::Pi()) x -= 2*TMath::Pi();
707  while (x < -TMath::Pi()) x += 2*TMath::Pi();
708  return x;
709  };
710 
711  std::array<double, 6> bds;
712  bds[0] = wrap_pi(best_fit - 3.0 * sigma_n); // -3σ
713  bds[1] = wrap_pi(best_fit - 2.0 * sigma_n); // -2σ
714  bds[2] = wrap_pi(best_fit - 1.0 * sigma_n); // -1σ
715  bds[3] = wrap_pi(best_fit + 1.0 * sigma_p); // +1σ
716  bds[4] = wrap_pi(best_fit + 2.0 * sigma_p); // +2σ
717  bds[5] = wrap_pi(best_fit + 3.0 * sigma_p); // +3σ
718 
719  constexpr double radius = 0.4;
720  constexpr double rad_to_deg = 180.0 / TMath::Pi();
721 
722  // ROOT expects TEllipse angles in degrees, counterclockwise from the x-axis.
723  // If phimax < phimin, ROOT draws counterclockwise across the full circle, causing overlaps.
724  // This ensures threesigA slice stays within the intended range.
725  auto normalize_angle = [](double rad) {
726  // If rad is negative, add 2*pi to wrap into [0, 2*pi)
727  if (rad < 0) rad += 2.0 * TMath::Pi();
728  return rad;
729  };
730 
731  TEllipse onesig (0.5, 0.5, radius, radius, bds[2] * rad_to_deg, bds[4] * rad_to_deg);
732  TEllipse twosigA (0.5, 0.5, radius, radius, bds[1] * rad_to_deg, bds[2] * rad_to_deg);
733  TEllipse twosigB (0.5, 0.5, radius, radius, bds[3] * rad_to_deg, bds[4] * rad_to_deg);
734 
735  // three sigma slices
736  TEllipse threesigA(0.5, 0.5, radius, radius, bds[0] * rad_to_deg, normalize_angle(bds[1]) * rad_to_deg);
737  TEllipse threesigB(0.5, 0.5, radius, radius, bds[4] * rad_to_deg, bds[5] * rad_to_deg);
738 
739  // Remaining slices
740  TEllipse rest(0.5, 0.5, radius, radius, bds[5]*rad_to_deg, bds[0]*rad_to_deg);
741  TEllipse restA(0.5, 0.5, radius, radius, bds[5]*rad_to_deg, 180.0);
742  TEllipse restB(0.5, 0.5, radius, radius, -180.0, bds[0]*rad_to_deg);
743 
744  onesig.SetFillColor(13);
745  twosigA.SetFillColor(12);
746  twosigB.SetFillColor(12);
747  threesigA.SetFillColor(11);
748  threesigB.SetFillColor(11);
749  TLine line1(0.5 - radius, 0.5, 0.5 + radius, 0.5);
750  line1.SetLineWidth(3);
751 
752  TLine line2(0.5, 0.5 - radius, 0.5, 0.5 + radius);
753  line2.SetLineWidth(3);
754 
755  TArrow bf(0.5, 0.5, 0.5 + radius * cos(best_fit),0.5 + radius * sin(best_fit),0.04, "|>");
756  bf.SetLineWidth(3);
757  bf.SetLineColor(kRed);
758  bf.SetFillColor(kRed);
759 
760  TCanvas canvas("canvas", "canvas", 0, 0, 1000, 1000);
761  onesig.Draw();
762  twosigA.Draw();
763  twosigB.Draw();
764  threesigA.Draw();
765  threesigB.Draw();
766 
767  // Check if the rest wraps around the circle
768  if (bds[5] > 0) {
769  // Single rest slice
770  rest.Draw();
771  } else {
772  // Split rest into two slices
773  restA.Draw();
774  restB.Draw();
775  }
776 
777  line1.Draw();
778  line2.Draw();
779  bf.Draw();
780 
781  TLegend leg(0.0, 0.8, 0.23, 0.95);
782  leg.AddEntry(&bf, "Best Fit", "L");
783  leg.AddEntry(&onesig, "1#sigma", "F");
784  leg.AddEntry(&twosigA, "2#sigma", "F");
785  leg.AddEntry(&threesigA, "3#sigma", "F");
786  leg.Draw();
787 
788  // KS: Simple lambda to avoid copy-pasting
789  auto draw_text = [](auto& txt, Color_t color = kBlack) {
790  txt.SetTextAlign(22);
791  txt.SetTextColor(color);
792  txt.SetTextFont(43);
793  txt.SetTextSize(40);
794  txt.SetTextAngle(0);
795  txt.Draw();
796  };
797 
798  //KS: If best fit point is somehow very close text we simply not plot it
799  // Define a threshold for "too close"
800  constexpr double too_close_threshold = 0.1;
801 
802  // Position of tbf
803  const double tbf_x = 0.5 + (radius + 0.02) * cos(best_fit);
804  const double tbf_y = 0.5 + (radius + 0.02) * sin(best_fit);
805 
806  // Function to calculate distance between two points
807  auto distance = [](double x1, double y1, double x2, double y2) {
808  return std::sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2));
809  };
810 
811  // Check and draw right 0
812  constexpr double t0_x = 0.5 + radius + 0.02;
813  constexpr double t0_y = 0.5;
814  TText t0(t0_x, t0_y, "0");
815  if (distance(tbf_x, tbf_y, t0_x, t0_y) > too_close_threshold) {
816  draw_text(t0);
817  }
818 
819  // Check and draw left pi
820  constexpr double tp_x = 0.5 - radius - 0.02;
821  constexpr double tp_y = 0.5;
822  TLatex tp(tp_x, tp_y, "#pi");
823  if (distance(tbf_x, tbf_y, tp_x, tp_y) > too_close_threshold) {
824  draw_text(tp);
825  }
826 
827  // Check and draw top pi/2
828  constexpr double tp2_x = 0.5;
829  constexpr double tp2_y = 0.5 + radius + 0.04;
830  TLatex tp2(tp2_x, tp2_y, "#frac{#pi}{2}");
831  if (distance(tbf_x, tbf_y, tp2_x, tp2_y) > too_close_threshold) {
832  draw_text(tp2);
833  }
834 
835  // Check and draw bottom -pi/2
836  constexpr double tmp2_x = 0.5;
837  constexpr double tmp2_y = 0.5 - radius - 0.04;
838  TLatex tmp2(tmp2_x, tmp2_y, "-#frac{#pi}{2}");
839  if (distance(tbf_x, tbf_y, tmp2_x, tmp2_y) > too_close_threshold) {
840  draw_text(tmp2);
841  }
842 
843  TLatex tbf(0.5 + (radius + 0.02) * cos(best_fit),
844  0.5 + (radius + 0.02) * sin(best_fit),
845  fmt::format("{:.2f}", best_fit).c_str());
846  draw_text(tbf, kRed);
847 
848  canvas.Print(CanvasName);
849 }
850 
851 // ***************
852 // MP: Calculate PMNS matrix elements
853 std::array<std::array<TComplex, 3>, 3> OscProcessor::CalculatePMNSElements(const double s2th13, const double s2th23, const double s2th12, const double dcp) const {
854 // ***************
855  const double s13 = std::sqrt(s2th13);
856  const double s23 = std::sqrt(s2th23);
857  const double s12 = std::sqrt(s2th12);
858 
859  const double sdcp = std::sin(dcp);
860  const double cdcp = std::cos(dcp);
861 
862  const double c13 = std::sqrt(1.-s2th13);
863  const double c12 = std::sqrt(1.-s2th12);
864  const double c23 = std::sqrt(1.-s2th23);
865 
866  double real_ue[3];
867  double imag_ue[3];
868  real_ue[0] = c12*c13;
869  imag_ue[0] = 0.;
870  real_ue[1] = s12*c13;
871  imag_ue[1] = 0.;
872  real_ue[2] = s13*cdcp;
873  imag_ue[2] = -s13*sdcp;
874 
875  double real_umu[3];
876  double imag_umu[3];
877  real_umu[0] = -s12*c23 - c12*s23*s13*cdcp;
878  imag_umu[0] = -c12*s23*s13*sdcp;
879  real_umu[1] = c12*c23 - s12*s23*s13*cdcp;
880  imag_umu[1] = -s12*s23*s13*sdcp;
881  real_umu[2] = s23*c13;
882  imag_umu[2] = 0.;
883 
884  double real_utau[3];
885  double imag_utau[3];
886  real_utau[0] = s12*s23 - c12*c23*s13*cdcp;
887  imag_utau[0] = -c12*c23*s13*sdcp;
888  real_utau[1] = -c12*s23 - s12*c23*s13*cdcp;
889  imag_utau[1] = -s12*c23*s13*sdcp;
890  real_utau[2] = c23*c13;
891  imag_utau[2] = 0.;
892 
893  TComplex U[3][3];
894  for (int i = 0; i < 3; i++) {
895  U[0][i] = TComplex(real_ue[i], imag_ue[i]);
896  U[1][i] = TComplex(real_umu[i], imag_umu[i]);
897  U[2][i] = TComplex(real_utau[i], imag_utau[i]);
898  }
899 
900  return {{{U[0][0], U[0][1], U[0][2]}, {U[1][0], U[1][1], U[1][2]}, {U[2][0], U[2][1], U[2][2]}}};
901 }
902 
903 // ***************
904 // MP: Produce PMNS matrix elements
906 // ***************
907  if(!OscEnabled ||
912  {
913  MACH3LOG_WARN("Will not {}, as oscillation parameters are missing", __func__);
914  return;
915  }
916  MACH3LOG_INFO("Starting {}", __func__);
917 
918  double s2th13, s2th23, s2th12, dcp = M3::_BAD_DOUBLE_;
919  double weight = 1.0;
920 
921  TDirectory *PMNSElementsDir = OutputFile->mkdir("PMNSElements");
922  PMNSElementsDir->cd();
923 
924  unsigned int step = 0;
925  Chain->SetBranchStatus("*", false);
926 
927  Chain->SetBranchStatus(Sin2Theta13Name.c_str(), true);
928  Chain->SetBranchAddress(Sin2Theta13Name.c_str(), &s2th13);
929 
930  Chain->SetBranchStatus(Sin2Theta23Name.c_str(), true);
931  Chain->SetBranchAddress(Sin2Theta23Name.c_str(), &s2th23);
932 
933  Chain->SetBranchStatus(Sin2Theta12Name.c_str(), true);
934  Chain->SetBranchAddress(Sin2Theta12Name.c_str(), &s2th12);
935 
936  Chain->SetBranchStatus(DeltaCPName.c_str(), true);
937  Chain->SetBranchAddress(DeltaCPName.c_str(), &dcp);
938 
939  Chain->SetBranchStatus("step", true);
940  Chain->SetBranchAddress("step", &step);
941 
942  if (Chain->GetBranch("Weight")) {
943  Chain->SetBranchStatus("Weight", true);
944  Chain->SetBranchAddress("Weight", &weight);
945  } else {
946  MACH3LOG_WARN("No Weight branch found — defaulting to 1.0");
947  weight = 1.0;
948  }
949  constexpr int n_bins = 1000;
950 
951  std::unique_ptr<TH1D> h_ue[3];
952  std::unique_ptr<TH1D> h_umu[3];
953  std::unique_ptr<TH1D> h_utau[3];
954 
955  std::unique_ptr<TH1D> h_ue_real[3];
956  std::unique_ptr<TH1D> h_umu_real[3];
957  std::unique_ptr<TH1D> h_utau_real[3];
958 
959  std::unique_ptr<TH1D> h_ue_imag[3];
960  std::unique_ptr<TH1D> h_umu_imag[3];
961  std::unique_ptr<TH1D> h_utau_imag[3];
962 
963  std::unique_ptr<TH2D> h_ue_s2th12[3];
964  std::unique_ptr<TH2D> h_umu_s2th12[3];
965  std::unique_ptr<TH2D> h_utau_s2th12[3];
966 
967  std::unique_ptr<TH2D> h_ue_s2th13[3];
968  std::unique_ptr<TH2D> h_umu_s2th13[3];
969  std::unique_ptr<TH2D> h_utau_s2th13[3];
970 
971  std::unique_ptr<TH2D> h_ue_s2th23[3];
972  std::unique_ptr<TH2D> h_umu_s2th23[3];
973  std::unique_ptr<TH2D> h_utau_s2th23[3];
974 
975  std::unique_ptr<TH2D> h_ue_dcp[3];
976  std::unique_ptr<TH2D> h_umu_dcp[3];
977  std::unique_ptr<TH2D> h_utau_dcp[3];
978 
979  for (int iU = 0; iU < 3; iU++) {
980  h_ue[iU] = std::make_unique<TH1D>(Form("h_ue%d", iU+1), Form(";|U_{e%d}|", iU+1), n_bins, 0., 1.);
981  h_ue[iU]->SetDirectory(nullptr);
982  h_umu[iU] = std::make_unique<TH1D>(Form("h_umu%d", iU+1), Form(";|U_{#mu%d}|", iU+1), n_bins, 0., 1.);
983  h_umu[iU]->SetDirectory(nullptr);
984  h_utau[iU] = std::make_unique<TH1D>(Form("h_utau%d", iU+1), Form(";|U_{#tau%d}|", iU+1), n_bins, 0., 1.);
985  h_utau[iU]->SetDirectory(nullptr);
986 
987  h_ue_real[iU] = std::make_unique<TH1D>(Form("h_ue%d_real", iU+1), Form(";Re(U_{e%d})", iU+1), n_bins, -1., 1.);
988  h_ue_real[iU]->SetDirectory(nullptr);
989  h_umu_real[iU] = std::make_unique<TH1D>(Form("h_umu%d_real", iU+1), Form(";Re(U_{#mu%d})", iU+1), n_bins, -1., 1.);
990  h_umu_real[iU]->SetDirectory(nullptr);
991  h_utau_real[iU] = std::make_unique<TH1D>(Form("h_utau%d_real", iU+1), Form(";Re(U_{#tau%d})", iU+1), n_bins, -1., 1.);
992  h_utau_real[iU]->SetDirectory(nullptr);
993 
994  h_ue_imag[iU] = std::make_unique<TH1D>(Form("h_ue%d_imag", iU+1), Form(";Im(U_{e%d})", iU+1), n_bins, -1., 1.);
995  h_ue_imag[iU]->SetDirectory(nullptr);
996  h_umu_imag[iU] = std::make_unique<TH1D>(Form("h_umu%d_imag", iU+1), Form(";Im(U_{#mu%d})", iU+1), n_bins, -1., 1.);
997  h_umu_imag[iU]->SetDirectory(nullptr);
998  h_utau_imag[iU] = std::make_unique<TH1D>(Form("h_utau%d_imag", iU+1), Form(";Im(U_{#tau%d})", iU+1), n_bins, -1., 1.);
999  h_utau_imag[iU]->SetDirectory(nullptr);
1000 
1001  h_ue_s2th12[iU] = std::make_unique<TH2D>(Form("h_ue%d_s2th12", iU+1), Form(";|U_{e%d}|;sin^{2}(#theta_{12})", iU+1),
1002  n_bins, 0., 1., n_bins, 0.15, 0.5);
1003  h_ue_s2th12[iU]->SetDirectory(nullptr);
1004  h_umu_s2th12[iU] = std::make_unique<TH2D>(Form("h_umu%d_s2th12", iU+1), Form(";|U_{#mu%d}|;sin^{2}(#theta_{12})", iU+1),
1005  n_bins, 0., 1., n_bins, 0.15, 0.5);
1006  h_umu_s2th12[iU]->SetDirectory(nullptr);
1007  h_utau_s2th12[iU] = std::make_unique<TH2D>(Form("h_utau%d_s2th12", iU+1), Form(";|U_{#tau%d}|;sin^{2}(#theta_{12})", iU+1),
1008  n_bins, 0., 1., n_bins, 0.15, 0.5);
1009  h_utau_s2th12[iU]->SetDirectory(nullptr);
1010 
1011 
1012  h_ue_s2th13[iU] = std::make_unique<TH2D>(Form("h_ue%d_s2th13", iU+1), Form(";|U_{e%d}|;sin^{2}(#theta_{13})", iU+1),
1013  n_bins, 0., 1., n_bins, 0., 0.1);
1014  h_ue_s2th13[iU]->SetDirectory(nullptr);
1015  h_umu_s2th13[iU] = std::make_unique<TH2D>(Form("h_umu%d_s2th13", iU+1), Form(";|U_{#mu%d}|;sin^{2}(#theta_{13})", iU+1),
1016  n_bins, 0., 1., n_bins, 0., 0.1);
1017  h_umu_s2th13[iU]->SetDirectory(nullptr);
1018  h_utau_s2th13[iU] = std::make_unique<TH2D>(Form("h_utau%d_s2th13", iU+1), Form(";|U_{#tau%d}|;sin^{2}(#theta_{13})", iU+1),
1019  n_bins, 0., 1., n_bins, 0., 0.1);
1020  h_utau_s2th13[iU]->SetDirectory(nullptr);
1021 
1022 
1023  h_ue_s2th23[iU] = std::make_unique<TH2D>(Form("h_ue%d_s2th23", iU+1), Form(";|U_{e%d}|;sin^{2}(#theta_{23})", iU+1),
1024  n_bins, 0., 1., n_bins, 0.3, 0.8);
1025  h_ue_s2th23[iU]->SetDirectory(nullptr);
1026  h_umu_s2th23[iU] = std::make_unique<TH2D>(Form("h_umu%d_s2th23", iU+1), Form(";|U_{#mu%d}|;sin^{2}(#theta_{23})", iU+1),
1027  n_bins, 0., 1., n_bins, 0.3, 0.8);
1028  h_umu_s2th23[iU]->SetDirectory(nullptr);
1029  h_utau_s2th23[iU] = std::make_unique<TH2D>(Form("h_utau%d_s2th23", iU+1), Form(";|U_{#tau%d}|;sin^{2}(#theta_{23})", iU+1),
1030  n_bins, 0., 1., n_bins, 0.3, 0.8);
1031  h_utau_s2th23[iU]->SetDirectory(nullptr);
1032 
1033 
1034  h_ue_dcp[iU] = std::make_unique<TH2D>(Form("h_ue%d_dcp", iU+1), Form(";|U_{e%d}|;#delta_{CP}", iU+1),
1035  n_bins, 0., 1., n_bins, -TMath::Pi(), TMath::Pi());
1036  h_ue_dcp[iU]->SetDirectory(nullptr);
1037  h_umu_dcp[iU] = std::make_unique<TH2D>(Form("h_umu%d_dcp", iU+1), Form(";|U_{#mu%d}|;#delta_{CP}", iU+1),
1038  n_bins, 0., 1., n_bins, -TMath::Pi(), TMath::Pi());
1039  h_umu_dcp[iU]->SetDirectory(nullptr);
1040  h_utau_dcp[iU] = std::make_unique<TH2D>(Form("h_utau%d_dcp", iU+1), Form(";|U_{#tau%d}|;#delta_{CP}", iU+1),
1041  n_bins, 0., 1., n_bins, -TMath::Pi(), TMath::Pi());
1042  h_utau_dcp[iU]->SetDirectory(nullptr);
1043  }
1044 
1045  // MP: 2D histograms for all unique pairs of |U_{alpha i}| vs |U_{beta j}| (no self-pairs)
1046  std::vector<std::unique_ptr<TH2D>> h_UU;
1047  std::string U_names[9] = {"ue1", "ue2", "ue3", "umu1", "umu2", "umu3", "utau1", "utau2", "utau3"};
1048  std::string U_tex[9] = {"|U_{e1}|", "|U_{e2}|", "|U_{e3}|", "|U_{#mu1}|", "|U_{#mu2}|", "|U_{#mu3}|", "|U_{#tau1}|", "|U_{#tau2}|", "|U_{#tau3}|"};
1049  for(int i=0; i < 9; ++i){
1050  for(int j = i+1; j < 9; ++j){
1051  std::string name = "h_" + U_names[i] + "_" + U_names[j];
1052  std::string title = ";" + U_tex[i] + ";" + U_tex[j];
1053  h_UU.push_back(std::make_unique<TH2D>(name.c_str(), title.c_str(), n_bins, 0., 1., n_bins, 0., 1.));
1054  h_UU.back()->SetDirectory(nullptr);
1055  }
1056  }
1057  const Long64_t countwidth = nEntries/5;
1058  for(int i = 0; i < nEntries; i++) {
1059  if (i % countwidth == 0) {
1062  } else {
1063  Chain->GetEntry(i);
1064  }
1065 
1066  if(step < BurnInCut) continue; // burn-in cut
1067 
1068  const std::array<std::array<TComplex, 3>, 3> U = CalculatePMNSElements(s2th13, s2th23, s2th12, dcp);
1069 
1070  for(int iU = 0; iU < 3; iU++){
1071  h_ue[iU]->Fill(TComplex::Abs(U[0][iU]), weight);
1072  h_umu[iU]->Fill(TComplex::Abs(U[1][iU]), weight);
1073  h_utau[iU]->Fill(TComplex::Abs(U[2][iU]), weight);
1074 
1075  h_ue_real[iU]->Fill(U[0][iU].Re(), weight);
1076  h_umu_real[iU]->Fill(U[1][iU].Re(), weight);
1077  h_utau_real[iU]->Fill(U[2][iU].Re(), weight);
1078 
1079  h_ue_imag[iU]->Fill(U[0][iU].Im(), weight);
1080  h_umu_imag[iU]->Fill(U[1][iU].Im(), weight);
1081  h_utau_imag[iU]->Fill(U[2][iU].Im(), weight);
1082 
1083  h_ue_s2th12[iU]->Fill(TComplex::Abs(U[0][iU]), s2th12, weight);
1084  h_umu_s2th12[iU]->Fill(TComplex::Abs(U[1][iU]), s2th12, weight);
1085  h_utau_s2th12[iU]->Fill(TComplex::Abs(U[2][iU]), s2th12, weight);
1086 
1087  h_ue_s2th13[iU]->Fill(TComplex::Abs(U[0][iU]), s2th13, weight);
1088  h_umu_s2th13[iU]->Fill(TComplex::Abs(U[1][iU]), s2th13, weight);
1089  h_utau_s2th13[iU]->Fill(TComplex::Abs(U[2][iU]), s2th13, weight);
1090 
1091  h_ue_s2th23[iU]->Fill(TComplex::Abs(U[0][iU]), s2th23, weight);
1092  h_umu_s2th23[iU]->Fill(TComplex::Abs(U[1][iU]), s2th23, weight);
1093  h_utau_s2th23[iU]->Fill(TComplex::Abs(U[2][iU]), s2th23, weight);
1094 
1095  h_ue_dcp[iU]->Fill(TComplex::Abs(U[0][iU]), dcp, weight);
1096  h_umu_dcp[iU]->Fill(TComplex::Abs(U[1][iU]), dcp, weight);
1097  h_utau_dcp[iU]->Fill(TComplex::Abs(U[2][iU]), dcp, weight);
1098  }
1099 
1100  // MP: Store in an array for easy access
1101  double U_mod[9] = {TComplex::Abs(U[0][0]), TComplex::Abs(U[0][1]), TComplex::Abs(U[0][2]), TComplex::Abs(U[1][0]), TComplex::Abs(U[1][1]), TComplex::Abs(U[1][2]), TComplex::Abs(U[2][0]), TComplex::Abs(U[2][1]), TComplex::Abs(U[2][2])};
1102  int idx = 0;
1103  for(int ix = 0; ix < 9; ++ix){
1104  for(int iy = ix+1; iy<9; ++iy){
1105  h_UU[idx]->Fill(U_mod[ix], U_mod[iy], weight);
1106  ++idx;
1107  }
1108  }
1109  } // end loop over steps
1110 
1111  // Now we save
1112  PMNSElementsDir->cd();
1113 
1114  for (int iU = 0; iU < 3; iU++) {
1115  h_ue[iU]->Write(Form("h_ue%d", iU+1));
1116  h_umu[iU]->Write(Form("h_umu%d", iU+1));
1117  h_utau[iU]->Write(Form("h_utau%d", iU+1));
1118 
1119  h_ue_real[iU]->Write(Form("h_ue%d_real", iU+1));
1120  h_umu_real[iU]->Write(Form("h_umu%d_real", iU+1));
1121  h_utau_real[iU]->Write(Form("h_utau%d_real", iU+1));
1122 
1123  h_ue_imag[iU]->Write(Form("h_ue%d_imag", iU+1));
1124  h_umu_imag[iU]->Write(Form("h_umu%d_imag", iU+1));
1125  h_utau_imag[iU]->Write(Form("h_utau%d_imag", iU+1));
1126 
1127  h_ue_s2th12[iU]->Write(Form("h_ue%d_s2th12", iU+1));
1128  h_umu_s2th12[iU]->Write(Form("h_umu%d_s2th12", iU+1));
1129  h_utau_s2th12[iU]->Write(Form("h_utau%d_s2th12", iU+1));
1130 
1131  h_ue_s2th13[iU]->Write(Form("h_ue%d_s2th13", iU+1));
1132  h_umu_s2th13[iU]->Write(Form("h_umu%d_s2th13", iU+1));
1133  h_utau_s2th13[iU]->Write(Form("h_utau%d_s2th13", iU+1));
1134 
1135  h_ue_s2th23[iU]->Write(Form("h_ue%d_s2th23", iU+1));
1136  h_umu_s2th23[iU]->Write(Form("h_umu%d_s2th23", iU+1));
1137  h_utau_s2th23[iU]->Write(Form("h_utau%d_s2th23", iU+1));
1138 
1139  h_ue_dcp[iU]->Write(Form("h_ue%d_dcp", iU+1));
1140  h_umu_dcp[iU]->Write(Form("h_umu%d_dcp", iU+1));
1141  h_utau_dcp[iU]->Write(Form("h_utau%d_dcp", iU+1));
1142  }
1143 
1144  // MP: Write all |U_{alpha i}| vs |U_{beta j}| 2D histograms (no self-pairs)
1145  for(size_t iPlot = 0; iPlot < h_UU.size(); ++iPlot){
1146  h_UU[iPlot]->Write();
1147  }
1148 
1149  PMNSElementsDir->Close();
1150  delete PMNSElementsDir;
1151 
1152  Chain->SetBranchStatus("*", true);
1153  OutputFile->cd();
1154 }
1155 
1156 // ***************
1157 // MP: Produce unitarity triangles from PMNS matrix elements
1159 // ***************
1160  if(!OscEnabled ||
1165  {
1166  MACH3LOG_WARN("Will not {}, as oscillation parameters are missing", __func__);
1167  return;
1168  }
1169  MACH3LOG_INFO("Starting {}", __func__);
1170 
1171  double s2th13, s2th23, s2th12, dcp = M3::_BAD_DOUBLE_;
1172  double weight = 1.0;
1173 
1174  TComplex tr_emu_num, tr_etau_num, tr_mutau_num, tr_12_num, tr_13_num, tr_23_num;
1175  TComplex tr_emu_denom, tr_etau_denom, tr_mutau_denom, tr_12_denom, tr_13_denom, tr_23_denom;
1176  TComplex tr_emu, tr_etau, tr_mutau, tr_12, tr_13, tr_23;
1177 
1178  TDirectory *UnitarityTrianglesDir = OutputFile->mkdir("UnitarityTriangles");
1179  UnitarityTrianglesDir->cd();
1180 
1181  unsigned int step = 0;
1182  Chain->SetBranchStatus("*", false);
1183 
1184  Chain->SetBranchStatus(Sin2Theta13Name.c_str(), true);
1185  Chain->SetBranchAddress(Sin2Theta13Name.c_str(), &s2th13);
1186 
1187  Chain->SetBranchStatus(Sin2Theta23Name.c_str(), true);
1188  Chain->SetBranchAddress(Sin2Theta23Name.c_str(), &s2th23);
1189 
1190  Chain->SetBranchStatus(Sin2Theta12Name.c_str(), true);
1191  Chain->SetBranchAddress(Sin2Theta12Name.c_str(), &s2th12);
1192 
1193  Chain->SetBranchStatus(DeltaCPName.c_str(), true);
1194  Chain->SetBranchAddress(DeltaCPName.c_str(), &dcp);
1195 
1196  Chain->SetBranchStatus("step", true);
1197  Chain->SetBranchAddress("step", &step);
1198 
1199  if (Chain->GetBranch("Weight")) {
1200  Chain->SetBranchStatus("Weight", true);
1201  Chain->SetBranchAddress("Weight", &weight);
1202  } else {
1203  MACH3LOG_WARN("No Weight branch found — defaulting to 1.0");
1204  weight = 1.0;
1205  }
1206  constexpr int n_bins = 1000;
1207 
1208  std::unique_ptr<TH2D> h_tr_emu;
1209  std::unique_ptr<TH2D> h_tr_etau;
1210  std::unique_ptr<TH2D> h_tr_mutau;
1211 
1212  std::unique_ptr<TH2D> h_tr_12;
1213  std::unique_ptr<TH2D> h_tr_13;
1214  std::unique_ptr<TH2D> h_tr_23;
1215 
1216  h_tr_emu = std::make_unique<TH2D>("h_tr_emu",";#rho_{e#mu};#eta_{e#mu}",n_bins,-5,5,n_bins,-5,5);
1217  h_tr_emu->SetDirectory(nullptr);
1218  h_tr_etau = std::make_unique<TH2D>("h_tr_etau",";#rho_{e#tau};#eta_{e#tau}",n_bins,0,2,n_bins,-1,1);
1219  h_tr_etau->SetDirectory(nullptr);
1220  h_tr_mutau = std::make_unique<TH2D>("h_tr_mutau",";#rho_{#mu#tau};#eta_{#mu#tau}",n_bins,0,2,n_bins,-1,1);
1221  h_tr_mutau->SetDirectory(nullptr);
1222 
1223  h_tr_12 = std::make_unique<TH2D>("h_tr_12",";#rho_{12};#eta_{12}",n_bins,0,3.2,n_bins,-1.5,1.5);
1224  h_tr_12->SetDirectory(nullptr);
1225  h_tr_13 = std::make_unique<TH2D>("h_tr_13",";#rho_{13};#eta_{13}",n_bins,0,2,n_bins,-1,1);
1226  h_tr_13->SetDirectory(nullptr);
1227  h_tr_23 = std::make_unique<TH2D>("h_tr_23",";#rho_{23};#eta_{23}",n_bins,-8,8,n_bins,-8,8);
1228  h_tr_23->SetDirectory(nullptr);
1229 
1230  const Long64_t countwidth = nEntries/5;
1231  for(int i = 0; i < nEntries; i++) {
1232  if (i % countwidth == 0) {
1235  } else {
1236  Chain->GetEntry(i);
1237  }
1238 
1239  if(step < BurnInCut) continue; // burn-in cut
1240 
1241  const std::array<std::array<TComplex, 3>, 3> U = CalculatePMNSElements(s2th13, s2th23, s2th12, dcp);
1242 
1243  // Calculate the sides of the unitarity triangles
1244  tr_emu_num = U[0][0].operator*(TComplex::Conjugate(U[1][0]));
1245  tr_emu_denom = U[0][2].operator*(TComplex::Conjugate(U[1][2]));
1246  tr_emu = - tr_emu_num.operator/(tr_emu_denom);
1247 
1248  tr_etau_num = U[0][1].operator*(TComplex::Conjugate(U[2][1]));
1249  tr_etau_denom = U[0][0].operator*(TComplex::Conjugate(U[2][0]));
1250  tr_etau = - tr_etau_num.operator/(tr_etau_denom);
1251 
1252  tr_mutau_num = U[1][2].operator*(TComplex::Conjugate(U[2][2]));
1253  tr_mutau_denom = U[1][1].operator*(TComplex::Conjugate(U[2][1]));
1254  tr_mutau = - tr_mutau_num.operator/(tr_mutau_denom);
1255 
1256  tr_12_num = U[0][0].operator*(TComplex::Conjugate(U[0][1]));
1257  tr_12_denom = U[1][0].operator*(TComplex::Conjugate(U[1][1]));
1258  tr_12 = - tr_12_num.operator/(tr_12_denom);
1259 
1260  tr_13_num = U[1][0].operator*(TComplex::Conjugate(U[1][2]));
1261  tr_13_denom = U[2][0].operator*(TComplex::Conjugate(U[2][2]));
1262  tr_13 = - tr_13_num.operator/(tr_13_denom);
1263 
1264  tr_23_num = U[2][1].operator*(TComplex::Conjugate(U[2][2]));
1265  tr_23_denom = U[0][1].operator*(TComplex::Conjugate(U[0][2]));
1266  tr_23 = - tr_23_num.operator/(tr_23_denom);
1267 
1268  h_tr_emu->Fill(tr_emu.Re(), tr_emu.Im(), weight);
1269  h_tr_etau->Fill(tr_etau.Re(), tr_etau.Im(), weight);
1270  h_tr_mutau->Fill(tr_mutau.Re(), tr_mutau.Im(), weight);
1271 
1272  h_tr_12->Fill(tr_12.Re(), tr_12.Im(), weight);
1273  h_tr_13->Fill(tr_13.Re(), tr_13.Im(), weight);
1274  h_tr_23->Fill(tr_23.Re(), tr_23.Im(), weight);
1275  } // end loop over steps
1276 
1277  // Now we save
1278  UnitarityTrianglesDir->cd();
1279 
1280  h_tr_emu->Write("h_tr_emu");
1281  h_tr_etau->Write("h_tr_etau");
1282  h_tr_mutau->Write("h_tr_mutau");
1283 
1284  h_tr_12->Write("h_tr_12");
1285  h_tr_13->Write("h_tr_13");
1286  h_tr_23->Write("h_tr_23");
1287 
1288  UnitarityTrianglesDir->Close();
1289  delete UnitarityTrianglesDir;
1290 
1291  Chain->SetBranchStatus("*", true);
1292  OutputFile->cd();
1293 }
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
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:61
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.
void ProduceUnitarityTriangles()
MP: Produce unitarity triangles from PMNS matrix elements.
int DeltaCPIndex
Index of in the parameter list.
Definition: OscProcessor.h:101
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:84
std::string Sin2Theta12Name
Name of the parameter representing .
Definition: OscProcessor.h:86
int Sin2Theta12Index
Index of in the parameter list.
Definition: OscProcessor.h:97
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:103
std::string DeltaCPName
Name of the parameter representing (the CP-violating phase).
Definition: OscProcessor.h:90
bool PlotJarlskog
Will plot Jarlskog Invariant using information in the chain.
Definition: OscProcessor.h:78
void ProducePMNSElements()
MP: Produce PMNS matrix elements.
int Sin2Theta13Index
Index of in the parameter list.
Definition: OscProcessor.h:95
std::string Sin2Theta23Name
Name of the parameter representing .
Definition: OscProcessor.h:88
std::array< std::array< TComplex, 3 >, 3 > CalculatePMNSElements(const double s2th13, const double s2th23, const double s2th12, const double dcp) const
MP: Calculate PMNS matrix elements.
void Get1DReactorConstraintInfo(std::pair< double, double > &Sin13_NewPrior, bool &DoReweight) const
Extract 1D reactor constraint information from an MCMC file.
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:99
virtual ~OscProcessor()
Destroys the OscProcessor object.
std::string DeltaM2_23Name
Name of the parameter representing (mass-squared difference).
Definition: OscProcessor.h:92
void LoadAdditionalInfo() final
Read the Osc cov file and get the input central values and errors Here we allow Jarlskog Shenanigans.
bool OscEnabled
Will plot Jarlskog Invariant using information in the chain.
Definition: OscProcessor.h:81
void MakePiePlot()
Make fancy Pie plot for delta CP.
void PerformJarlskogAnalysis()
Perform Several Jarlskog Plotting.
void EstimateDataTransferRate(TChain *chain, const Long64_t entry)
KS: Check what CPU you are using.
Definition: Monitor.cpp:212
void PrintProgressBar(const Long64_t Done, const Long64_t All)
KS: Simply print progress bar.
Definition: Monitor.cpp:229
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