MaCh3  2.4.2
Reference Guide
Functions
PlotMCMCDiag.cpp File Reference

KS: This script is used to analyse output form DiagMCMC. More...

#include "Fitters/MCMCProcessor.h"
#include "Samples/HistogramUtils.h"
#include "THStack.h"
#include "TGraphAsymmErrors.h"
#include "TLegend.h"
#include <TSystem.h>
Include dependency graph for PlotMCMCDiag.cpp:

Go to the source code of this file.

Functions

double GetMinimumInRange (TH1D *hist, const double minRange, const double maxRange)
 KS: function which looks for minimum in given range. More...
 
bool IsHistogramAllOnes (TH1D *hist, double tolerance=0.001, int max_failures=100)
 HW: Check if histogram is flat within a given tolerance. More...
 
void MakeDiagPlot (std::vector< TString > input_files, std::vector< TString > hist_labels)
 function which loops over Diag MCMC output and prints everything to PDF More...
 
void PlotAutoCorr (const std::vector< TString > &fname)
 Plot autocorrelation for all params into single plot with color codding helping whether on average it is ok or not. More...
 
std::pair< std::unique_ptr< TGraph >, std::unique_ptr< TGraph > > CreateMinMaxBand (TH1D *hist, Color_t color)
 HW: Create a band of minimum and maximum values from a histogram. More...
 
std::unique_ptr< TGraphAsymmErrors > CalculateMinMaxBand (const std::vector< std::unique_ptr< TH1D >> &histograms, Color_t color)
 
void ProcessAutoCorrelationDirectory (TDirectoryFile *autocor_dir, std::unique_ptr< TH1D > &average_hist, int &parameter_count, std::vector< std::unique_ptr< TH1D >> &histograms)
 Loop over directory get histograms into vector and add to averaged hist. More...
 
void ProcessDiagnosticFile (const TString &file_path, std::unique_ptr< TH1D > &average_hist, int &parameter_count, std::vector< std::unique_ptr< TH1D >> &histograms)
 
std::unique_ptr< TH1D > AutocorrProcessInputs (const TString &input_file, std::vector< std::unique_ptr< TH1D >> &histograms)
 
void CompareAverageAC (const std::vector< std::vector< std::unique_ptr< TH1D >>> &histograms, const std::vector< std::unique_ptr< TH1D >> &averages, const std::vector< TString > &hist_labels, const TString &output_name, bool draw_min_max=true, bool draw_all=false, bool draw_errors=true)
 
void PlotAverageACMult (std::vector< TString > input_files, std::vector< TString > hist_labels, const TString &output_name, bool draw_min_max=true, bool draw_all=false, bool draw_errors=true)
 Calculate mean AC based on all parameters with error band. Great for comparing AC between different chains. More...
 
int main (int argc, char *argv[])
 

Detailed Description

KS: This script is used to analyse output form DiagMCMC.

Warning
This script support comparing up to 4 files, there is easy way to expand it up to five or six,
Todo:
this need serious refactor
Author
Henry Wallace
Kamil Skwarczynski

Definition in file PlotMCMCDiag.cpp.

Function Documentation

◆ AutocorrProcessInputs()

std::unique_ptr<TH1D> AutocorrProcessInputs ( const TString &  input_file,
std::vector< std::unique_ptr< TH1D >> &  histograms 
)

Definition at line 356 of file PlotMCMCDiag.cpp.

