MaCh3  2.2.3
Reference Guide
Functions | Variables
PlotLLH.cpp File Reference
#include <algorithm>
#include <iomanip>
#include "plottingUtils/plottingManager.h"
Include dependency graph for PlotLLH.cpp:

Go to the source code of this file.

Functions

void getSplitSampleStack (int fileIdx, std::string parameterName, TH1D LLH_allSams, std::vector< float > &cumSums, std::vector< bool > &drawLabel, THStack *sampleStack, TLegend *splitSamplesLegend, float baselineLLH_main=0.00001)
 
void drawRatioStack (THStack *ratioCompStack)
 
void makeLLHScanComparisons (const std::string &paramName, const std::string &LLHType, const std::string &outputFileName, const std::unique_ptr< TCanvas > &canv, const std::unique_ptr< TPad > &LLHPad, const std::unique_ptr< TPad > &ratioPad)
 
void makeSplitSampleLLHScanComparisons (const std::string &paramName, const std::string &outputFileName, const std::unique_ptr< TCanvas > &canv, const std::unique_ptr< TPad > &LLHPad, const std::unique_ptr< TPad > &ratioPad)
 
int PlotLLH ()
 
int main (int argc, char **argv)
 

Variables

double ratioPlotSplit
 
double yTitleOffset
 
double sampleLabelThreshold
 
int lineWidth
 
bool totalOnSplitPlots
 
bool sameAxis
 
double ratioLabelScaling
 
MaCh3Plotting::PlottingManagerman
 

Detailed Description

Author
Ewan Miller

Definition in file PlotLLH.cpp.

Function Documentation

◆ drawRatioStack()

void drawRatioStack ( THStack *  ratioCompStack)

Definition at line 94 of file PlotLLH.cpp.

94  {
95  // do this so 1.0 is in the middle of the plot vertically
96  double stackMax = ratioCompStack->GetMaximum("NOSTACK");
97  double stackMin = ratioCompStack->GetMinimum("NOSTACK");
98 
99  double stackLim = std::max(std::abs(1.0 - stackMax), std::abs(1.0 - stackMin));
100 
101  ratioCompStack->SetMinimum(1.0 - 1.05 * stackLim);
102  ratioCompStack->SetMaximum(1.0 + 1.05 * stackLim);
103 
104  // draw it
105  ratioCompStack->Draw(Form("NOSTACK%s", man->getDrawOptions().c_str()));
106 
107  // make it look a bit nicer
108  ratioCompStack->GetXaxis()->SetLabelSize(ratioLabelScaling *
109  ratioCompStack->GetXaxis()->GetLabelSize());
110  ratioCompStack->GetXaxis()->SetTitleSize(ratioLabelScaling *
111  ratioCompStack->GetXaxis()->GetLabelSize());
112  ratioCompStack->GetXaxis()->SetTitle("Parameter Variation");
113 
114  ratioCompStack->GetYaxis()->SetLabelSize(ratioLabelScaling *
115  ratioCompStack->GetYaxis()->GetLabelSize());
116  ratioCompStack->GetYaxis()->SetTitleSize(ratioLabelScaling *
117  ratioCompStack->GetYaxis()->GetLabelSize());
118  ratioCompStack->GetYaxis()->SetTitleOffset(yTitleOffset);
119  ratioCompStack->GetYaxis()->SetNdivisions(5, 2, 0);
120 }
double ratioLabelScaling
Definition: PlotLLH.cpp:23
double yTitleOffset
Definition: PlotLLH.cpp:17
MaCh3Plotting::PlottingManager * man
Definition: PlotLLH.cpp:25
const std::string getDrawOptions() const

◆ getSplitSampleStack()

void getSplitSampleStack ( int  fileIdx,
std::string  parameterName,
TH1D  LLH_allSams,
std::vector< float > &  cumSums,
std::vector< bool > &  drawLabel,
THStack *  sampleStack,
TLegend *  splitSamplesLegend,
float  baselineLLH_main = 0.00001 
)

Definition at line 27 of file PlotLLH.cpp.

