MaCh3  2.2.3
Reference Guide
MCMCProcessor.h
Go to the documentation of this file.
1 #pragma once
2 
3 // C++ includes
4 #include <complex>
5 #include <cstdio>
6 
7 // MaCh3 includes
11 
13 // ROOT includes
14 #include "TFile.h"
15 #include "TBranch.h"
16 #include "TCanvas.h"
17 #include "TLine.h"
18 #include "TLegend.h"
19 #include "TString.h"
20 #include "TGraphErrors.h"
21 #include "TVectorD.h"
22 #include "TColor.h"
23 #include "TStyle.h"
24 #include "TStopwatch.h"
25 #include "TText.h"
26 #include "TGaxis.h"
27 #include "TTree.h"
28 #include "TROOT.h"
29 #include "TKey.h"
30 #include "TRandom3.h"
31 #include "TGraphPolar.h"
32 #include "TCandle.h"
33 #include "TMath.h"
34 #include "TMatrixDSymEigen.h"
35 #include "TVirtualFFT.h"
37 
38 
39 //KS: Joy of forward declaration https://gieseanw.wordpress.com/2018/02/25/the-joys-of-forward-declarations-results-from-the-real-world/
40 class TChain;
41 class TF1;
42 
46 
49  kXSecPar = 0,
50  kNDPar = 1,
51  kFDDetPar = 2,
52 
53  kNParameterEnum = 3 //KS: keep it at the end to keep track of all parameters
54 };
55 
62  public:
65  MCMCProcessor(const std::string &InputFile);
67  virtual ~MCMCProcessor();
68 
70  void Initialise();
72  void MakePostfit(const std::map<std::string, std::pair<double, double>>& Edges = {});
75  void MakeCovariance();
77  void CacheSteps();
80  void MakeCovariance_MP(const bool Mute = false);
84  void MakeSubOptimality(const int NIntervals = 10);
85 
87  void ResetHistograms();
88 
90  void DrawPostfit();
92  void MakeViolin();
97  void MakeCredibleIntervals(const std::vector<double>& CredibleIntervals = {0.99, 0.90, 0.68 },
98  const std::vector<Color_t>& CredibleIntervalsColours = {kCyan+4, kCyan-2, kCyan-10},
99  const bool CredibleInSigmas = false
100  );
102  void DrawCovariance();
104  void MakeCovarianceYAML(const std::string& OutputYAMLFile, const std::string& MeansMethod) const;
110  void MakeCredibleRegions(const std::vector<double>& CredibleRegions = {0.99, 0.90, 0.68},
111  const std::vector<Style_t>& CredibleRegionStyle = {kDashed, kSolid, kDotted},
112  const std::vector<Color_t>& CredibleRegionColor = {kGreen-3, kGreen-10, kGreen},
113  const bool CredibleInSigmas = false,
114  const bool Draw2DPosterior = true,
115  const bool DrawBestFit = true
116  );
124  void MakeTrianglePlot(const std::vector<std::string>& ParNames,
125  // 1D
126  const std::vector<double>& CredibleIntervals = {0.99, 0.90, 0.68 },
127  const std::vector<Color_t>& CredibleIntervalsColours = {kCyan+4, kCyan-2, kCyan-10},
128  //2D
129  const std::vector<double>& CredibleRegions = {0.99, 0.90, 0.68},
130  const std::vector<Style_t>& CredibleRegionStyle = {kDashed, kSolid, kDotted},
131  const std::vector<Color_t>& CredibleRegionColor = {kGreen-3, kGreen-10, kGreen},
132  // Other
133  const bool CredibleInSigmas = false
134  );
135 
140  void CheckCredibleIntervalsOrder(const std::vector<double>& CredibleIntervals, const std::vector<Color_t>& CredibleIntervalsColours) const;
141 
147  void CheckCredibleRegionsOrder(const std::vector<double>& CredibleRegions,
148  const std::vector<Style_t>& CredibleRegionStyle,
149  const std::vector<Color_t>& CredibleRegionColor);
150 
153  void GetPolarPlot(const std::vector<std::string>& ParNames);
154 
160  void GetBayesFactor(const std::vector<std::string>& ParName,
161  const std::vector<std::vector<double>>& Model1Bounds,
162  const std::vector<std::vector<double>>& Model2Bounds,
163  const std::vector<std::vector<std::string>>& ModelNames);
165  void GetSavageDickey(const std::vector<std::string>& ParName,
166  const std::vector<double>& EvaluationPoint,
167  const std::vector<std::vector<double>>& Bounds);
168 
172  void SavageDickeyPlot(std::unique_ptr<TH1D>& PriorHist,
173  std::unique_ptr<TH1D>& PosteriorHist,
174  const std::string& Title,
175  const double EvaluationPoint) const ;
176 
181  void ReweightPrior(const std::vector<std::string>& Names,
182  const std::vector<double>& NewCentral,
183  const std::vector<double>& NewError);
184 
185 
190  void SmearChain(const std::vector<std::string>& Names,
191  const std::vector<double>& NewCentral,
192  const bool& SaveBranch);
193 
197  void ParameterEvolution(const std::vector<std::string>& Names,
198  const std::vector<int>& NIntervals);
199 
203  inline void ThinMCMC(const int ThinningCut) const { ThinningMCMC(MCMCFile+".root", ThinningCut); };
204 
206  void DiagMCMC();
207 
208  // Get the number of parameters
210  inline int GetNParams() const { return nDraw; };
211  inline int GetNXSec() const { return nParam[kXSecPar]; };
212  inline int GetNND() const { return nParam[kNDPar]; };
213  inline int GetNFD() const { return nParam[kFDDetPar]; };
214 
216  YAML::Node GetCovConfig(const int i) const {return CovConfig.at(i); }
217 
219  int GetGroup(const std::string& name) const;
220 
223  inline TH1D* GetHpost(const int i) const { return hpost[i]; };
227  inline TH2D* GetHpost2D(const int i, const int j) const { return hpost2D[i][j]; };
229  inline TH2D* GetViolin() const { return hviolin.get(); };
231  inline TH2D* GetViolinPrior() const { return hviolin_prior.get(); };
232 
233  //Covariance getters
234  inline std::vector<std::string> GetXSecCov() const { return CovPos[kXSecPar]; };
235  inline std::string GetNDCov() const { return CovPos[kNDPar].back(); };
236  inline std::string GetFDCov() const { return CovPos[kFDDetPar].back(); };
237 
239  void GetPostfit(TVectorD *&Central, TVectorD *&Errors, TVectorD *&Central_Gauss, TVectorD *&Errors_Gauss, TVectorD *&Peaks);
243  void GetCovariance(TMatrixDSym *&Cov, TMatrixDSym *&Corr);
245  void GetPostfit_Ind(TVectorD *&Central, TVectorD *&Errors, TVectorD *&Peaks, ParameterEnum kParam);
246 
248  const std::vector<TString>& GetBranchNames() const { return BranchNames;};
250  void GetNthParameter(const int param, double &Prior, double &PriorError, TString &Title) const;
252  int GetParamIndexFromName(const std::string& Name);
254  inline Long64_t GetnEntries(){return nEntries;};
256  inline Long64_t GetnSteps(){return nSteps;};
258  inline void SetNBins(const int NewBins) {nBins = NewBins;};
259 
262  inline void SetEntries(const int NewEntries) {
263  if (NewEntries > nEntries) {
264  MACH3LOG_ERROR("Cannot increase entries from {} to {}. Only decreasing is allowed.", nEntries, NewEntries);
265  throw MaCh3Exception(__FILE__, __LINE__);
266  }
267  if (NewEntries <= 0) {
268  MACH3LOG_ERROR("Entries cannot be below 0, but {} was passed.", NewEntries);
269  throw MaCh3Exception(__FILE__, __LINE__);
270  }
271 
272  if (static_cast<int>(BurnInCut) > NewEntries) {
273  MACH3LOG_ERROR("BurnInCut ({}) is larger than NewEntries ({})", BurnInCut, NewEntries);
274  throw MaCh3Exception(__FILE__, __LINE__);
275  }
276 
277  MACH3LOG_INFO("Setting entries to {} from {}.", NewEntries, nEntries);
278  MACH3LOG_WARN("This may behave not as expected when using merged multiple chains");
279  nEntries = NewEntries;
280  }
283  void SetStepCut(const std::string& Cuts);
286  void SetStepCut(const int Cuts);
288  void CheckStepCut() const;
289 
292  inline void SetPlotRelativeToPrior(const bool PlotOrNot){plotRelativeToPrior = PlotOrNot; };
293  inline void SetPrintToPDF(const bool PlotOrNot){printToPDF = PlotOrNot; };
295  inline void SetPlotErrorForFlatPrior(const bool PlotOrNot){PlotFlatPrior = PlotOrNot; };
296  inline void SetPlotBinValue(const bool PlotOrNot){plotBinValue = PlotOrNot; };
297  inline void SetFancyNames(const bool PlotOrNot){FancyPlotNames = PlotOrNot; };
299  inline void SetSmoothing(const bool PlotOrNot){ApplySmoothing = PlotOrNot; };
302  inline void SetPost2DPlotThreshold(const double Threshold){Post2DPlotThreshold = Threshold; };
304  inline void SetUseFFTAutoCorrelation(const bool useFFT){useFFTAutoCorrelation = useFFT; };
305 
308  inline void SetExcludedTypes(std::vector<std::string> Name){ExcludedTypes = Name; };
309  inline void SetExcludedNames(std::vector<std::string> Name){ExcludedNames = Name; };
310  inline void SetExcludedGroups(std::vector<std::string> Name){ExcludedGroups = Name; };
311 
314  inline void SetnBatches(const int Batches){nBatches = Batches; };
315  inline void SetnLags(const int nLags){AutoCorrLag = nLags; };
316 
318  inline void SetOutputSuffix(const std::string Suffix){OutputSuffix = Suffix; };
320  inline void SetPosterior1DCut(const std::string Cut){Posterior1DCut = Cut; };
321  protected:
323  inline std::unique_ptr<TH1D> MakePrefit();
325  inline void MakeOutputFile();
327  inline void DrawCorrelations1D();
331  inline void DrawCorrelationsGroup(const std::unique_ptr<TH2D>& CorrMatrix) const;
333  inline void ReadInputCov();
335  inline void ReadInputCovLegacy();
337  inline void FindInputFiles();
339  inline void FindInputFilesLegacy();
341  inline void ReadModelFile();
343  virtual void LoadAdditionalInfo() {};
345  inline void ReadNDFile();
347  inline void ReadFDFile();
349  inline void PrintInfo() const;
350 
352  inline void ScanInput();
354  inline void ScanParameterOrder();
356  inline void SetupOutput();
357 
358  // MCMC Diagnostic
360  inline void PrepareDiagMCMC();
362  inline void ParamTraces();
364  inline void AutoCorrelation();
367  inline void AutoCorrelation_FFT();
378  inline void CalculateESS(const int nLags, const std::vector<std::vector<double>>& LagL);
382  inline void BatchedAnalysis();
384  inline void BatchedMeans();
388  inline void GewekeDiagnostic();
390  inline void AcceptanceProbabilities();
394  inline void PowerSpectrumAnalysis();
395 
397  std::vector<double> GetMargins(const std::unique_ptr<TCanvas>& Canv) const;
399  void SetMargins(std::unique_ptr<TCanvas>& Canv, const std::vector<double>& margins);
405  void SetTLineStyle(TLine* Line, const Color_t Colour, const Width_t Width, const ELineStyle Style) const;
409  void SetLegendStyle(TLegend* Legend, const double size) const;
410 
412  std::string MCMCFile;
414  std::string OutputSuffix;
416  std::vector<std::vector<std::string>> CovPos;
418  std::vector<std::string> CovNamePos;
420  std::vector<YAML::Node> CovConfig;
421 
423  TChain *Chain;
425  std::string StepCut;
427  std::string Posterior1DCut;
429  unsigned int UpperCut;
431  unsigned int BurnInCut;
435  int nEntries;
437  int nSteps;
439  int nSamples;
441  int nSysts;
443  int nDraw;
444 
445  //Name of all branches as well as branches we don't want to include in the analysis
446  std::vector<TString> BranchNames;
447  std::vector<std::string> ExcludedTypes;
448  std::vector<std::string> ExcludedNames;
449  std::vector<std::string> ExcludedGroups;
450 
452  std::vector<bool> IamVaried;
454  std::vector<std::vector<TString>> ParamNames;
456  std::vector<std::vector<double>> ParamCentral;
457  std::vector<std::vector<double>> ParamNom;
459  std::vector<std::vector<double>> ParamErrors;
461  std::vector<std::vector<bool>> ParamFlat;
463  std::vector<int> nParam;
465  std::vector<ParameterEnum> ParamType;
468  std::vector<int> ParamTypeStartPos;
469  // KS: For example flux or detector within matrix
470  std::vector<std::string> ParameterGroup;
471 
473  std::vector<TString> SampleName_v;
475  std::vector<TString> SystName_v;
476 
478  std::string OutputName;
480  TString CanvasName;
481 
484 
485  //Even more flags
502 
503  std::vector<int> NDSamplesBins;
504  std::vector<std::string> NDSamplesNames;
505 
507  std::unique_ptr<TF1> Gauss;
508 
510  TFile *OutputFile;
511 
513  std::unique_ptr<TCanvas> Posterior;
514 
515  //Vector of best fit points and errors obtained with different methods
517  TVectorD *Central_Value;
519  TVectorD *Means;
521  TVectorD *Errors;
523  TVectorD *Means_Gauss;
525  TVectorD *Errors_Gauss;
527  TVectorD *Means_HPD;
529  TVectorD *Errors_HPD;
534 
536  TMatrixDSym *Covariance;
538  TMatrixDSym *Correlation;
539 
541  std::vector<TH1D*> hpost;
543  std::vector<std::vector<TH2D*>> hpost2D;
545  std::unique_ptr<TH2D> hviolin;
547  std::unique_ptr<TH2D> hviolin_prior;
548 
550  double** ParStep;
552  unsigned int* StepNumber;
553 
555  int nBins;
557  double DrawRange;
558 
560  bool CacheMCMC;
563 
564  //Number of batches and LagL used in MCMC diagnostic
566  int nBatches;
569 
571  double *ParamSums;
573  double **BatchedAverages;
574 
576  double **SampleValues;
578  double **SystValues;
579 
581  double *AccProbValues;
584 
588  std::string ReweightName;
590  double* WeightValue;
591 
592  //Only if GPU is enabled
593  #ifdef MaCh3_CUDA
595  inline void PrepareGPU_AutoCorr(const int nLags);
596 
598  float* ParStep_cpu;
599  float* NumeratorSum_cpu;
600  float* ParamSums_cpu;
601  float* DenomSum_cpu;
602 
604  float* ParStep_gpu;
605  float* NumeratorSum_gpu;
606  float* ParamSums_gpu;
607  float* DenomSum_gpu;
608  #endif
609 };
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:109
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:120
int size
ParameterEnum
KS: Enum for different covariance classes.
Definition: MCMCProcessor.h:48
@ kNDPar
Definition: MCMCProcessor.h:50
@ kXSecPar
Definition: MCMCProcessor.h:49
@ kNParameterEnum
Definition: MCMCProcessor.h:53
@ kFDDetPar
Definition: MCMCProcessor.h:51
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:25
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:26
constexpr ELineStyle Style[NVars]
void ThinningMCMC(const std::string &FilePath, const int ThinningCut)
Thin MCMC Chain, to save space and maintain low autocorrelations.
Utility functions for statistical interpretations in MaCh3.
Class responsible for processing MCMC chains, performing diagnostics, generating plots,...
Definition: MCMCProcessor.h:61
void ScanParameterOrder()
Scan order of params from a different groups.
int nBatches
Number of batches for Batched Mean.
void CheckStepCut() const
Check if step cut isn't larger than highest values of step in a chain.
void GewekeDiagnostic()
Geweke Diagnostic based on the methods described by Fang (2014) and Karlsbakk (2011)....
TMatrixDSym * Correlation
Posterior Correlation Matrix.
void MakeViolin()
Make and Draw Violin.
void GetNthParameter(const int param, double &Prior, double &PriorError, TString &Title) const
Get properties of parameter by passing it number.
void SetExcludedGroups(std::vector< std::string > Name)
double ** ParStep
Array holding values for all parameters.
void SetExcludedNames(std::vector< std::string > Name)
void ThinMCMC(const int ThinningCut) const
Thin MCMC Chain, to save space and maintain low autocorrelations.
const std::vector< TString > & GetBranchNames() const
Get the vector of branch names from root file.
void PrintInfo() const
Print info like how many params have been loaded etc.
void SetnLags(const int nLags)
void MakeCredibleIntervals(const std::vector< double > &CredibleIntervals={0.99, 0.90, 0.68 }, const std::vector< Color_t > &CredibleIntervalsColours={kCyan+4, kCyan-2, kCyan-10}, const bool CredibleInSigmas=false)
Make and Draw Credible intervals.
void ReadModelFile()
Read the xsec file and get the input central values and errors.
std::vector< double > GetMargins(const std::unique_ptr< TCanvas > &Canv) const
Get TCanvas margins, to be able to reset them if particular function need different margins.
std::unique_ptr< TF1 > Gauss
Gaussian fitter.
void Initialise()
Scan chain, what parameters we have and load information from covariance matrices.
MCMCProcessor(const std::string &InputFile)
Constructs an MCMCProcessor object with the specified input file and options.
double Post2DPlotThreshold
KS: Set Threshold when to plot 2D posterior as by default we get a LOT of plots.
void AcceptanceProbabilities()
Acceptance Probability.
double ** SampleValues
Holds the sample values.
void BatchedAnalysis()
Get the batched means variance estimation and variable indicating if number of batches is sensible .
void SetOutputSuffix(const std::string Suffix)
Sett output suffix, this way jobs using the same file will have different names.
void AutoCorrelation()
KS: Calculate autocorrelations supports both OpenMP and CUDA :)
TVectorD * Means_HPD
Vector with mean values using Highest Posterior Density.
void ResetHistograms()
Reset 2D posteriors, in case we would like to calculate in again with different BurnInCut.
std::vector< TString > SampleName_v
Vector of each systematic.
void GetBayesFactor(const std::vector< std::string > &ParName, const std::vector< std::vector< double >> &Model1Bounds, const std::vector< std::vector< double >> &Model2Bounds, const std::vector< std::vector< std::string >> &ModelNames)
Calculate Bayes factor for vector of params, and model boundaries.
double * WeightValue
Stores value of weight for each step.
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.
void MakeTrianglePlot(const std::vector< std::string > &ParNames, const std::vector< double > &CredibleIntervals={0.99, 0.90, 0.68 }, const std::vector< Color_t > &CredibleIntervalsColours={kCyan+4, kCyan-2, kCyan-10}, const std::vector< double > &CredibleRegions={0.99, 0.90, 0.68}, const std::vector< Style_t > &CredibleRegionStyle={kDashed, kSolid, kDotted}, const std::vector< Color_t > &CredibleRegionColor={kGreen-3, kGreen-10, kGreen}, const bool CredibleInSigmas=false)
Make fancy triangle plot for selected parameters.
std::string OutputName
Name of output files.
std::unique_ptr< TH2D > hviolin_prior
Holds prior violin plot for all dials,.
int nSamples
Number of sample PDF objects.
std::vector< int > nParam
Number of parameters per type.
void SetExcludedTypes(std::vector< std::string > Name)
Setter related what parameters we want to exclude from analysis, for example if cross-section paramet...
int GetGroup(const std::string &name) const
Number of params from a given group, for example flux.
std::vector< std::string > CovNamePos
Covariance matrix name position.
double DrawRange
Drawrange for SetMaximum.
std::vector< std::vector< bool > > ParamFlat
Whether Param has flat prior or not.
void SetSmoothing(const bool PlotOrNot)
Set whether want to use smoothing for histograms using ROOT algorithm.
std::vector< YAML::Node > CovConfig
Covariance matrix config.
std::vector< std::vector< double > > ParamNom
double ** SystValues
Holds the systs values.
virtual ~MCMCProcessor()
Destroys the MCMCProcessor object.
TVectorD * Errors
Vector with errors values using RMS.
double * AccProbBatchedAverages
Holds all accProb in batches.
void MakePostfit(const std::map< std::string, std::pair< double, double >> &Edges={})
Make 1D projection for each parameter and prepare structure.
bool useFFTAutoCorrelation
MJR: Use FFT-based autocorrelation algorithm (save time & resources)?
void MakeCredibleRegions(const std::vector< double > &CredibleRegions={0.99, 0.90, 0.68}, const std::vector< Style_t > &CredibleRegionStyle={kDashed, kSolid, kDotted}, const std::vector< Color_t > &CredibleRegionColor={kGreen-3, kGreen-10, kGreen}, const bool CredibleInSigmas=false, const bool Draw2DPosterior=true, const bool DrawBestFit=true)
Make and Draw Credible Regions.
std::string StepCut
BurnIn Cuts.
void GetPostfit_Ind(TVectorD *&Central, TVectorD *&Errors, TVectorD *&Peaks, ParameterEnum kParam)
Or the individual post-fits.
int GetParamIndexFromName(const std::string &Name)
Get parameter number based on name.
void DrawCorrelations1D()
Draw 1D correlations which might be more helpful than looking at huge 2D Corr matrix.
void SetnBatches(const int Batches)
Set value of Nbatches used for batched mean, this need to be done earlier as batches are made when re...
int GetNFD() const
void SetEntries(const int NewEntries)
Set number of entries to make potentially MCMC Processing faster.
void GetPostfit(TVectorD *&Central, TVectorD *&Errors, TVectorD *&Central_Gauss, TVectorD *&Errors_Gauss, TVectorD *&Peaks)
Get the post-fit results (arithmetic and Gaussian)
TVectorD * Means_Gauss
Vector with mean values using Gaussian fit.
unsigned int * StepNumber
Step number for step, important if chains were merged.
int AutoCorrLag
LagL used in AutoCorrelation.
std::unique_ptr< TCanvas > Posterior
Fancy canvas used for our beautiful plots.
TFile * OutputFile
The output file.
void SetPosterior1DCut(const std::string Cut)
Allow to set addtional cuts based on ROOT TBrowser cut, for to only affect one mass ordering.
void ReadNDFile()
Read the ND cov file and get the input central values and errors.
TH1D * GetHpost(const int i) const
Get 1D posterior for a given parameter.
bool ApplySmoothing
Apply smoothing for 2D histos using root algorithm.
unsigned int UpperCut
KS: Used only for SubOptimality.
int nBins
Number of bins.
TChain * Chain
Main chain storing all steps etc.
std::string MCMCFile
Name of MCMC file.
void SetPlotRelativeToPrior(const bool PlotOrNot)
You can set relative to prior or relative to generated. It is advised to use relate to prior.
bool ReweightPosterior
Whether to apply reweighting weight or not.
TH2D * GetViolin() const
Get Violin plot for all parameters with posterior values.
void SetLegendStyle(TLegend *Legend, const double size) const
Configures the style of a TLegend object.
int GetNParams() const
Get total number of used parameters.
std::vector< std::string > ExcludedNames
std::unique_ptr< TH1D > MakePrefit()
Prepare prefit histogram for parameter overlay plot.
std::vector< TH1D * > hpost
Holds 1D Posterior Distributions.
YAML::Node GetCovConfig(const int i) const
Get Yaml config obtained from a Chain.
TVectorD * Errors_HPD_Negative
Vector with negative error (left hand side) values using Highest Posterior Density.
std::vector< std::vector< std::string > > CovPos
Covariance matrix file name position.
void DrawCorrelationsGroup(const std::unique_ptr< TH2D > &CorrMatrix) const
Produces correlation matrix but instead of giving name for each param it only give name for param gro...
double * ParamSums
Total parameter sum for each param.
Long64_t GetnEntries()
Get Number of entries that Chain has, for merged chains will not be the same Nsteps.
void SetFancyNames(const bool PlotOrNot)
void DiagMCMC()
KS: Perform MCMC diagnostic including Autocorrelation, Trace etc.
std::vector< std::string > ExcludedGroups
std::string Posterior1DCut
Cut used when making 1D Posterior distribution.
double * AccProbValues
Holds all accProb.
void FindInputFilesLegacy()
void DrawCovariance()
Draw the post-fit covariances.
void SetMargins(std::unique_ptr< TCanvas > &Canv, const std::vector< double > &margins)
Set TCanvas margins to specified values.
void ParamTraces()
CW: Draw trace plots of the parameters i.e. parameter vs step.
void PrepareDiagMCMC()
CW: Prepare branches etc. for DiagMCMC.
std::vector< bool > IamVaried
Is the ith parameter varied.
TH2D * GetHpost2D(const int i, const int j) const
Get 2D posterior for a given parameter combination.
std::unique_ptr< TH2D > hviolin
Holds violin plot for all dials.
void SetPrintToPDF(const bool PlotOrNot)
Long64_t GetnSteps()
Get Number of Steps that Chain has, for merged chains will not be the same nEntries.
int nDraw
Number of all parameters used in the analysis.
std::string OutputSuffix
Output file suffix useful when running over same file with different settings.
int GetNXSec() const
void SetNBins(const int NewBins)
Modify number of bins used for 1D and 2D Histograms.
void SetupOutput()
Prepare all objects used for output.
void MakeOutputFile()
prepare output root file and canvas to which we will save EVERYTHING
bool plotBinValue
If true it will print value on each bin of covariance matrix.
TVectorD * Errors_Gauss
Vector with error values using Gaussian fit.
void CheckCredibleIntervalsOrder(const std::vector< double > &CredibleIntervals, const std::vector< Color_t > &CredibleIntervalsColours) const
Checks the order and size consistency of the CredibleIntervals and CredibleIntervalsColours vectors.
void CalculateESS(const int nLags, const std::vector< std::vector< double >> &LagL)
KS: calc Effective Sample Size.
std::vector< ParameterEnum > ParamType
Make an enum for which class this parameter belongs to so we don't have to keep string comparing.
std::vector< std::string > NDSamplesNames
virtual void LoadAdditionalInfo()
allow loading additional info for example used for oscillation parameters
TVectorD * Central_Value
Vector with central value for each parameter.
std::vector< std::string > ExcludedTypes
int nSteps
KS: For merged chains number of entries will be different from nSteps.
TH2D * GetViolinPrior() const
Get Violin plot for all parameters with prior values.
void CheckCredibleRegionsOrder(const std::vector< double > &CredibleRegions, const std::vector< Style_t > &CredibleRegionStyle, const std::vector< Color_t > &CredibleRegionColor)
Checks the order and size consistency of the CredibleRegions, CredibleRegionStyle,...
std::vector< int > NDSamplesBins
void MakeCovariance_MP(const bool Mute=false)
Calculate covariance by making 2D projection of each combination of parameters using multithreading.
double ** BatchedAverages
Values of batched average for every param and batch.
void DrawPostfit()
Draw the post-fit comparisons.
TString CanvasName
Name of canvas which help to save to the sample pdf.
void SetPost2DPlotThreshold(const double Threshold)
Code will only plot 2D posteriors if Correlation are larger than defined threshold.
std::vector< std::string > ParameterGroup
std::string GetFDCov() const
void SetUseFFTAutoCorrelation(const bool useFFT)
Toggle using the FFT-based autocorrelation calculator.
TVectorD * Errors_HPD
Vector with error values using Highest Posterior Density.
void AutoCorrelation_FFT()
MJR: Autocorrelation function using FFT algorithm for extra speed.
void SetTLineStyle(TLine *Line, const Color_t Colour, const Width_t Width, const ELineStyle Style) const
Configures a TLine object with the specified style parameters.
void ReadInputCovLegacy()
void GetPolarPlot(const std::vector< std::string > &ParNames)
Make funny polar plot.
std::vector< std::vector< TH2D * > > hpost2D
Holds 2D Posterior Distributions.
void ParameterEvolution(const std::vector< std::string > &Names, const std::vector< int > &NIntervals)
Make .gif of parameter evolution.
void GetCovariance(TMatrixDSym *&Cov, TMatrixDSym *&Corr)
Get the post-fit covariances and correlations.
TVectorD * Means
Vector with mean values using Arithmetic Mean.
std::vector< TString > SystName_v
Vector of each sample PDF object.
void CacheSteps()
KS:By caching each step we use multithreading.
std::vector< TString > BranchNames
bool PlotFlatPrior
Whether we plot flat prior or not, we usually provide error even for flat prior params.
bool FancyPlotNames
Whether we want fancy plot names or not.
std::string GetNDCov() const
void ReweightPrior(const std::vector< std::string > &Names, const std::vector< double > &NewCentral, const std::vector< double > &NewError)
Reweight Prior by giving new central value and new error.
std::string ReweightName
Name of branch used for chain reweighting.
void SetStepCut(const std::string &Cuts)
Set the step cutting by string.
bool printToPDF
Will plot all plot to PDF not only to root file.
void SmearChain(const std::vector< std::string > &Names, const std::vector< double > &NewCentral, const bool &SaveBranch)
Smear chain contours.
int GetNND() const
void GetSavageDickey(const std::vector< std::string > &ParName, const std::vector< double > &EvaluationPoint, const std::vector< std::vector< double >> &Bounds)
Calculate Bayes factor for point like hypothesis using SavageDickey.
void SetPlotErrorForFlatPrior(const bool PlotOrNot)
Set whether you want to plot error for parameters which have flat prior.
std::vector< std::vector< TString > > ParamNames
Name of parameters which we are going to analyse.
bool doDiagMCMC
Doing MCMC Diagnostic.
int nSysts
Number of covariance objects.
void ReadFDFile()
Read the FD cov file and get the input central values and errors.
std::vector< std::string > GetXSecCov() const
void FindInputFiles()
Read the output MCMC file and find what inputs were used.
bool CacheMCMC
MCMC Chain has been cached.
std::vector< int > ParamTypeStartPos
bool MadePostfit
Sanity check if Postfit is already done to not make several times.
void SetPlotBinValue(const bool PlotOrNot)
void PowerSpectrumAnalysis()
RC: Perform spectral analysis of MCMC .
void ScanInput()
Scan Input etc.
void BatchedMeans()
CW: Batched means, literally read from an array and chuck into TH1D.
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.
int nBranches
Number of branches in a TTree.
TVectorD * Errors_HPD_Positive
Vector with positive error (right hand side) values using Highest Posterior Density.
TMatrixDSym * Covariance
Posterior Covariance Matrix.
void MakeSubOptimality(const int NIntervals=10)
Make and Draw SubOptimality .
void MakeCovarianceYAML(const std::string &OutputYAMLFile, const std::string &MeansMethod) const
Make YAML file from post-fit covariance.
bool plotRelativeToPrior
Whether we plot relative to prior or nominal, in most cases is prior.
void MakeCovariance()
Calculate covariance by making 2D projection of each combination of parameters.
unsigned int BurnInCut
Value of burn in cut.
void ReadInputCov()
CW: Read the input Covariance matrix entries. Get stuff like parameter input errors,...
Custom exception class for MaCh3 errors.