11 #pragma GCC diagnostic ignored "-Wfloat-conversion"
47 for (
int i = 0; i <
nDraw; ++i)
52 const std::string CurrentName =
ParamNames[ParamEnum][ParamNo].Data();
55 if (CurrentName ==
"sin2th_13") {
58 }
else if (CurrentName ==
"sin2th_12") {
61 }
else if (CurrentName ==
"sin2th_23") {
64 }
else if (CurrentName ==
"delta_cp") {
67 }
else if (CurrentName ==
"delm2_23") {
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)");
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);
108 const double j = s13*c13*c13*s12*c12*s23*c23*sdcp;
117 double Prior = 1.0, PriorError = 1.0;
118 bool FlatPrior =
false;
125 FlatPrior =
ParamFlat[ParType][ParamTemp];
128 return randGen->Uniform(FlatBounds[0], FlatBounds[1]);
131 return randGen->Gaus(Prior, PriorError);
146 MACH3LOG_WARN(
"Will not {}, as oscillation parameters are missing", __func__);
151 bool DoReweight =
false;
155 std::pair<double, double> Sin13_NewPrior;
158 TFile *TempFile =
M3::Open((
MCMCFile +
".root"),
"open", __FILE__, __LINE__);
161 TMacro *Config = TempFile->Get<TMacro>(
"Reweight_Config");
163 if (Config !=
nullptr) {
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__);
178 MACH3LOG_INFO(
"No valid reweighting configuration (1D Gaussian on sin2th_13 only) found for Jarlskog analysis");
181 MACH3LOG_INFO(
"No reweighting configuration found for Jarlskog analysis");
188 TDirectory *JarlskogDir =
OutputFile->mkdir(
"Jarlskog");
191 unsigned int step = 0;
192 Chain->SetBranchStatus(
"*",
false);
209 Chain->SetBranchStatus(
"step",
true);
210 Chain->SetBranchAddress(
"step", &step);
213 Chain->SetBranchStatus(
"Weight",
true);
214 Chain->SetBranchAddress(
"Weight", &weight);
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);
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");
232 auto jarl_IH =
M3::Clone(jarl.get(),
"jarl_IH");
233 auto jarl_NH =
M3::Clone(jarl.get(),
"jarl_NH");
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");
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");
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");
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");
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;
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");
260 auto prior3 = std::make_unique<TF1>(
"prior3",
"TMath::Abs(TMath::Cos(x))");
263 auto randGen = std::make_unique<TRandom3>(0);
264 const Long64_t countwidth =
nEntries/5;
267 if (i % countwidth == 0) {
276 const double j =
CalcJarlskog(s2th13, s2th23, s2th12, dcp);
277 const double prior_weight = prior3->Eval(dcp);
279 jarl->Fill(j, weight);
280 jarl_th23->Fill(j, s2th23, weight);
281 jarl_dcp->Fill(j, dcp, weight);
283 jarl_flatsindcp->Fill(j, prior_weight*weight);
284 jarl_th23_flatsindcp->Fill(j, s2th23, prior_weight*weight);
291 const double prior_sindcp = randGen->Uniform(-1., 1.);
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;
303 jarl_prior->Fill(prior_j);
304 jarl_prior_flatsindcp->Fill(prior_flatsindcp_j);
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);
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);
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);
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);
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");
343 jarl_dcp->Write(
"jarlskog_dcp_both");
344 jarl_dcp_NH->Write(
"jarlskog_dcp_NH");
345 jarl_dcp_IH->Write(
"jarlskog_dcp_IH");
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");
355 jarl_prior->Write(
"jarl_prior");
356 jarl_prior_flatsindcp->Write(
"jarl_prior_flatsindcp");
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");
364 jarl_NH, jarl_NH_flatsindcp,
365 jarl_IH, jarl_IH_flatsindcp);
370 SavageDickeyPlot(jarl_flatsindcp, jarl_wRC_prior_flatsindcp,
"Jarlskog flat sin#delta_{CP}", 0);
373 SavageDickeyPlot(jarl_flatsindcp, jarl_prior_flatsindcp,
"Jarlskog flat sin#delta_{CP}", 0);
376 JarlskogDir->Close();
379 Chain->SetBranchStatus(
"*",
true);
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) {
393 int originalErrorLevel = gErrorIgnoreLevel;
394 gErrorIgnoreLevel = kFatal;
397 for(
int hierarchy = -1; hierarchy <= 1; hierarchy++)
399 std::unique_ptr<TH1D> j_hist;
400 std::unique_ptr<TH1D> j_hist_sdcp;
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) {
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) {
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");
414 MACH3LOG_ERROR(
"Invalid hierarchy option. 1 for NH, 0 for both, -1 for IH");
419 j_hist_sdcp->Rebin(7);
421 j_hist->SetLineColor(kAzure-2);
422 j_hist_sdcp->SetLineColor(kOrange+1);
423 j_hist->SetLineWidth(2);
424 j_hist_sdcp->SetLineWidth(2);
426 auto StyleAxis = [](TH1* h) {
427 auto xAxis = h->GetXaxis();
428 auto yAxis = h->GetYaxis();
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);
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);
447 StyleAxis(j_hist.get());
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());
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");
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.;
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");
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.;
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.;
487 integral = j_hist_copy->Integral();
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.;
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.;
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.;
518 integral = j_hist_sdcp_copy->Integral();
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.;
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.;
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.;
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);
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);
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);
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);
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);
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);
587 double arrowLength = 0.003;
588 double arrowHeight = vertUp;
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);
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());
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());
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);
608 auto CopyLineStyle = [](
const TH1D* src, TLine* dst) {
609 dst->SetLineColor(src->GetLineColor());
610 dst->SetLineStyle(src->GetLineStyle());
611 dst->SetLineWidth(src->GetLineWidth());
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());
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());
628 auto leg = std::make_unique<TLegend>(0.45, 0.60, 0.75, 0.90);
629 leg->SetTextSize(0.05);
630 leg->SetFillStyle(0);
632 leg->SetTextFont(132);
633 leg->SetBorderSize(0);
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");
641 j_hist->GetYaxis()->SetRangeUser(0., j_hist->GetMaximum()*1.15);
643 j_hist_sdcp->Draw(
"same h");
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");
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();
660 auto ttext = std::make_unique<TText>();
662 ttext->SetTextSize(0.03);
663 ttext->SetTextAlign(13);
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");
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");
680 gErrorIgnoreLevel = originalErrorLevel;
688 MACH3LOG_WARN(
"Will not {}, as oscillation parameters are missing", __func__);
696 const double sigma_p = (*Errors_HPD_Positive)(
DeltaCPIndex);
697 const double sigma_n = (*Errors_HPD_Negative)(
DeltaCPIndex);
699 auto wrap_pi = [](
double x) {
700 while (x > TMath::Pi()) x -= 2*TMath::Pi();
701 while (x < -TMath::Pi()) x += 2*TMath::Pi();
705 std::array<double, 6> bds;
706 bds[0] = wrap_pi(best_fit - 3.0 * sigma_n);
707 bds[1] = wrap_pi(best_fit - 2.0 * sigma_n);
708 bds[2] = wrap_pi(best_fit - 1.0 * sigma_n);
709 bds[3] = wrap_pi(best_fit + 1.0 * sigma_p);
710 bds[4] = wrap_pi(best_fit + 2.0 * sigma_p);
711 bds[5] = wrap_pi(best_fit + 3.0 * sigma_p);
713 constexpr
double radius = 0.4;
714 constexpr
double rad_to_deg = 180.0 / TMath::Pi();
719 auto normalize_angle = [](
double rad) {
721 if (rad < 0) rad += 2.0 * TMath::Pi();
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);
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);
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);
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);
746 TLine line2(0.5, 0.5 - radius, 0.5, 0.5 + radius);
747 line2.SetLineWidth(3);
749 TArrow bf(0.5, 0.5, 0.5 + radius * cos(best_fit),0.5 + radius * sin(best_fit),0.04,
"|>");
751 bf.SetLineColor(kRed);
752 bf.SetFillColor(kRed);
754 TCanvas canvas(
"canvas",
"canvas", 0, 0, 1000, 1000);
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");
783 auto draw_text = [](
auto& txt, Color_t color = kBlack) {
784 txt.SetTextAlign(22);
785 txt.SetTextColor(color);
794 constexpr
double too_close_threshold = 0.1;
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);
801 auto distance = [](
double x1,
double y1,
double x2,
double y2) {
802 return std::sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2));
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) {
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) {
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) {
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) {
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);
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
YAML::Node TMacroToYAML(const TMacro ¯o)
KS: Convert a ROOT TMacro object to a YAML node.
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Class responsible for processing MCMC chains, performing diagnostics, generating plots,...
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.
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 .
std::string Sin2Theta12Name
Name of the parameter representing .
int Sin2Theta12Index
Index of in the parameter list.
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.
std::string DeltaCPName
Name of the parameter representing (the CP-violating phase).
bool PlotJarlskog
Will plot Jarlskog Invariant using information in the chain.
int Sin2Theta13Index
Index of in the parameter list.
std::string Sin2Theta23Name
Name of the parameter representing .
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.
virtual ~OscProcessor()
Destroys the OscProcessor object.
std::string DeltaM2_23Name
Name of the parameter representing (mass-squared difference).
bool OscEnabled
Will plot Jarlskog Invariant using information in the chain.
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.
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.
void PrintProgressBar(const Long64_t Done, const Long64_t All)
KS: Simply print progress bar.
void EstimateDataTransferRate(TChain *chain, const Long64_t entry)
KS: Check what CPU you are using.