357 {
358  std::unique_ptr<TH1D> average_hist = nullptr;
359  int parameter_count = 0;
360 
361  try {
362  ProcessDiagnosticFile(input_file, average_hist, parameter_count, histograms);
363  } catch (const std::exception &e) {
364  MACH3LOG_ERROR("Error processing file {} : {}", input_file, e.what());
365  }
366 
367  if (average_hist && parameter_count > 0)
368  {
369  MACH3LOG_INFO("Processed {} parameters from {} files", parameter_count, input_file);
370  average_hist->Scale(1.0 / parameter_count);
371  }
372 
373  return average_hist;
374 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
void ProcessDiagnosticFile(const TString &file_path, std::unique_ptr< TH1D > &average_hist, int &parameter_count, std::vector< std::unique_ptr< TH1D >> &histograms)

◆ CalculateMinMaxBand()

std::unique_ptr<TGraphAsymmErrors> CalculateMinMaxBand ( const std::vector< std::unique_ptr< TH1D >> &  histograms,
Color_t  color 
)

Definition at line 267 of file PlotMCMCDiag.cpp.

268 {
269  if (histograms.empty()) {
270  throw MaCh3Exception(__FILE__, __LINE__, "Empty histogram vector provided");
271  }
272 
273  int nBins = histograms[0]->GetNbinsX();
274  std::vector<double> x(nBins), ymin(nBins, 1e10), ymax(nBins, -1e10);
275 
276  for (int bin = 0; bin < nBins; bin++)
277  {
278  x[bin] = histograms[0]->GetBinCenter(bin + 1);
279 
280  for (const auto &hist : histograms)
281  {
282  double content = hist->GetBinContent(bin + 1);
283  ymin[bin] = std::min(ymin[bin], content);
284  ymax[bin] = std::max(ymax[bin], content);
285  }
286  }
287 
288  // Create a band using TGraphAsymmErrors for the shaded region
289  auto band = std::make_unique<TGraphAsymmErrors>(nBins);
290  for (int i = 0; i < nBins; ++i)
291  {
292  band->SetPoint(i, x[i], (ymax[i] + ymin[i]) / 2.0);
293  band->SetPointError(i, 0, 0, (ymax[i] - ymin[i]) / 2.0, (ymax[i] - ymin[i]) / 2.0);
294  }
295 
296  band->SetFillColorAlpha(color, 0.2f);
297  band->SetLineWidth(1);
298 
299  return band; // Return band and min graph (min graph can be used for legend)
300 }
Custom exception class used throughout MaCh3.

◆ CompareAverageAC()

void CompareAverageAC ( const std::vector< std::vector< std::unique_ptr< TH1D >>> &  histograms,
const std::vector< std::unique_ptr< TH1D >> &  averages,
const std::vector< TString > &  hist_labels,
const TString &  output_name,
bool  draw_min_max = true,
bool  draw_all = false,
bool  draw_errors = true 
)

Definition at line 376 of file PlotMCMCDiag.cpp.

383 {
384  auto canvas = std::make_unique<TCanvas>("AverageAC", "Average Auto Correlation", 800, 600);
385  canvas->SetGrid();
386 
387  // Setup colour palette
388  constexpr Int_t nb = 255.0;
389  Double_t stops[8] = {0.0, 0.25, 0.5, 0.75};
390  Double_t red[8] = {0.83203125, 0.796875, 0.0, 0.9375};
391  Double_t green[8] = {0.3671875, 0.47265625, 0.4453125, 0.890625};
392  Double_t blue[8] = {0.0, 0.65234375, 0.6953125, 0.2578125};
393  TColor::CreateGradientColorTable(8, stops, red, green, blue, nb);
394 
395  auto leg = std::make_unique<TLegend>(0.5, 0.7, 0.9, 0.9);
396  leg->SetFillColorAlpha(kWhite, 0.7f);
397  if (draw_min_max)
398  {
399  for (size_t i = 0; i < histograms.size(); ++i)
400  {
401  auto colour = static_cast<Color_t>(TColor::GetColorPalette(static_cast<Int_t>(i * nb / histograms.size())));
402  auto band = CalculateMinMaxBand(histograms[i], colour);
403  band->SetLineWidth(0);
404  band->SetTitle("Average Auto Correlation");
405  band->GetXaxis()->SetTitle("lag");
406  band->GetYaxis()->SetTitle("Autocorrelation Function");
407  if (i == 0) {
408  band->Draw("A3");
409  } else {
410  band->Draw("3 SAME");
411  }
412  }
413  }
414  std::vector<std::unique_ptr<TH1D>> error_hist(averages.size());
415  for (size_t i = 0; i < averages.size(); ++i)
416  {
417  auto colour = static_cast<Color_t>(TColor::GetColorPalette(static_cast<Int_t>(i * nb / histograms.size())));
418 
419  if (draw_errors)
420  {
421  error_hist[i] = M3::Clone(averages[i].get());
422  error_hist[i]->SetFillColorAlpha(colour, 0.3f);
423  error_hist[i]->SetLineWidth(0);
424  error_hist[i]->SetTitle("Average Auto Correlation");
425  error_hist[i]->GetXaxis()->SetTitle("lag");
426  error_hist[i]->GetYaxis()->SetTitle("Autocorrelation Function");
427  error_hist[i]->Draw("E3 SAME");
428  }
429  if (draw_all)
430  {
431  for (const auto &hist : histograms[i])
432  {
433  hist->SetLineColorAlpha(colour, 0.05f);
434  hist->SetLineWidth(1);
435  hist->Draw("HIST SAME");
436  }
437  }
438  averages[i]->SetLineColor(colour);
439  averages[i]->SetLineWidth(2);
440  averages[i]->SetTitle("Average Auto Correlation");
441  averages[i]->GetXaxis()->SetTitle("lag");
442  averages[i]->GetYaxis()->SetTitle("Autocorrelation Function");
443  averages[i]->Draw("HIST SAME");
444 
445  TString int_str = TString::Format("Average integrated autocorrelation %.2f", averages[i]->Integral());
446  leg->AddEntry(averages[i].get(), hist_labels[i] + int_str, "l");
447  }
448 
449  if (averages.size() > 1)
450  {
451  leg->Draw();
452  }
453 
454  canvas->SaveAs(output_name + ".pdf");
455 }
std::unique_ptr< TGraphAsymmErrors > CalculateMinMaxBand(const std::vector< std::unique_ptr< TH1D >> &histograms, Color_t color)
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.

◆ CreateMinMaxBand()

std::pair<std::unique_ptr<TGraph>, std::unique_ptr<TGraph> > CreateMinMaxBand ( TH1D *  hist,
Color_t  color 
)

HW: Create a band of minimum and maximum values from a histogram.

Definition at line 244 of file PlotMCMCDiag.cpp.

245 {
246  int nBins = hist->GetNbinsX();
247  std::vector<double> x(nBins), ymin(nBins), ymax(nBins);
248 
249  for (int i = 0; i < nBins; ++i)
250  {
251  x[i] = hist->GetBinCenter(i + 1);
252  ymin[i] = hist->GetBinContent(i + 1);
253  ymax[i] = hist->GetBinContent(i + 1);
254  }
255 
256  auto minGraph = std::make_unique<TGraph>(nBins, x.data(), ymin.data());
257  auto maxGraph = std::make_unique<TGraph>(nBins, x.data(), ymax.data());
258 
259  minGraph->SetLineColor(color);
260  maxGraph->SetLineColor(color);
261  minGraph->SetFillColorAlpha(color, 0.3f);
262  maxGraph->SetFillColorAlpha(color, 0.3f);
263 
264  return {std::move(minGraph), std::move(maxGraph)};
265 }

◆ GetMinimumInRange()

double GetMinimumInRange ( TH1D *  hist,
const double  minRange,
const double  maxRange 
)

KS: function which looks for minimum in given range.

Definition at line 20 of file PlotMCMCDiag.cpp.

21 {
22  double MinVale = 1234567890.;
23  TAxis *xaxis = hist->GetXaxis();
24  for (int x = 1; x <= hist->GetNbinsX(); x++)
25  {
26  if (xaxis->GetBinLowEdge(x) > minRange && xaxis->GetBinUpEdge(x) < maxRange)
27  {
28  if (MinVale > hist->GetBinContent(x))
29  MinVale = hist->GetBinContent(x);
30  }
31  }
32  return MinVale;
33 }

◆ IsHistogramAllOnes()

bool IsHistogramAllOnes ( TH1D *  hist,
double  tolerance = 0.001,
int  max_failures = 100 
)

HW: Check if histogram is flat within a given tolerance.

Definition at line 36 of file PlotMCMCDiag.cpp.

37 {
38  int failure_count = 0;
39 
40  for (int bin = 2; bin <= hist->GetNbinsX(); ++bin)
41  {
42  if (std::fabs(hist->GetBinContent(bin) - 1.0) > tolerance)
43  {
44  if (++failure_count > max_failures)
45  {
46  return false;
47  }
48  }
49  }
50  return true;
51 }

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 488 of file PlotMCMCDiag.cpp.

489 {
491  // Require at least one file-title pair
492  if (argc < 3 || (argc - 1) % 2 != 0)
493  {
494  MACH3LOG_ERROR("Wrong number of arguments ({}) provided", argc);
495  MACH3LOG_ERROR("Usage: {} file1 title1 [file2 title2 ... fileN titleN]", argv[0]);
496  MACH3LOG_ERROR("Example: {} DiagMCMC_Output.root PlotName", argv[0]);
497  throw MaCh3Exception(__FILE__, __LINE__);
498  }
499  std::vector<TString> filename;
500  std::vector<TString> title;
501 
502  std::vector<TString> filenames;
503  std::vector<TString> titles;
504 
505  // Parse all pairs (filename, title)
506  for (int i = 1; i < argc; i += 2)
507  {
508  filenames.emplace_back(argv[i]);
509  titles.emplace_back(argv[i + 1]);
510  }
511 
512  // Log how many pairs we found
513  MACH3LOG_INFO("Processing {} file(s)", filenames.size());
514 
515  // Generate plots
516  MakeDiagPlot(filenames, titles);
517  PlotAutoCorr(filenames);
518  PlotAverageACMult(filenames, titles, Form("%s_Average_Auto_Corr", argv[2]), true);
519  return 0;
520 }
void SetMaCh3LoggerFormat()
Set messaging format of the logger.
Definition: MaCh3Logger.h:61
void PlotAverageACMult(std::vector< TString > input_files, std::vector< TString > hist_labels, const TString &output_name, bool draw_min_max=true, bool draw_all=false, bool draw_errors=true)
Calculate mean AC based on all parameters with error band. Great for comparing AC between different c...
void PlotAutoCorr(const std::vector< TString > &fname)
Plot autocorrelation for all params into single plot with color codding helping whether on average it...
void MakeDiagPlot(std::vector< TString > input_files, std::vector< TString > hist_labels)
function which loops over Diag MCMC output and prints everything to PDF

◆ MakeDiagPlot()

void MakeDiagPlot ( std::vector< TString >  input_files,
std::vector< TString >  hist_labels 
)

function which loops over Diag MCMC output and prints everything to PDF

Definition at line 54 of file PlotMCMCDiag.cpp.

56 {
57  auto c1 = std::make_unique<TCanvas>("c1", " ", 0, 0, 800, 630);
58  gStyle->SetOptStat(0); // Set 0 to disable statistic box
59  // To avoid TCanvas::Print> messages
60  gErrorIgnoreLevel = kWarning;
61 
62  TKey *key = nullptr;
63  std::vector<TFile*> infiles;
64  bool add_legend = false;
65 
66  for (const auto& fname : input_files)
67  {
68  TFile* f = M3::Open(fname.Data(), "UPDATE", __FILE__, __LINE__);
69  MACH3LOG_INFO("Adding file {}", fname.Data());
70  infiles.push_back(f);
71  add_legend = (input_files.size() > 1);
72  }
73 
74  constexpr int NFiles_hard = 4;
75  constexpr Color_t Colour[NFiles_hard] = {kRed, kBlue, kGreen, kOrange};
76  constexpr Style_t Style[NFiles_hard] = {kSolid, kDashed, kDotted, kDashDotted};
77  if(input_files.size() > NFiles_hard){
78  MACH3LOG_ERROR("Passed {}, while I have nasty hardcoding for max {} filse", input_files.size(), NFiles_hard);
79  throw MaCh3Exception(__FILE__, __LINE__);
80  }
81  TIter next(infiles[0]->GetListOfKeys());
82  while ((key = static_cast<TKey *>(next())))
83  {
84  std::string dirname = std::string(key->GetName());
85  if (std::string(key->GetClassName()) != "TDirectoryFile")
86  continue;
87  // KS: Script will work with LogL and Batched_means, you can comment it if you are interested in it
88  if ((dirname == "LogL") || (dirname == "Batched_means") || (dirname == "PowerSpectrum"))
89  continue;
90  // KS: Trace wo longer chains is super big, the way to avoid is to plot as png but I don't like png,
91  // keep possibility to skip it
92  // if( (dirname == "Trace") ) continue;
93  infiles[0]->cd(dirname.c_str());
94  TIter nextsub(gDirectory->GetListOfKeys());
95  c1->Print(Form("%s_%s.pdf[", hist_labels[0].Data(), dirname.c_str()), "pdf");
96  TKey *subkey;
97  while ((subkey = static_cast<TKey *>(nextsub())))
98  {
99  auto leg = std::make_unique<TLegend>(0.7, 0.7, 0.9, 0.9);
100  leg->SetFillColorAlpha(kWhite, 0.7f);
101  std::string name = std::string(subkey->GetName());
102  name = dirname + "/" + name;
103  MACH3LOG_INFO("{}", name);
104  if (std::string(subkey->GetClassName()) != "TH1D") {
105  continue;
106  } else {
107  MACH3LOG_WARN("continuing along my way for {}", dirname);
108  }
109  if (name == "AccProb/AcceptanceProbability") continue;
110 
111  std::vector<std::unique_ptr<TH1D>> blarb(infiles.size());
112  MACH3LOG_INFO("Looking for {} from file {}", name.c_str(), input_files[0].Data());
113  blarb[0] = M3::Clone(static_cast<TH1D *>(infiles[0]->Get(name.c_str())));
114  // KS: Some fixe params can go crazy
115  if (std::isnan(blarb[0]->GetBinContent(1)))
116  continue;
117 
118  RemoveFitter(blarb[0].get(), "Fitter");
119  if (dirname == "AccProb")
120  blarb[0]->GetYaxis()->SetRangeUser(0, 1.0);
121 
122  blarb[0]->SetLineStyle(Style[0]);
123  blarb[0]->SetLineColor(Colour[0]);
124  blarb[0]->Draw();
125  leg->AddEntry(blarb[0].get(), hist_labels[0].Data(), "l");
126 
127  for (size_t i = 1; i < infiles.size(); ++i)
128  {
129  blarb[i] = M3::Clone(static_cast<TH1D *>(infiles[i]->Get(name.c_str())));
130 
131  blarb[i]->SetLineStyle(Style[i]);
132  blarb[i]->SetLineColor(Colour[i]);
133 
134  blarb[i]->Draw(i == 0 ? "" : "same");
135  leg->AddEntry(blarb[i].get(), hist_labels[i].Data(), "l");
136  }
137  if (add_legend)
138  {
139  leg->Draw();
140  }
141  c1->Print(Form("%s_%s.pdf", hist_labels[0].Data(), dirname.c_str()), "pdf");
142  }
143  gDirectory->cd("..");
144  c1->Print(Form("%s_%s.pdf]", hist_labels[0].Data(), dirname.c_str()), "pdf");
145  }
146 
147  for (size_t i = 0; i < infiles.size(); ++i)
148  {
149  infiles[i]->Close();
150  delete infiles[i];
151  }
152 }
void RemoveFitter(TH1D *hist, const std::string &name)
KS: Remove fitted TF1 from hist to make comparison easier.
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
constexpr ELineStyle Style[NVars]
Type Get(const YAML::Node &node, const std::string File, const int Line)
Get content of config file.
Definition: YamlHelper.h:291
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.

◆ PlotAutoCorr()

void PlotAutoCorr ( const std::vector< TString > &  fname)

Plot autocorrelation for all params into single plot with color codding helping whether on average it is ok or not.

Definition at line 155 of file PlotMCMCDiag.cpp.

156 {
157  std::vector<TFile*> infile(fname.size());
158  int Nfiles = 0;
159  for(size_t i = 0; i < fname.size(); i++)
160  {
161  infile[i] = M3::Open(fname[i].Data(), "open", __FILE__, __LINE__);
162  Nfiles++;
163  }
164 
165  auto c1 = std::make_unique<TCanvas>("c1", " ", 0, 0, 800, 630);
166  gStyle->SetOptStat(0); // Set 0 to disable statistic box
167  // To avoid TCanvas::Print> messages
168  gErrorIgnoreLevel = kWarning;
169 
170  c1->Print("Auto_Corr_PerFile.pdf[", "pdf");
171  for (int ik = 0; ik < Nfiles; ik++)
172  {
173  TIter next(infile[ik]->GetListOfKeys());
174 
175  TKey *key;
176  // KS: Keep track of histos to delete them
177  std::vector<std::unique_ptr<TH1D>> histos;
178  while ((key = static_cast<TKey *>(next())))
179  {
180  std::string dirname = std::string(key->GetName());
181 
182  // KS: Script We are only interested in auto corr
183  if ((dirname != "Auto_corr"))
184  continue;
185 
186  infile[ik]->cd(dirname.c_str());
187  TIter nextsub(gDirectory->GetListOfKeys());
188 
189  TKey *subkey;
190  bool FirstTime = true;
191  while ((subkey = static_cast<TKey *>(nextsub())))
192  {
193  std::string name = std::string(subkey->GetName());
194  name = dirname + "/" + name;
195 
196  if (std::string(subkey->GetClassName()) != "TH1D")
197  continue;
198  MACH3LOG_DEBUG("{}", name.c_str());
199  histos.push_back(M3::Clone(static_cast<TH1D *>(infile[ik]->Get(name.c_str()))));
200  // KS: Some fixe pramas can go crazy
201  if (std::isnan(histos.back()->GetBinContent(1)))
202  continue;
203  // KS: This is unfortunately hardcoded, need to find better way to write this
204  // histos.back()->GetListOfFunctions()->ls();
205  RemoveFitter(histos.back().get(), "Fitter");
206  // KS: Very arbitrary check if autocorrelation dropped below threshold in half of consider LagL
207  double UpperEdgeBeforeLast = histos.back()->GetXaxis()->GetBinUpEdge(histos.back()->GetNbinsX()/2);
208  double MinValue = GetMinimumInRange(histos.back().get(), 0, UpperEdgeBeforeLast);
209 
210  if (MinValue >= 0.80)
211  histos.back()->SetLineColor(kRed);
212  else if (MinValue >= 0.40 && MinValue < 0.80)
213  histos.back()->SetLineColor(kOrange);
214  else if (MinValue > 0.20 && MinValue < 0.40)
215  histos.back()->SetLineColor(kYellow);
216  else if (MinValue <= 0.20)
217  histos.back()->SetLineColor(kGreen);
218  histos.back()->GetXaxis()->UnZoom();
219 
220  if (FirstTime)
221  histos.back()->SetTitle(Form("Auto_Corr_%s.pdf", fname[ik].Data()));
222 
223  histos.back()->SetLineStyle(kDashed);
224 
225  if (FirstTime)
226  histos.back()->Draw();
227 
228  TString integral = Form("Integral: %.2f", histos.back()->Integral());
229 
230  if (!FirstTime)
231  histos.back()->Draw("same");
232  FirstTime = false;
233  }
234  gDirectory->cd("..");
235  }
236  c1->Print("Auto_Corr_PerFile.pdf", "pdf");
237  }
238  c1->Print("Auto_Corr_PerFile.pdf]", "pdf");
239 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
double GetMinimumInRange(TH1D *hist, const double minRange, const double maxRange)
KS: function which looks for minimum in given range.

◆ PlotAverageACMult()

void PlotAverageACMult ( std::vector< TString >  input_files,
std::vector< TString >  hist_labels,
const TString &  output_name,
bool  draw_min_max = true,
bool  draw_all = false,
bool  draw_errors = true 
)

Calculate mean AC based on all parameters with error band. Great for comparing AC between different chains.

Definition at line 458 of file PlotMCMCDiag.cpp.

464 {
465  // Process first folder
466  std::vector<std::vector<std::unique_ptr<TH1D>>> histograms;
467  std::vector<std::unique_ptr<TH1D>> averages;
468 
469  for (int i = 0; i < static_cast<int>(input_files.size()); i++)
470  {
471  TString folder = input_files[i];
472  MACH3LOG_INFO("Folder : {}", folder);
473 
474  if (folder.IsNull())
475  {
476  MACH3LOG_WARN("Skipping empty or dummy folder: {}", folder.Data());
477  continue;
478  }
479 
480  std::vector<std::unique_ptr<TH1D>> histograms_i;
481  averages.emplace_back((AutocorrProcessInputs(folder, histograms_i)));
482  histograms.push_back(std::move(histograms_i));
483  }
484 
485  CompareAverageAC(histograms, averages, hist_labels, output_name, draw_min_max, draw_all, draw_errors);
486 }
std::unique_ptr< TH1D > AutocorrProcessInputs(const TString &input_file, std::vector< std::unique_ptr< TH1D >> &histograms)
void CompareAverageAC(const std::vector< std::vector< std::unique_ptr< TH1D >>> &histograms, const std::vector< std::unique_ptr< TH1D >> &averages, const std::vector< TString > &hist_labels, const TString &output_name, bool draw_min_max=true, bool draw_all=false, bool draw_errors=true)

◆ ProcessAutoCorrelationDirectory()

void ProcessAutoCorrelationDirectory ( TDirectoryFile *  autocor_dir,
std::unique_ptr< TH1D > &  average_hist,
int &  parameter_count,
std::vector< std::unique_ptr< TH1D >> &  histograms 
)

Loop over directory get histograms into vector and add to averaged hist.

Definition at line 303 of file PlotMCMCDiag.cpp.

307 {
308  TIter next(autocor_dir->GetListOfKeys());
309  TKey *key;
310 
311  while ((key = dynamic_cast<TKey *>(next())))
312  {
313  TH1D *current_hist = nullptr;
314  autocor_dir->GetObject(key->GetName(), current_hist);
315 
316  if (!current_hist ||
317  current_hist->GetMaximum() <= 0 ||
318  IsHistogramAllOnes(current_hist))
319  {
320  continue;
321  }
322 
323  histograms.push_back(M3::Clone(current_hist));
324  if (!average_hist) {
325  average_hist = M3::Clone(current_hist);
326  } else {
327  average_hist->Add(current_hist);
328  }
329  parameter_count++;
330  }
331 }
bool IsHistogramAllOnes(TH1D *hist, double tolerance=0.001, int max_failures=100)
HW: Check if histogram is flat within a given tolerance.

◆ ProcessDiagnosticFile()

void ProcessDiagnosticFile ( const TString &  file_path,
std::unique_ptr< TH1D > &  average_hist,
int &  parameter_count,
std::vector< std::unique_ptr< TH1D >> &  histograms 
)

Definition at line 333 of file PlotMCMCDiag.cpp.

337 {
338  std::unique_ptr<TFile> input_file(TFile::Open(file_path));
339  if (!input_file || input_file->IsZombie())
340  {
341  throw MaCh3Exception(__FILE__, __LINE__, "Could not open file: " + std::string(file_path.Data()));
342  }
343 
344  TDirectoryFile *autocor_dir = nullptr;
345  input_file->GetObject("Auto_corr", autocor_dir);
346 
347  if (!autocor_dir)
348  {
349  throw MaCh3Exception(__FILE__, __LINE__,
350  "Auto_corr directory not found in file: " + std::string(file_path.Data()));
351  }
352 
353  ProcessAutoCorrelationDirectory(autocor_dir, average_hist, parameter_count, histograms);
354 }
void ProcessAutoCorrelationDirectory(TDirectoryFile *autocor_dir, std::unique_ptr< TH1D > &average_hist, int &parameter_count, std::vector< std::unique_ptr< TH1D >> &histograms)
Loop over directory get histograms into vector and add to averaged hist.