31  {
32  std::vector<std::string> sampNames = man->input().getTaggedSamples(man->getOption<std::vector<std::string>>("sampleTags"));
33  size_t nSamples = sampNames.size();
34 
35  cumSums.resize(nSamples);
36  drawLabel.resize(nSamples);
37 
38  float LLH_main_integ = LLH_allSams.Integral();
39  float cumSum = 0.0;
40  int nBins = LLH_allSams.GetNbinsX();
41 
42  MACH3LOG_DEBUG("getting split sample THStack for {} known samples", nSamples);
43  for (uint i = 0; i < nSamples; i++)
44  {
45  std::string sampName = sampNames[i];
46 
47  MACH3LOG_DEBUG(" on sample {}/{}: {}", i, nSamples, sampName);
48 
49  TH1D *LLH_indivSam =
50  new TH1D(man->input().getSampleSpecificLLHScan_TH1D(fileIdx, parameterName, sampName));
51  LLH_indivSam->SetName(Form("%i_%s_%s", fileIdx, parameterName.c_str(), sampName.c_str()));
52  LLH_indivSam->SetBit(kCanDelete);
53 
54  // the hist for this sample didn't exist so move on
55  if (LLH_indivSam->GetNbinsX() == 1)
56  {
57  MACH3LOG_DEBUG(" sample hist had only 1 bin - assuming it doesn't exist");
58  delete LLH_indivSam;
59  drawLabel[i] = false;
60  cumSums[i] = -999.9;
61  continue;
62  }
63 
64  LLH_indivSam->SetStats(0);
65  LLH_indivSam->SetLineColor(TColor::GetColorPalette(
66  floor(static_cast<float>(i) * TColor::GetNumberOfColors() / static_cast<float>(nSamples))));
67  LLH_indivSam->SetFillColor(TColor::GetColorPalette(
68  floor(static_cast<float>(i) * TColor::GetNumberOfColors() / static_cast<float>(nSamples))));
69  sampleStack->Add(LLH_indivSam);
70  splitSamplesLegend->AddEntry(LLH_indivSam, man->style().prettifySampleName(sampName).c_str(), "lf");
71 
72  float lastBinLLH = LLH_indivSam->GetBinContent(nBins);
73  cumSum += lastBinLLH;
74 
75  MACH3LOG_DEBUG(" Last bin LLH = {} :: cumulative LLH = {}", lastBinLLH, cumSum);
76  MACH3LOG_DEBUG(" LLH fraction = {} / {} = {}", LLH_indivSam->Integral(), LLH_main_integ, LLH_indivSam->Integral() / LLH_main_integ);
77 
78  cumSums[i] = cumSum;
79 
80  // dont draw a label if the likelihood contribution is less than threshold%
81  if ((LLH_indivSam->Integral() / LLH_main_integ > sampleLabelThreshold) &&
82  (LLH_indivSam->Integral() / baselineLLH_main > sampleLabelThreshold))
83  {
84  drawLabel[i] = true;
85  } else {
86  drawLabel[i] = false;
87  }
88  MACH3LOG_DEBUG(" drawLabel = {}", drawLabel.back());
89  MACH3LOG_DEBUG("");
90  }
91 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:24
double sampleLabelThreshold
Definition: PlotLLH.cpp:18
std::vector< std::string > getTaggedSamples(const std::vector< std::string > &tags, std::string checkType="all") const
Get all samples which have some set of tags.
Definition: inputManager.h:493
TH1D getSampleSpecificLLHScan_TH1D(const int fileNum, const std::string &paramName, const std::string &sample) const
Definition: inputManager.h:405
const InputManager & input() const
Get the InputManager contained within this PlottingManager, for doing input related things.
T getOption(const std::string &option)
Get a specific option from the config for this executable.
const StyleManager & style() const
Get the StyleManager contained within this PlottingManager, for doing style related things.
std::string prettifySampleName(const std::string &origName) const
Convert hideous and vulgar internal sample name into a beautiful presentable name.
Definition: styleManager.h:48

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 453 of file PlotLLH.cpp.

453  {
456 
458  man->parseInputs(argc, argv);
459 
460  man->setExec("PlotLLH");
461 
462  ratioPlotSplit = man->getOption<double>("ratioPlotSplit");
463  yTitleOffset = man->getOption<double>("yTitleOffset");
464  sampleLabelThreshold = man->getOption<double>("sampleLabelThreshold");
465  lineWidth = man->getOption<int>("lineWidth");
466  totalOnSplitPlots = man->getOption<bool>("totalOnSplitPlots");
467  sameAxis = man->getOption<bool>("sameAxis");
468 
469  // scale the ratio plot labels by this much to make them same size as the normal plot
470  ratioLabelScaling = (1.0 / ratioPlotSplit - 1.0);
471 
472  return PlotLLH();
473 }
void SetMaCh3LoggerFormat()
Set messaging format of the logger.
Definition: MaCh3Logger.h:51
bool sameAxis
Definition: PlotLLH.cpp:21
double ratioPlotSplit
Definition: PlotLLH.cpp:16
int PlotLLH()
Definition: PlotLLH.cpp:386
bool totalOnSplitPlots
Definition: PlotLLH.cpp:20
int lineWidth
Definition: PlotLLH.cpp:19
The main class to be used in plotting scripts.
void setExec(const std::string &execName)
Internally set the name of the executable that manager is being used in.
void parseInputs(int argc, char *const *argv)
Parse command line arguments.
void MaCh3Welcome()
KS: Prints welcome message with MaCh3 logo.
Definition: Monitor.cpp:12

◆ makeLLHScanComparisons()

void makeLLHScanComparisons ( const std::string &  paramName,
const std::string &  LLHType,
const std::string &  outputFileName,
const std::unique_ptr< TCanvas > &  canv,
const std::unique_ptr< TPad > &  LLHPad,
const std::unique_ptr< TPad > &  ratioPad 
)

Definition at line 122 of file PlotLLH.cpp.

127  {
128  // will use these to make comparisons of likelihoods
129  auto compStack = std::make_unique<THStack>("LLH_Stack", "");
130  auto ratioCompStack = std::make_unique<THStack>("LLH_Ratio_Stack", "");
131  auto legend = std::make_unique<TLegend>(0.3, 0.6, 0.7, 0.8);
132 
133  // get the sample reweight hist from the main file
134  TH1D LLH_main = man->input().getLLHScan_TH1D(0, paramName, LLHType);
135  LLH_main.SetStats(0);
136 
137  LLH_main.SetLineColor(kBlack);
138  compStack->Add(&LLH_main);
139  legend->AddEntry(&LLH_main, man->getFileLabel(0).c_str(), "l");
140 
141  int nBins = LLH_main.GetNbinsX();
142 
143  // go through the other files
144  for (unsigned int extraFileIdx = 1; extraFileIdx < man->input().getNInputFiles(); extraFileIdx++)
145  {
146  int originalErrorLevel = gErrorIgnoreLevel;
147  gErrorIgnoreLevel = kFatal;
148 
149  TH1D *compHist = new TH1D(man->input().getLLHScan_TH1D(extraFileIdx, paramName, LLHType));
150  compHist->SetName(Form("LLHScan_%s_%s_%d", paramName.c_str(), LLHType.c_str(), extraFileIdx));
151  compHist->SetBit(kCanDelete); // <- will allow this to be deleted by root once done plotting
152  gErrorIgnoreLevel = originalErrorLevel;
153  if (compHist->GetNbinsX() == 0)
154  continue;
155 
156  // make them look different to each other
157  compHist->SetLineColor(
158  TColor::GetColorPalette(floor(static_cast<float>(extraFileIdx) * TColor::GetNumberOfColors() /
159  static_cast<float>(man->input().getNInputFiles()))));
160  compHist->SetLineStyle(2 + extraFileIdx % 9);
161  compHist->SetLineWidth(lineWidth);
162 
163  TH1D *divHist = static_cast<TH1D*>(compHist->Clone(Form("RatioHist_%i", extraFileIdx)));
164  divHist->SetBit(kCanDelete);
165 
166  if (man->getPlotRatios())
167  {
168  // get the ratio hist
169  divHist->Divide(compHist, &LLH_main);
170  ratioCompStack->Add(divHist);
171  }
172 
173  // add it to the comparisson hstack and legend
174  compStack->Add(compHist);
175  legend->AddEntry(compHist, man->getFileLabel(extraFileIdx).c_str(), "l");
176  }
177 
178  // draw the log likelihoods
179  canv->cd();
180  canv->Draw();
181  LLHPad->Draw();
182  LLHPad->cd();
183  if (man->getDrawGrid())
184  LLHPad->SetGrid();
185  compStack->Draw(Form("NOSTACK%s", man->getDrawOptions().c_str()));
186  if (!man->getPlotRatios())
187  compStack->GetXaxis()->SetTitle("Parameter Variation");
188  compStack->GetYaxis()->SetTitle(Form("-2LLH_{%s}", LLHType.c_str()));
189  compStack->SetTitle(man->style().prettifyParamName(paramName).c_str());
190  compStack->GetYaxis()->SetTitleOffset(yTitleOffset);
191  legend->Draw();
192 
193  // add the ratio plot if specified
194  if (man->getPlotRatios())
195  {
196  canv->cd();
197  ratioPad->Draw();
198  ratioPad->cd();
199  if (man->getDrawGrid())
200  ratioPad->SetGrid();
201 
202  drawRatioStack(ratioCompStack.get());
203 
204  // add horizontal line at 1 for reference
205  TLine line = TLine();
206  line.SetLineColor(kBlack);
207  line.SetLineWidth(lineWidth);
208  line.DrawLine(LLH_main.GetBinLowEdge(1), 1.0, LLH_main.GetBinLowEdge(nBins + 1), 1.0);
209  }
210 
211  // save to the output file
212  canv->SaveAs(outputFileName.c_str());
213 }
TCanvas * canv
void drawRatioStack(THStack *ratioCompStack)
Definition: PlotLLH.cpp:94
TH1D getLLHScan_TH1D(const int fileNum, const std::string &paramName, const std::string &LLHType) const
Definition: inputManager.h:368
size_t getNInputFiles() const
Definition: inputManager.h:473
const std::string getFileLabel(int i) const
std::string prettifyParamName(const std::string &origName) const
Convert hideous and vulgar internal parameter name into a beautiful presentable name.
Definition: styleManager.h:39

◆ makeSplitSampleLLHScanComparisons()

void makeSplitSampleLLHScanComparisons ( const std::string &  paramName,
const std::string &  outputFileName,
const std::unique_ptr< TCanvas > &  canv,
const std::unique_ptr< TPad > &  LLHPad,
const std::unique_ptr< TPad > &  ratioPad 
)

Definition at line 215 of file PlotLLH.cpp.

219  {
220  MACH3LOG_DEBUG(" Making split sample LLH comparison");
221  canv->Clear();
222  canv->Draw();
223 
224  canv->Divide(man->getNFiles());
225 
226  // get the sample hist from the main file
227  TH1D LLH_main = man->input().getLLHScan_TH1D(0, paramName, "sample");
228  if (LLH_main.GetNbinsX() == 1)
229  {
230  MACH3LOG_DEBUG(" Main LLH had only 1 bin, assuming it doesn't exist");
231  return;
232  }
233 
234  THStack *baseSplitSamplesStack = new THStack(
235  paramName.c_str(), Form("%s - %s", paramName.c_str(), man->getFileLabel(0).c_str()));
236  TLegend *baseSplitSamplesLegend = new TLegend(0.37, 0.475, 0.63, 0.9);
237 
238  if (man->getDrawGrid())
239  canv->cd(1)->SetGrid();
240  else
241  canv->cd(1);
242 
243  std::vector<float> cumSums;
244  std::vector<bool> drawLabel;
245 
246  getSplitSampleStack(0, paramName, LLH_main, cumSums, drawLabel, baseSplitSamplesStack,
247  baseSplitSamplesLegend);
248 
249  baseSplitSamplesStack->Draw(man->getDrawOptions().c_str());
250 
251  if (totalOnSplitPlots)
252  {
253  LLH_main.SetLineWidth(1); // undo SetLineWidth that was done above
254  LLH_main.Draw(Form("same%s", man->getDrawOptions().c_str()));
255  baseSplitSamplesLegend->AddEntry(&LLH_main, "All Samples", "l");
256  }
257 
258  baseSplitSamplesLegend->Draw();
259 
260  TLatex *label = new TLatex;
261  // format the label
262  label->SetTextAlign(11);
263  label->SetTextAngle(-55);
264  label->SetTextSize(0.012);
265 
266  // need to draw the labels after other stuff or they dont show up
267  for (uint i = 0; i < man->input().getTaggedSamples(man->getOption<std::vector<std::string>>("sampleTags")).size(); i++)
268  {
269  MACH3LOG_DEBUG(" Will I draw the label for sample {}??", i);
270  std::string sampName = man->input().getTaggedSamples(man->getOption<std::vector<std::string>>("sampleTags"))[i];
271  if (!drawLabel[i])
272  {
273  MACH3LOG_DEBUG(" - Not drawing label");
274  continue;
275  }
276  MACH3LOG_DEBUG(" - Drawing label");
277 
278  label->DrawLatex(LLH_main.GetBinLowEdge(LLH_main.GetNbinsX() + 1), cumSums[i],
279  Form("#leftarrow%s", man->style().prettifySampleName(sampName).c_str()));
280  MACH3LOG_DEBUG(" I drew the label!");
281  }
282 
283  // now we plot the comparisson file plots
284  for (unsigned int extraFileIdx = 1; extraFileIdx < man->getNFiles(); extraFileIdx++)
285  {
286  MACH3LOG_DEBUG(" - Adding plot for additional file {}", extraFileIdx);
287  canv->cd(1 + extraFileIdx);
288 
289  std::vector<float> extraCumSums;
290  std::vector<bool> extraDrawLabel;
291 
292  THStack *splitSamplesStack =
293  new THStack(paramName.c_str(),
294  Form("%s - %s", paramName.c_str(), man->getFileLabel(extraFileIdx).c_str()));
295  TLegend *splitSamplesLegend = new TLegend(0.37, 0.475, 0.63, 0.9);
296 
297  splitSamplesStack->SetBit(kCanDelete);
298  splitSamplesLegend->SetBit(kCanDelete);
299 
300  TH1D compLLH_main = man->input().getLLHScan_TH1D(extraFileIdx, paramName, "sample");
301  if (compLLH_main.GetNbinsX() == 1)
302  {
303  delete splitSamplesStack;
304  delete splitSamplesLegend;
305  continue;
306  }
307 
308  // if on the same y axis, also check that the contribution of the sample compared to the
309  // baseline LLH integral is above the threshold otherwise the labels might get very crowded if
310  // the comparisson LLH is much smaller than the baseline one
311  if (sameAxis)
312  getSplitSampleStack(extraFileIdx, paramName, compLLH_main, extraCumSums, extraDrawLabel,
313  splitSamplesStack, splitSamplesLegend, LLH_main.Integral());
314  else
315  getSplitSampleStack(extraFileIdx, paramName, compLLH_main, extraCumSums, extraDrawLabel,
316  splitSamplesStack, splitSamplesLegend);
317 
318  // go to the pad for the histograms
319  if (man->getDrawGrid())
320  LLHPad->SetGrid();
321  LLHPad->Draw();
322  LLHPad->cd();
323 
324  splitSamplesStack->Draw(man->getDrawOptions().c_str());
325  splitSamplesStack->GetXaxis()->SetTitle("Parameter Variation");
326  splitSamplesStack->GetYaxis()->SetTitle("-2LLH_{sam}");
327  if (sameAxis)
328  splitSamplesStack->SetMaximum(baseSplitSamplesStack->GetMaximum());
329 
330  if (totalOnSplitPlots)
331  {
332  compLLH_main.Draw(Form("same%s", man->getDrawOptions().c_str()));
333  compLLH_main.SetLineColor(kBlack);
334  splitSamplesLegend->AddEntry(&compLLH_main, "All Samples", "l");
335  }
336  splitSamplesLegend->Draw();
337 
338  // need to draw the labels after other stuff or they dont show up
339  for (uint i = 0; i < man->input().getTaggedSamples(man->getOption<std::vector<std::string>>("sampleTags")).size(); i++)
340  {
341  std::string sampName = man->input().getTaggedSamples(man->getOption<std::vector<std::string>>("sampleTags"))[i];
342  if (!drawLabel[i])
343  continue;
344  label->DrawLatex(compLLH_main.GetBinLowEdge(compLLH_main.GetNbinsX() + 1), extraCumSums[i],
345  Form("#leftarrow%s", man->style().prettifySampleName(sampName).c_str()));
346  }
347 
348  if (man->getPlotRatios())
349  {
350  THStack *splitSamplesStackRatios = new THStack(paramName.c_str(), "");
351 
352  TList *baselineHistList = baseSplitSamplesStack->GetHists();
353  TList *compHistList = splitSamplesStack->GetHists();
354 
355  for (uint sampleIdx = 0; sampleIdx < man->input().getTaggedSamples(man->getOption<std::vector<std::string>>("sampleTags")).size(); sampleIdx++)
356  {
357  TH1D *divHist = new TH1D(
358  Form("%s_%s_splitDiv_%i", paramName.c_str(), man->getFileLabel(extraFileIdx).c_str(),
359  sampleIdx),
360  Form("%s_%s_splitDiv", paramName.c_str(), man->getFileLabel(extraFileIdx).c_str()),
361  compLLH_main.GetNbinsX(), compLLH_main.GetBinLowEdge(1),
362  compLLH_main.GetBinLowEdge(compLLH_main.GetNbinsX() + 1));
363  divHist->Divide(static_cast<TH1D*>(compHistList->At(sampleIdx)),
364  static_cast<TH1D*>(baselineHistList->At(sampleIdx)));
365  splitSamplesStackRatios->Add(divHist);
366  divHist->SetLineColor((static_cast<TH1D*>(compHistList->At(sampleIdx))->GetLineColor()));
367  }
368 
369  canv->cd(2 + extraFileIdx);
370  if (man->getDrawGrid())
371  ratioPad->SetGrid();
372  ratioPad->Draw();
373  ratioPad->cd();
374 
375  drawRatioStack(splitSamplesStackRatios);
376  }
377  }
378 
379  canv->SaveAs(outputFileName.c_str());
380 
381  delete label;
382  delete baseSplitSamplesStack;
383  delete baseSplitSamplesLegend;
384 }
void getSplitSampleStack(int fileIdx, std::string parameterName, TH1D LLH_allSams, std::vector< float > &cumSums, std::vector< bool > &drawLabel, THStack *sampleStack, TLegend *splitSamplesLegend, float baselineLLH_main=0.00001)
Definition: PlotLLH.cpp:27

◆ PlotLLH()

int PlotLLH ( )

Definition at line 386 of file PlotLLH.cpp.

386  {
387  man->style().setPalette(man->getOption<std::string>("colorPalette"));
388 
389  // make the canvas and other plotting stuff
390  auto canv = std::make_unique<TCanvas>("canv", "", 1024, 1024);
391  gStyle->SetOptTitle(2);
392  // split up the canvas so can have side by side plots, one for each file
393  auto splitSamplesCanv = std::make_unique<TCanvas>("splitSampCanv", "", 4096 * man->getNFiles(), 4096);
394 
395  // split the canvas if plotting ratios
396  std::unique_ptr<TPad> LLHPad, ratioPad;
397  if (man->getPlotRatios())
398  {
399  LLHPad = std::make_unique<TPad>("LLHPad", "LLHPad", 0.0, ratioPlotSplit, 1.0, 1.0);
400  LLHPad->SetBottomMargin(0.0);
401  ratioPad = std::make_unique<TPad>("ratioPad", "ratioPad", 0.0, 0.0, 1.0, ratioPlotSplit);
402  ratioPad->SetTopMargin(0.0);
403  ratioPad->SetBottomMargin(0.3);
404  } else {
405  LLHPad = std::make_unique<TPad>("AllSampPad", "AllSampPad", 0.0, 0.0, 1.0, 1.0);
406  ratioPad = std::make_unique<TPad>("AllSampRatioPad", "AllSampRatioPad", 0.0, 0.0, 0.0, 0.0);
407  }
408 
409  canv->SaveAs((man->getOutputName("_Sample") + "[").c_str());
410  canv->SaveAs((man->getOutputName("_Penalty") + "[").c_str());
411  canv->SaveAs((man->getOutputName("_Total") + "[").c_str());
412 
413  if (man->getSplitBySample())
414  canv->SaveAs((man->getOutputName("_bySample") + "[").c_str());
415 
416  for( std::string par: man->getOption<std::vector<std::string>>("parameterTags")) std::cout << par << ", ";
417  // loop over the spline parameters
418  for (std::string paramName : man->input().getTaggedParameters(man->getOption<std::vector<std::string>>("parameterTags")))
419  {
420  MACH3LOG_DEBUG("working on parameter {}", paramName);
421  // ###############################################################
422  // First lets do just the straight up likelihoods from all samples
423  // ###############################################################
424 
425  makeLLHScanComparisons(paramName, "sample", man->getOutputName("_Sample"), canv, LLHPad,
426  ratioPad);
427  makeLLHScanComparisons(paramName, "penalty", man->getOutputName("_Penalty"), canv, LLHPad,
428  ratioPad);
429  makeLLHScanComparisons(paramName, "total", man->getOutputName("_Total"), canv, LLHPad,
430  ratioPad);
431 
432  // #########################################
433  // ## now lets make plots split by sample ##
434  // #########################################
435  if (man->getSplitBySample())
436  {
437  makeSplitSampleLLHScanComparisons(paramName, man->getOutputName("_bySample"),
438  splitSamplesCanv, LLHPad, ratioPad);
439  }
440  }
441 
442  canv->SaveAs((man->getOutputName("_Sample") + "]").c_str());
443  canv->SaveAs((man->getOutputName("_Penalty") + "]").c_str());
444  canv->SaveAs((man->getOutputName("_Total") + "]").c_str());
445  if (man->getSplitBySample())
446  canv->SaveAs((man->getOutputName("_bySample") + "]").c_str());
447 
448  if(man != nullptr) delete man;
449 
450  return 0;
451 }
void makeSplitSampleLLHScanComparisons(const std::string &paramName, const std::string &outputFileName, const std::unique_ptr< TCanvas > &canv, const std::unique_ptr< TPad > &LLHPad, const std::unique_ptr< TPad > &ratioPad)
Definition: PlotLLH.cpp:215
void makeLLHScanComparisons(const std::string &paramName, const std::string &LLHType, const std::string &outputFileName, const std::unique_ptr< TCanvas > &canv, const std::unique_ptr< TPad > &LLHPad, const std::unique_ptr< TPad > &ratioPad)
Definition: PlotLLH.cpp:122
std::vector< std::string > getTaggedParameters(const std::vector< std::string > &tags, std::string checkType="all") const
Get all parameters which have some set of tags.
Definition: inputManager.h:482
const std::string getOutputName() const
Get the straight up output file name with no bells or whistles, just the file extension.
void setPalette(const int rootPlotStyle) const
Set the root colour palette to one of the default root pallettes as defined in (root docs)[https://ro...

Variable Documentation

◆ lineWidth

int lineWidth

Definition at line 19 of file PlotLLH.cpp.

◆ man

Definition at line 25 of file PlotLLH.cpp.

◆ ratioLabelScaling

double ratioLabelScaling

Definition at line 23 of file PlotLLH.cpp.

◆ ratioPlotSplit

double ratioPlotSplit

Definition at line 16 of file PlotLLH.cpp.

◆ sameAxis

bool sameAxis

Definition at line 21 of file PlotLLH.cpp.

◆ sampleLabelThreshold

double sampleLabelThreshold

Definition at line 18 of file PlotLLH.cpp.

◆ totalOnSplitPlots

bool totalOnSplitPlots

Definition at line 20 of file PlotLLH.cpp.

◆ yTitleOffset

double yTitleOffset

Definition at line 17 of file PlotLLH.cpp.