MaCh3  2.2.3
Reference Guide
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
OscProcessor Class Reference

This class extends MCMC and allow specialised for Oscillation parameters analysis which require specialised hardcoding. More...

#include <Fitters/OscProcessor.h>

Inheritance diagram for OscProcessor:
[legend]
Collaboration diagram for OscProcessor:
[legend]

Public Member Functions

 OscProcessor (const std::string &InputFile)
 Constructs an OscProcessor object with the specified input file and options. More...
 
virtual ~OscProcessor ()
 Destroys the OscProcessor object. More...
 
void PerformJarlskogAnalysis ()
 Perform Several Jarlskog Plotting. More...
 
- Public Member Functions inherited from MCMCProcessor
 MCMCProcessor (const std::string &InputFile)
 Constructs an MCMCProcessor object with the specified input file and options. More...
 
virtual ~MCMCProcessor ()
 Destroys the MCMCProcessor object. More...
 
void Initialise ()
 Scan chain, what parameters we have and load information from covariance matrices. More...
 
void MakePostfit (const std::map< std::string, std::pair< double, double >> &Edges={})
 Make 1D projection for each parameter and prepare structure. More...
 
void MakeCovariance ()
 Calculate covariance by making 2D projection of each combination of parameters. More...
 
void CacheSteps ()
 KS:By caching each step we use multithreading. More...
 
void MakeCovariance_MP (const bool Mute=false)
 Calculate covariance by making 2D projection of each combination of parameters using multithreading. More...
 
void MakeSubOptimality (const int NIntervals=10)
 Make and Draw SubOptimality [27]. More...
 
void ResetHistograms ()
 Reset 2D posteriors, in case we would like to calculate in again with different BurnInCut. More...
 
void DrawPostfit ()
 Draw the post-fit comparisons. More...
 
void MakeViolin ()
 Make and Draw Violin. More...
 
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. More...
 
void DrawCovariance ()
 Draw the post-fit covariances. More...
 
void MakeCovarianceYAML (const std::string &OutputYAMLFile, const std::string &MeansMethod) const
 Make YAML file from post-fit covariance. More...
 
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. More...
 
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. More...
 
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. More...
 
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, and CredibleRegionColor vectors. More...
 
void GetPolarPlot (const std::vector< std::string > &ParNames)
 Make funny polar plot. More...
 
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. More...
 
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. More...
 
void SavageDickeyPlot (std::unique_ptr< TH1D > &PriorHist, std::unique_ptr< TH1D > &PosteriorHist, const std::string &Title, const double EvaluationPoint) const
 Produce Savage Dickey plot. More...
 
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. More...
 
void SmearChain (const std::vector< std::string > &Names, const std::vector< double > &NewCentral, const bool &SaveBranch)
 Smear chain contours. More...
 
void ParameterEvolution (const std::vector< std::string > &Names, const std::vector< int > &NIntervals)
 Make .gif of parameter evolution. More...
 
void ThinMCMC (const int ThinningCut) const
 Thin MCMC Chain, to save space and maintain low autocorrelations. More...
 
void DiagMCMC ()
 KS: Perform MCMC diagnostic including Autocorrelation, Trace etc. More...
 
int GetNParams () const
 Get total number of used parameters. More...
 
int GetNXSec () const
 
int GetNND () const
 
int GetNFD () const
 
YAML::Node GetCovConfig (const int i) const
 Get Yaml config obtained from a Chain. More...
 
int GetGroup (const std::string &name) const
 Number of params from a given group, for example flux. More...
 
TH1D * GetHpost (const int i) const
 Get 1D posterior for a given parameter. More...
 
TH2D * GetHpost2D (const int i, const int j) const
 Get 2D posterior for a given parameter combination. More...
 
TH2D * GetViolin () const
 Get Violin plot for all parameters with posterior values. More...
 
TH2D * GetViolinPrior () const
 Get Violin plot for all parameters with prior values. More...
 
std::vector< std::string > GetXSecCov () const
 
std::string GetNDCov () const
 
std::string GetFDCov () const
 
void GetPostfit (TVectorD *&Central, TVectorD *&Errors, TVectorD *&Central_Gauss, TVectorD *&Errors_Gauss, TVectorD *&Peaks)
 Get the post-fit results (arithmetic and Gaussian) More...
 
void GetCovariance (TMatrixDSym *&Cov, TMatrixDSym *&Corr)
 Get the post-fit covariances and correlations. More...
 
void GetPostfit_Ind (TVectorD *&Central, TVectorD *&Errors, TVectorD *&Peaks, ParameterEnum kParam)
 Or the individual post-fits. More...
 
const std::vector< TString > & GetBranchNames () const
 Get the vector of branch names from root file. More...
 
void GetNthParameter (const int param, double &Prior, double &PriorError, TString &Title) const
 Get properties of parameter by passing it number. More...
 
int GetParamIndexFromName (const std::string &Name)
 Get parameter number based on name. More...
 
Long64_t GetnEntries ()
 Get Number of entries that Chain has, for merged chains will not be the same Nsteps. More...
 
Long64_t GetnSteps ()
 Get Number of Steps that Chain has, for merged chains will not be the same nEntries. More...
 
void SetNBins (const int NewBins)
 Modify number of bins used for 1D and 2D Histograms. More...
 
void SetEntries (const int NewEntries)
 Set number of entries to make potentially MCMC Processing faster. More...
 
void SetStepCut (const std::string &Cuts)
 Set the step cutting by string. More...
 
void SetStepCut (const int Cuts)
 Set the step cutting by int. More...
 
void CheckStepCut () const
 Check if step cut isn't larger than highest values of step in a chain. More...
 
void SetPlotRelativeToPrior (const bool PlotOrNot)
 You can set relative to prior or relative to generated. It is advised to use relate to prior. More...
 
void SetPrintToPDF (const bool PlotOrNot)
 
void SetPlotErrorForFlatPrior (const bool PlotOrNot)
 Set whether you want to plot error for parameters which have flat prior. More...
 
void SetPlotBinValue (const bool PlotOrNot)
 
void SetFancyNames (const bool PlotOrNot)
 
void SetSmoothing (const bool PlotOrNot)
 Set whether want to use smoothing for histograms using ROOT algorithm. More...
 
void SetPost2DPlotThreshold (const double Threshold)
 Code will only plot 2D posteriors if Correlation are larger than defined threshold. More...
 
void SetUseFFTAutoCorrelation (const bool useFFT)
 Toggle using the FFT-based autocorrelation calculator. More...
 
void SetExcludedTypes (std::vector< std::string > Name)
 Setter related what parameters we want to exclude from analysis, for example if cross-section parameters look like xsec_, then passing "xsec_" will. More...
 
void SetExcludedNames (std::vector< std::string > Name)
 
void SetExcludedGroups (std::vector< std::string > Name)
 
void SetnBatches (const int Batches)
 Set value of Nbatches used for batched mean, this need to be done earlier as batches are made when reading tree. More...
 
void SetnLags (const int nLags)
 
void SetOutputSuffix (const std::string Suffix)
 Sett output suffix, this way jobs using the same file will have different names. More...
 
void SetPosterior1DCut (const std::string Cut)
 Allow to set addtional cuts based on ROOT TBrowser cut, for to only affect one mass ordering. More...
 

Protected Member Functions

void LoadAdditionalInfo () override
 Read the Osc cov file and get the input central values and errors Here we allow Jarlskog Shenanigans. More...
 
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. More...
 
double CalcJarlskog (const double s2th13, const double s2th23, const double s2th12, const double dcp) const
 Calculate Jarlskog Invariant using oscillation parameters. More...
 
double SamplePriorForParam (const int paramIndex, const std::unique_ptr< TRandom3 > &randGen, const std::vector< double > &FlatBounds) const
 Draw Prior value. More...
 
- Protected Member Functions inherited from MCMCProcessor
std::unique_ptr< TH1D > MakePrefit ()
 Prepare prefit histogram for parameter overlay plot. More...
 
void MakeOutputFile ()
 prepare output root file and canvas to which we will save EVERYTHING More...
 
void DrawCorrelations1D ()
 Draw 1D correlations which might be more helpful than looking at huge 2D Corr matrix. More...
 
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 group. More...
 
void ReadInputCov ()
 CW: Read the input Covariance matrix entries. Get stuff like parameter input errors, names, and so on. More...
 
void ReadInputCovLegacy ()
 
void FindInputFiles ()
 Read the output MCMC file and find what inputs were used. More...
 
void FindInputFilesLegacy ()
 
void ReadModelFile ()
 Read the xsec file and get the input central values and errors. More...
 
void ReadNDFile ()
 Read the ND cov file and get the input central values and errors. More...
 
void ReadFDFile ()
 Read the FD cov file and get the input central values and errors. More...
 
void PrintInfo () const
 Print info like how many params have been loaded etc. More...
 
void ScanInput ()
 Scan Input etc. More...
 
void ScanParameterOrder ()
 Scan order of params from a different groups. More...
 
void SetupOutput ()
 Prepare all objects used for output. More...
 
void PrepareDiagMCMC ()
 CW: Prepare branches etc. for DiagMCMC. More...
 
void ParamTraces ()
 CW: Draw trace plots of the parameters i.e. parameter vs step. More...
 
void AutoCorrelation ()
 KS: Calculate autocorrelations supports both OpenMP and CUDA :) More...
 
void AutoCorrelation_FFT ()
 MJR: Autocorrelation function using FFT algorithm for extra speed. More...
 
void CalculateESS (const int nLags, const std::vector< std::vector< double >> &LagL)
 KS: calc Effective Sample Size. More...
 
void BatchedAnalysis ()
 Get the batched means variance estimation and variable indicating if number of batches is sensible [4] [28]. More...
 
void BatchedMeans ()
 CW: Batched means, literally read from an array and chuck into TH1D. More...
 
void GewekeDiagnostic ()
 Geweke Diagnostic based on the methods described by Fang (2014) and Karlsbakk (2011). [8] [20]. More...
 
void AcceptanceProbabilities ()
 Acceptance Probability. More...
 
void PowerSpectrumAnalysis ()
 RC: Perform spectral analysis of MCMC [7]. More...
 
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. More...
 
void SetMargins (std::unique_ptr< TCanvas > &Canv, const std::vector< double > &margins)
 Set TCanvas margins to specified values. More...
 
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. More...
 
void SetLegendStyle (TLegend *Legend, const double size) const
 Configures the style of a TLegend object. More...
 

Protected Attributes

bool PlotJarlskog
 Will plot Jarlskog Invariant using information in the chain. More...
 
bool OscEnabled
 Will plot Jarlskog Invariant using information in the chain. More...
 
std::string Sin2Theta13Name
 Name of the parameter representing \(\sin^2\theta_{13}\). More...
 
std::string Sin2Theta12Name
 Name of the parameter representing \(\sin^2\theta_{12}\). More...
 
std::string Sin2Theta23Name
 Name of the parameter representing \(\sin^2\theta_{23}\). More...
 
std::string DeltaCPName
 Name of the parameter representing \(\delta_{\mathrm{CP}}\) (the CP-violating phase). More...
 
std::string DeltaM2_23Name
 Name of the parameter representing \(\Delta m^2_{32}\) (mass-squared difference). More...
 
int Sin2Theta13Index
 Index of \(\sin^2\theta_{13}\) in the parameter list. More...
 
int Sin2Theta12Index
 Index of \(\sin^2\theta_{12}\) in the parameter list. More...
 
int Sin2Theta23Index
 Index of \(\sin^2\theta_{23}\) in the parameter list. More...
 
int DeltaCPIndex
 Index of \(\delta_{\mathrm{CP}}\) in the parameter list. More...
 
int DeltaM2_23Index
 Index of \(\Delta m^2_{32}\) in the parameter list. More...
 
- Protected Attributes inherited from MCMCProcessor
std::string MCMCFile
 Name of MCMC file. More...
 
std::string OutputSuffix
 Output file suffix useful when running over same file with different settings. More...
 
std::vector< std::vector< std::string > > CovPos
 Covariance matrix file name position. More...
 
std::vector< std::string > CovNamePos
 Covariance matrix name position. More...
 
std::vector< YAML::Node > CovConfig
 Covariance matrix config. More...
 
TChain * Chain
 Main chain storing all steps etc. More...
 
std::string StepCut
 BurnIn Cuts. More...
 
std::string Posterior1DCut
 Cut used when making 1D Posterior distribution. More...
 
unsigned int UpperCut
 KS: Used only for SubOptimality. More...
 
unsigned int BurnInCut
 Value of burn in cut. More...
 
int nBranches
 Number of branches in a TTree. More...
 
int nEntries
 KS: For merged chains number of entries will be different from nSteps. More...
 
int nSteps
 KS: For merged chains number of entries will be different from nSteps. More...
 
int nSamples
 Number of sample PDF objects. More...
 
int nSysts
 Number of covariance objects. More...
 
int nDraw
 Number of all parameters used in the analysis. More...
 
std::vector< TString > BranchNames
 
std::vector< std::string > ExcludedTypes
 
std::vector< std::string > ExcludedNames
 
std::vector< std::string > ExcludedGroups
 
std::vector< bool > IamVaried
 Is the ith parameter varied. More...
 
std::vector< std::vector< TString > > ParamNames
 Name of parameters which we are going to analyse. More...
 
std::vector< std::vector< double > > ParamCentral
 Parameters central values which we are going to analyse. More...
 
std::vector< std::vector< double > > ParamNom
 
std::vector< std::vector< double > > ParamErrors
 Uncertainty on a single parameter. More...
 
std::vector< std::vector< bool > > ParamFlat
 Whether Param has flat prior or not. More...
 
std::vector< int > nParam
 Number of parameters per type. More...
 
std::vector< ParameterEnumParamType
 Make an enum for which class this parameter belongs to so we don't have to keep string comparing. More...
 
std::vector< int > ParamTypeStartPos
 
std::vector< std::string > ParameterGroup
 
std::vector< TString > SampleName_v
 Vector of each systematic. More...
 
std::vector< TString > SystName_v
 Vector of each sample PDF object. More...
 
std::string OutputName
 Name of output files. More...
 
TString CanvasName
 Name of canvas which help to save to the sample pdf. More...
 
bool PlotFlatPrior
 Whether we plot flat prior or not, we usually provide error even for flat prior params. More...
 
bool plotRelativeToPrior
 Whether we plot relative to prior or nominal, in most cases is prior. More...
 
bool MadePostfit
 Sanity check if Postfit is already done to not make several times. More...
 
bool printToPDF
 Will plot all plot to PDF not only to root file. More...
 
bool FancyPlotNames
 Whether we want fancy plot names or not. More...
 
bool plotBinValue
 If true it will print value on each bin of covariance matrix. More...
 
bool ApplySmoothing
 Apply smoothing for 2D histos using root algorithm. More...
 
double Post2DPlotThreshold
 KS: Set Threshold when to plot 2D posterior as by default we get a LOT of plots. More...
 
bool useFFTAutoCorrelation
 MJR: Use FFT-based autocorrelation algorithm (save time & resources)? More...
 
std::vector< int > NDSamplesBins
 
std::vector< std::string > NDSamplesNames
 
std::unique_ptr< TF1 > Gauss
 Gaussian fitter. More...
 
TFile * OutputFile
 The output file. More...
 
std::unique_ptr< TCanvas > Posterior
 Fancy canvas used for our beautiful plots. More...
 
TVectorD * Central_Value
 Vector with central value for each parameter. More...
 
TVectorD * Means
 Vector with mean values using Arithmetic Mean. More...
 
TVectorD * Errors
 Vector with errors values using RMS. More...
 
TVectorD * Means_Gauss
 Vector with mean values using Gaussian fit. More...
 
TVectorD * Errors_Gauss
 Vector with error values using Gaussian fit. More...
 
TVectorD * Means_HPD
 Vector with mean values using Highest Posterior Density. More...
 
TVectorD * Errors_HPD
 Vector with error values using Highest Posterior Density. More...
 
TVectorD * Errors_HPD_Positive
 Vector with positive error (right hand side) values using Highest Posterior Density. More...
 
TVectorD * Errors_HPD_Negative
 Vector with negative error (left hand side) values using Highest Posterior Density. More...
 
TMatrixDSym * Covariance
 Posterior Covariance Matrix. More...
 
TMatrixDSym * Correlation
 Posterior Correlation Matrix. More...
 
std::vector< TH1D * > hpost
 Holds 1D Posterior Distributions. More...
 
std::vector< std::vector< TH2D * > > hpost2D
 Holds 2D Posterior Distributions. More...
 
std::unique_ptr< TH2D > hviolin
 Holds violin plot for all dials. More...
 
std::unique_ptr< TH2D > hviolin_prior
 Holds prior violin plot for all dials,. More...
 
double ** ParStep
 Array holding values for all parameters. More...
 
unsigned int * StepNumber
 Step number for step, important if chains were merged. More...
 
int nBins
 Number of bins. More...
 
double DrawRange
 Drawrange for SetMaximum. More...
 
bool CacheMCMC
 MCMC Chain has been cached. More...
 
bool doDiagMCMC
 Doing MCMC Diagnostic. More...
 
int nBatches
 Number of batches for Batched Mean. More...
 
int AutoCorrLag
 LagL used in AutoCorrelation. More...
 
double * ParamSums
 Total parameter sum for each param. More...
 
double ** BatchedAverages
 Values of batched average for every param and batch. More...
 
double ** SampleValues
 Holds the sample values. More...
 
double ** SystValues
 Holds the systs values. More...
 
double * AccProbValues
 Holds all accProb. More...
 
double * AccProbBatchedAverages
 Holds all accProb in batches. More...
 
bool ReweightPosterior
 Whether to apply reweighting weight or not. More...
 
std::string ReweightName
 Name of branch used for chain reweighting. More...
 
double * WeightValue
 Stores value of weight for each step. More...
 

Detailed Description

This class extends MCMC and allow specialised for Oscillation parameters analysis which require specialised hardcoding.

Author
Clarence Wret
Kamil Skwarczynski

Definition at line 10 of file OscProcessor.h.

Constructor & Destructor Documentation

◆ OscProcessor()

_MaCh3_Safe_Include_Start_ _MaCh3_Safe_Include_End_ OscProcessor::OscProcessor ( const std::string &  InputFile)

Constructs an OscProcessor object with the specified input file and options.

Parameters
InputFileThe path to the input file containing MCMC data.
Todo:
Here where we should add all unitarity triangles, fancy Jarlskog studies and other hacky things that only make sense for oscitations

Definition at line 11 of file OscProcessor.cpp.

11  : MCMCProcessor(InputFile) {
12 // ****************************
13  //KS: WARNING this only work when you project from Chain, will nor work when you try SetBranchAddress etc. Turn it on only if you know how to use it
14  PlotJarlskog = false;
15 
22 }
MCMCProcessor(const std::string &InputFile)
Constructs an MCMCProcessor object with the specified input file and options.
int DeltaCPIndex
Index of in the parameter list.
Definition: OscProcessor.h:75
int Sin2Theta12Index
Index of in the parameter list.
Definition: OscProcessor.h:71
int DeltaM2_23Index
Index of in the parameter list.
Definition: OscProcessor.h:77
bool PlotJarlskog
Will plot Jarlskog Invariant using information in the chain.
Definition: OscProcessor.h:52
int Sin2Theta13Index
Index of in the parameter list.
Definition: OscProcessor.h:69
int Sin2Theta23Index
Index of in the parameter list.
Definition: OscProcessor.h:73
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:48

◆ ~OscProcessor()

OscProcessor::~OscProcessor ( )
virtual

Destroys the OscProcessor object.

Definition at line 26 of file OscProcessor.cpp.

26  {
27 // ****************************
28 
29 }

Member Function Documentation

◆ CalcJarlskog()

double OscProcessor::CalcJarlskog ( const double  s2th13,
const double  s2th23,
const double  s2th12,
const double  dcp 
) const
protected

Calculate Jarlskog Invariant using oscillation parameters.

Parameters
s2th13Value of \( \sin^2\theta_{13} \)
s2th23Value of \( \sin^2\theta_{23} \)
s2th12Value of \( \sin^2\theta_{12} \)
dcpCP-violating phase \( \delta_{\text{CP}} \) (in radians)
Returns
The value of the Jarlskog invariant \( J_{\text{CP}} \) [18]

Definition at line 97 of file OscProcessor.cpp.

97  {
98 // ***************
99  const double s13 = std::sqrt(s2th13);
100  const double s23 = std::sqrt(s2th23);
101  const double s12 = std::sqrt(s2th12);
102  const double sdcp = std::sin(dcp);
103  const double c13 = std::sqrt(1.-s2th13);
104  const double c12 = std::sqrt(1.-s2th12);
105  const double c23 = std::sqrt(1.-s2th23);
106 
107  const double j = s13*c13*c13*s12*c12*s23*c23*sdcp;
108 
109  return j;
110 }

◆ LoadAdditionalInfo()

void OscProcessor::LoadAdditionalInfo ( )
overrideprotectedvirtual

Read the Osc cov file and get the input central values and errors Here we allow Jarlskog Shenanigans.

Todo:
remove this hardcoding (e.g., use a map or enum-to-name function)
Todo:
we should actually calculate central value and prior error but leave it for now...

Reimplemented from MCMCProcessor.

Definition at line 33 of file OscProcessor.cpp.

33  {
34 // ***************
35  // KS: Check if OscParams were enabled, in future we will also get
36  for(size_t i = 0; i < ParameterGroup.size(); i++) {
37  if(ParameterGroup[i] == "Osc"){
38  OscEnabled = true;
39  break;
40  }
41  }
42 
43  if(OscEnabled)
44  {
45  for (int i = 0; i < nDraw; ++i)
46  {
47  //Those keep which parameter type we run currently and relative number
48  const int ParamEnum = ParamType[i];
49  const int ParamNo = i - ParamTypeStartPos[ParameterEnum(ParamEnum)];
50  const std::string CurrentName = ParamNames[ParamEnum][ParamNo].Data();
51 
53  if (CurrentName == "sin2th_13") {
54  Sin2Theta13Index = i;
55  Sin2Theta13Name = CurrentName;
56  } else if (CurrentName == "sin2th_12") {
57  Sin2Theta12Index = i;
58  Sin2Theta12Name = CurrentName;
59  } else if (CurrentName == "sin2th_23") {
60  Sin2Theta23Index = i;
61  Sin2Theta23Name = CurrentName;
62  } else if (CurrentName == "delta_cp") {
63  DeltaCPIndex = i;
64  DeltaCPName = CurrentName;
65  } else if (CurrentName == "delm2_23") {
66  DeltaM2_23Index = i;
67  DeltaM2_23Name = CurrentName;
68  }
69  }
70  } else{
71  MACH3LOG_WARN("Didn't find oscillation parameters");
72  }
73 
75  {
76  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)");
77  BranchNames.push_back("J_cp");
78  ParamType.push_back(kXSecPar);
79  nParam[kXSecPar]++;
80  nDraw++;
81 
83  ParamNom[kXSecPar].push_back( 0. );
84  ParamCentral[kXSecPar].push_back( 0. );
85  ParamErrors[kXSecPar].push_back( 1. );
86  // Push back the name
87  ParamNames[kXSecPar].push_back("J_cp");
88  ParamFlat[kXSecPar].push_back( false );
89  } else if(PlotJarlskog && !OscEnabled) {
90  MACH3LOG_ERROR("Trying to enable Jarlskog without oscillations");
91  throw MaCh3Exception(__FILE__,__LINE__);
92  }
93 }
ParameterEnum
KS: Enum for different covariance classes.
Definition: MCMCProcessor.h:48
@ kXSecPar
Definition: MCMCProcessor.h:49
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:26
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::vector< std::vector< double > > ParamNom
TChain * Chain
Main chain storing all steps etc.
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.
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
Custom exception class for MaCh3 errors.
std::string Sin2Theta13Name
Name of the parameter representing .
Definition: OscProcessor.h:58
std::string Sin2Theta12Name
Name of the parameter representing .
Definition: OscProcessor.h:60
std::string DeltaCPName
Name of the parameter representing (the CP-violating phase).
Definition: OscProcessor.h:64
std::string Sin2Theta23Name
Name of the parameter representing .
Definition: OscProcessor.h:62
std::string DeltaM2_23Name
Name of the parameter representing (mass-squared difference).
Definition: OscProcessor.h:66
bool OscEnabled
Will plot Jarlskog Invariant using information in the chain.
Definition: OscProcessor.h:55

◆ MakeJarlskogPlot()

void OscProcessor::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 
)
protected

Perform Jarlskog Plotting.

Author
Kevin Wood
Note
based on drawJarl_dcpPriorComparison.C

Definition at line 371 of file OscProcessor.cpp.

376  {
377 // ***************
378  MACH3LOG_INFO("Starting {}", __func__);
379  int originalErrorLevel = gErrorIgnoreLevel;
380  gErrorIgnoreLevel = kFatal;
381 
382  // 1-->NH, 0-->both, -1-->IH
383  for(int hierarchy = -1; hierarchy <= 1; hierarchy++)
384  {
385  std::unique_ptr<TH1D> j_hist;
386  std::unique_ptr<TH1D> j_hist_sdcp;
387  if(hierarchy == 1) {
388  j_hist = M3::Clone(jarl_NH.get(), "");
389  j_hist_sdcp = M3::Clone(jarl_NH_flatsindcp.get(), "");
390  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
391  } else if(hierarchy == 0) {
392  j_hist = M3::Clone(jarl.get(), "");
393  j_hist_sdcp = M3::Clone(jarl_flatsindcp.get(), "");
394  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
395  } else if(hierarchy == -1) {
396  j_hist = M3::Clone(jarl_IH.get(), "");
397  j_hist_sdcp = M3::Clone(jarl_IH_flatsindcp.get(), "");
398  j_hist->SetTitle(";J_{CP} #equiv s_{13}c^{2}_{13}s_{12}c_{12}s_{23}c_{23}sin#delta_{CP};Posterior probability");
399  } else {
400  MACH3LOG_ERROR("Invalid hierarchy option. 1 for NH, 0 for both, -1 for IH");
401  throw MaCh3Exception(__FILE__ , __LINE__ );
402  }
403 
404  j_hist->Rebin(7);
405  j_hist_sdcp->Rebin(7);
406 
407  j_hist->SetLineColor(kAzure-2);
408  j_hist_sdcp->SetLineColor(kOrange+1);
409  j_hist->SetLineWidth(2);
410  j_hist_sdcp->SetLineWidth(2);
411 
412  auto StyleAxis = [](TH1* h) {
413  auto xAxis = h->GetXaxis();
414  auto yAxis = h->GetYaxis();
415 
416  xAxis->SetLabelSize(0.04);
417  xAxis->SetLabelFont(132);
418  xAxis->SetTitleSize(0.04);
419  xAxis->SetTitleOffset(0.80);
420  xAxis->SetTitleFont(132);
421  xAxis->SetNdivisions(505);
422  xAxis->SetTickSize(0.04);
423 
424  yAxis->SetLabelSize(0.04);
425  yAxis->SetLabelFont(132);
426  yAxis->SetTitleSize(0.04);
427  yAxis->SetTitleOffset(1.2);
428  yAxis->SetTitleFont(132);
429  yAxis->SetNdivisions(505);
430  yAxis->SetTickSize(0.04);
431  };
432 
433  StyleAxis(j_hist.get());
434 
435  j_hist->GetXaxis()->SetRangeUser(-0.04,0.04);
436  j_hist->Scale(1./j_hist->Integral());
437  j_hist_sdcp->Scale(1./j_hist_sdcp->Integral());
438 
439  std::unique_ptr<TH1D> j_hist_copy = M3::Clone(j_hist.get(), "j_hist_copy");
440  std::unique_ptr<TH1D> j_hist_1sig = M3::Clone(j_hist.get(), "j_hist_1sig");
441  std::unique_ptr<TH1D> j_hist_2sig = M3::Clone(j_hist.get(), "j_hist_2sig");
442  std::unique_ptr<TH1D> j_hist_3sig = M3::Clone(j_hist.get(), "j_hist_3sig");
443 
444  //upper and lower edges
445  double j_bf = j_hist_copy->GetXaxis()->GetBinCenter(j_hist_copy->GetMaximumBin());
446  double j_1sig_low = 9999999.;
447  double j_1sig_up = -9999999.;
448  double j_2sig_low = 9999999.;;
449  double j_2sig_up = -9999999.;
450  double j_3sig_low = 9999999.;;
451  double j_3sig_up = -9999999.;
452 
453 
454  std::unique_ptr<TH1D> j_hist_sdcp_copy = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_copy");
455  std::unique_ptr<TH1D> j_hist_sdcp_1sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_1sig");
456  std::unique_ptr<TH1D> j_hist_sdcp_2sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_2sig");
457  std::unique_ptr<TH1D> j_hist_sdcp_3sig = M3::Clone(j_hist_sdcp.get(), "j_hist_sdcp_3sig");
458 
459  //upper and lower edges
460  double j_sdcp_1sig_low = 9999999.;
461  double j_sdcp_1sig_up = -9999999.;
462  double j_sdcp_2sig_low = 9999999.;;
463  double j_sdcp_2sig_up = -9999999.;
464  double j_sdcp_3sig_low = 9999999.;;
465  double j_sdcp_3sig_up = -9999999.;
466 
467  double contlevel1 = 0.68;
468  double contlevel2 = 0.90;
469  double contlevel4 = 0.99;
470  double contlevel5 = 0.9973;
471  double integral, tsum = 0.;
472 
473  integral = j_hist_copy->Integral();
474 
475  while((tsum/integral)<contlevel5) {
476  double tmax = j_hist_copy->GetMaximum();
477  int bin = j_hist_copy->GetMaximumBin();
478  double xval = j_hist_copy->GetXaxis()->GetBinCenter(bin);
479  double xwidth = j_hist_copy->GetXaxis()->GetBinWidth(bin);
480  if((tsum/integral)<contlevel1) {
481  j_hist_copy->SetBinContent(bin,-1.0);
482  j_hist_1sig->SetBinContent(bin,0.);
483  j_hist_2sig->SetBinContent(bin,0.);
484  j_hist_3sig->SetBinContent(bin,0.);
485  if(xval<j_1sig_low && xval<j_bf) j_1sig_low = xval - xwidth/2.;
486  if(xval>j_1sig_up && xval>j_bf) j_1sig_up = xval + xwidth/2.;
487  }
488  if((tsum/integral)<contlevel2 && (tsum / integral > contlevel1) ) {
489  j_hist_copy->SetBinContent(bin,-5.0);
490  j_hist_2sig->SetBinContent(bin,0.);
491  j_hist_3sig->SetBinContent(bin,0.);
492  if(xval<j_2sig_low && xval<j_bf) j_2sig_low = xval - xwidth/2.;
493  if(xval>j_2sig_up && xval>j_bf) j_2sig_up = xval + xwidth/2.;
494  }
495  if((tsum/integral)<contlevel4 && (tsum / integral > contlevel1) ) {
496  j_hist_copy->SetBinContent(bin,-9.0);
497  j_hist_3sig->SetBinContent(bin,0.);
498  if(xval < j_3sig_low && xval <j_bf) j_3sig_low = xval - xwidth/2.;
499  if(xval > j_3sig_up && xval > j_bf) j_3sig_up = xval + xwidth/2.;
500  }
501  tsum+=tmax;
502  }
503 
504  integral = j_hist_sdcp_copy->Integral();
505  tsum = 0.;
506 
507  while((tsum/integral)<contlevel5) {
508  double tmax = j_hist_sdcp_copy->GetMaximum();
509  int bin = j_hist_sdcp_copy->GetMaximumBin();
510  double xval = j_hist_sdcp_copy->GetXaxis()->GetBinCenter(bin);
511  double xwidth = j_hist_sdcp_copy->GetXaxis()->GetBinWidth(bin);
512  if((tsum/integral)<contlevel1) {
513  j_hist_sdcp_copy->SetBinContent(bin,-1.0);
514  j_hist_sdcp_1sig->SetBinContent(bin,0.);
515  j_hist_sdcp_2sig->SetBinContent(bin,0.);
516  j_hist_sdcp_3sig->SetBinContent(bin,0.);
517  if(xval<j_sdcp_1sig_low && xval<j_bf) j_sdcp_1sig_low = xval - xwidth/2.;
518  if(xval>j_sdcp_1sig_up && xval>j_bf) j_sdcp_1sig_up = xval + xwidth/2.;
519  }
520  if((tsum/integral)<contlevel2 && (tsum / integral > contlevel1) ) {
521  j_hist_sdcp_copy->SetBinContent(bin,-5.0);
522  j_hist_sdcp_2sig->SetBinContent(bin,0.);
523  j_hist_sdcp_3sig->SetBinContent(bin,0.);
524  if(xval<j_sdcp_2sig_low && xval<j_bf) j_sdcp_2sig_low = xval - xwidth/2.;
525  if(xval>j_sdcp_2sig_up && xval>j_bf) j_sdcp_2sig_up = xval + xwidth/2.;
526  }
527  if((tsum/integral)<contlevel4 && (tsum / integral > contlevel1) ) {
528  j_hist_sdcp_copy->SetBinContent(bin,-9.0);
529  j_hist_sdcp_3sig->SetBinContent(bin,0.);
530  if(xval<j_sdcp_3sig_low && xval<j_bf) j_sdcp_3sig_low = xval - xwidth/2.;
531  if(xval>j_sdcp_3sig_up && xval>j_bf) j_sdcp_3sig_up = xval + xwidth/2.;
532  }
533  tsum+=tmax;
534  }
535 
536  j_hist_1sig->SetLineStyle(9);
537  j_hist_sdcp_1sig->SetLineStyle(9);
538  j_hist_2sig->SetLineStyle(7);
539  j_hist_sdcp_2sig->SetLineStyle(7);
540  j_hist_3sig->SetLineStyle(2);
541  j_hist_sdcp_3sig->SetLineStyle(2);
542 
543  auto ldash = std::make_unique<TH1D>("ldash", "solid", 10, -0.04, 0.04);
544  auto sdash = std::make_unique<TH1D>("sdash", "dashed", 10, -0.04, 0.04);
545  auto fdash = std::make_unique<TH1D>("fdash", "fdashed",10, -0.04, 0.04);
546  ldash->SetLineColor(kBlack);
547  sdash->SetLineColor(kBlack);
548  fdash->SetLineColor(kBlack);
549  ldash->SetLineWidth(2);
550  sdash->SetLineWidth(2);
551  fdash->SetLineWidth(2);
552  ldash->SetLineStyle(9);
553  sdash->SetLineStyle(7);
554  fdash->SetLineStyle(2);
555 
556  double vertUp = 0.5 * j_hist->GetMaximum();
557  auto jline_1sig_low = std::make_unique<TLine>(j_1sig_low, 0., j_1sig_low, vertUp);
558  auto jline_2sig_low = std::make_unique<TLine>(j_2sig_low, 0., j_2sig_low, vertUp);
559  auto jline_3sig_low = std::make_unique<TLine>(j_3sig_low, 0., j_3sig_low, vertUp);
560 
561  auto jline_1sig_up = std::make_unique<TLine>(j_1sig_up, 0., j_1sig_up,vertUp);
562  auto jline_2sig_up = std::make_unique<TLine>(j_2sig_up, 0., j_2sig_up,vertUp);
563  auto jline_3sig_up = std::make_unique<TLine>(j_3sig_up, 0., j_3sig_up,vertUp);
564 
565  auto jline_sdcp_1sig_low = std::make_unique<TLine>(j_sdcp_1sig_low, 0., j_sdcp_1sig_low, vertUp);
566  auto jline_sdcp_2sig_low = std::make_unique<TLine>(j_sdcp_2sig_low, 0., j_sdcp_2sig_low, vertUp);
567  auto jline_sdcp_3sig_low = std::make_unique<TLine>(j_sdcp_3sig_low, 0., j_sdcp_3sig_low, vertUp);
568 
569  auto jline_sdcp_1sig_up = std::make_unique<TLine>(j_sdcp_1sig_up, 0., j_sdcp_1sig_up, vertUp);
570  auto jline_sdcp_2sig_up = std::make_unique<TLine>(j_sdcp_2sig_up, 0., j_sdcp_2sig_up, vertUp);
571  auto jline_sdcp_3sig_up = std::make_unique<TLine>(j_sdcp_3sig_up, 0., j_sdcp_3sig_up, vertUp);
572 
573  double arrowLength = 0.003;
574  double arrowHeight = vertUp;
575 
576  auto MakeArrow = [&](double x, Color_t color, Width_t width) -> std::unique_ptr<TArrow> {
577  auto arrow = std::make_unique<TArrow>(x, arrowHeight, x - arrowLength, arrowHeight, 0.02, ">");
578  arrow->SetLineColor(color);
579  arrow->SetLineWidth(width);
580  return arrow;
581  };
582 
583  auto j_arrow_1sig_up = MakeArrow(j_1sig_up, j_hist_1sig->GetLineColor(), j_hist_1sig->GetLineWidth());
584  auto j_arrow_2sig_up = MakeArrow(j_2sig_up, j_hist_2sig->GetLineColor(), j_hist_2sig->GetLineWidth());
585  auto j_arrow_3sig_up = MakeArrow(j_3sig_up, j_hist_3sig->GetLineColor(), j_hist_3sig->GetLineWidth());
586 
587  auto j_sdcp_arrow_1sig_up = MakeArrow(j_sdcp_1sig_up, j_hist_sdcp_1sig->GetLineColor(), j_hist_sdcp_1sig->GetLineWidth());
588  auto j_sdcp_arrow_2sig_up = MakeArrow(j_sdcp_2sig_up, j_hist_sdcp_2sig->GetLineColor(), j_hist_sdcp_2sig->GetLineWidth());
589  auto j_sdcp_arrow_3sig_up = MakeArrow(j_sdcp_3sig_up, j_hist_sdcp_3sig->GetLineColor(), j_hist_sdcp_3sig->GetLineWidth());
590 
591  MACH3LOG_DEBUG("j_1sig_low = {:.4f}, j_2sig_low = {:.4f}, j_3sig_low = {:.4f}", j_1sig_low, j_2sig_low, j_3sig_low);
592  MACH3LOG_DEBUG("j_1sig_up = {:.4f}, j_2sig_up = {:.4f}, j_3sig_up = {:.4f}", j_1sig_up, j_2sig_up, j_3sig_up);
593 
594  auto CopyLineStyle = [](const TH1D* src, TLine* dst) {
595  dst->SetLineColor(src->GetLineColor());
596  dst->SetLineStyle(src->GetLineStyle());
597  dst->SetLineWidth(src->GetLineWidth());
598  };
599 
600  CopyLineStyle(j_hist_1sig.get(), jline_1sig_low.get());
601  CopyLineStyle(j_hist_1sig.get(), jline_1sig_up.get());
602  CopyLineStyle(j_hist_2sig.get(), jline_2sig_low.get());
603  CopyLineStyle(j_hist_2sig.get(), jline_2sig_up.get());
604  CopyLineStyle(j_hist_3sig.get(), jline_3sig_low.get());
605  CopyLineStyle(j_hist_3sig.get(), jline_3sig_up.get());
606 
607  CopyLineStyle(j_hist_sdcp_1sig.get(), jline_sdcp_1sig_low.get());
608  CopyLineStyle(j_hist_sdcp_1sig.get(), jline_sdcp_1sig_up.get());
609  CopyLineStyle(j_hist_sdcp_2sig.get(), jline_sdcp_2sig_low.get());
610  CopyLineStyle(j_hist_sdcp_2sig.get(), jline_sdcp_2sig_up.get());
611  CopyLineStyle(j_hist_sdcp_3sig.get(), jline_sdcp_3sig_low.get());
612  CopyLineStyle(j_hist_sdcp_3sig.get(), jline_sdcp_3sig_up.get());
613 
614  auto leg = std::make_unique<TLegend>(0.45, 0.60, 0.75, 0.90);
615  leg->SetTextSize(0.05);
616  leg->SetFillStyle(0);
617  leg->SetNColumns(1);
618  leg->SetTextFont(132);
619  leg->SetBorderSize(0);
620 
621  leg->AddEntry(j_hist.get(), "Prior flat in #delta_{CP}", "l");
622  leg->AddEntry(j_hist_sdcp.get(), "Prior flat in sin#delta_{CP}", "l");
623  leg->AddEntry(ldash.get(), "68% CI", "l");
624  leg->AddEntry(sdash.get(), "90% CI", "l");
625  leg->AddEntry(fdash.get(), "99% CI", "l");
626 
627  j_hist->GetYaxis()->SetRangeUser(0., j_hist->GetMaximum()*1.15);
628  j_hist->Draw("h");
629  j_hist_sdcp->Draw("same h");
630 
631  jline_sdcp_1sig_up->Draw("same");
632  jline_sdcp_2sig_up->Draw("same");
633  jline_sdcp_3sig_up->Draw("same");
634  jline_1sig_up->Draw("same");
635  jline_2sig_up->Draw("same");
636  jline_3sig_up->Draw("same");
637 
638  j_arrow_1sig_up->Draw();
639  j_arrow_2sig_up->Draw();
640  j_arrow_3sig_up->Draw();
641  j_sdcp_arrow_1sig_up->Draw();
642  j_sdcp_arrow_2sig_up->Draw();
643  j_sdcp_arrow_3sig_up->Draw();
644  leg->Draw("same");
645 
646  auto ttext = std::make_unique<TText>();
647  ttext->SetNDC(); // Use normalized device coordinates
648  ttext->SetTextSize(0.03); // Adjust size as needed
649  ttext->SetTextAlign(13); // Align left-top
650 
651  if (hierarchy == 1) ttext->DrawText(0.15, 0.85, "Normal Ordering");
652  else if (hierarchy == 0) ttext->DrawText(0.15, 0.85, "Both Orderings");
653  else if (hierarchy == -1) ttext->DrawText(0.15, 0.85, "Inverted Ordering");
654 
655  gPad->RedrawAxis();
656  Posterior->Update();
657  gPad->Update();
658 
659  Posterior->Print(CanvasName);
660 
661  if(hierarchy == 1) Posterior->Write("jarl1D_NH_comp");
662  else if(hierarchy == 0) Posterior->Write("jarl1D_both_comp");
663  else if(hierarchy == -1) Posterior->Write("jarl1D_IH_comp");
664  }
665 
666  gErrorIgnoreLevel = originalErrorLevel;
667 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:24
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:25
std::unique_ptr< TCanvas > Posterior
Fancy canvas used for our beautiful plots.
TString CanvasName
Name of canvas which help to save to the sample pdf.
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.

◆ PerformJarlskogAnalysis()

void OscProcessor::PerformJarlskogAnalysis ( )

Perform Several Jarlskog Plotting.

Author
Kevin Wood
Note
based on makeJarlskog.C

Definition at line 136 of file OscProcessor.cpp.

136  {
137 // ***************
138  if(!OscEnabled ||
144  {
145  MACH3LOG_WARN("Will not {}, as oscillation parameters are missing", __func__);
146  return;
147  }
148  MACH3LOG_INFO("Starting {}", __func__);
149 
150  bool DoReweight = false;
151 
152  double s2th13, s2th23, s2th12, dcp, dm2 = M3::_BAD_DOUBLE_;
153  double weight = 1.0;
154  std::pair<double, double> Sin13_NewPrior;
155 
156  // Now read the MCMC file
157  TFile *TempFile = new TFile((MCMCFile + ".root").c_str(), "open");
158 
159  // Get the settings for the MCMC
160  TMacro *Config = TempFile->Get<TMacro>("Reweight_Config");
161 
162  if (Config != nullptr) {
163  YAML::Node Settings = TMacroToYAML(*Config);
164  if(CheckNodeExists(Settings, "Weight", Sin2Theta13Name)) {
165  Sin13_NewPrior = Get<std::pair<double, double>>(Settings["Weight"][Sin2Theta13Name], __FILE__, __LINE__);
166  MACH3LOG_INFO("Found Weight in chain, using RC reweighting with new priors {} +- {}", Sin13_NewPrior.first, Sin13_NewPrior.second);
167  DoReweight = true;
168  }
169  }
170 
171  TempFile->Close();
172  delete TempFile;
173 
174  TDirectory *JarlskogDir = OutputFile->mkdir("Jarlskog");
175  JarlskogDir->cd();
176 
177  unsigned int step = 0;
178  Chain->SetBranchStatus("*", false);
179 
180  Chain->SetBranchStatus(Sin2Theta13Name.c_str(), true);
181  Chain->SetBranchAddress(Sin2Theta13Name.c_str(), &s2th13);
182 
183  Chain->SetBranchStatus(Sin2Theta23Name.c_str(), true);
184  Chain->SetBranchAddress(Sin2Theta23Name.c_str(), &s2th23);
185 
186  Chain->SetBranchStatus(Sin2Theta12Name.c_str(), true);
187  Chain->SetBranchAddress(Sin2Theta12Name.c_str(), &s2th12);
188 
189  Chain->SetBranchStatus(DeltaCPName.c_str(), true);
190  Chain->SetBranchAddress(DeltaCPName.c_str(), &dcp);
191 
192  Chain->SetBranchStatus(DeltaM2_23Name.c_str(), true);
193  Chain->SetBranchAddress(DeltaM2_23Name.c_str(), &dm2);
194 
195  Chain->SetBranchStatus("step", true);
196  Chain->SetBranchAddress("step", &step);
197 
198  if(DoReweight) {
199  Chain->SetBranchStatus("Weight", true);
200  Chain->SetBranchAddress("Weight", &weight);
201  } else {
202  MACH3LOG_WARN("Not applying reweighting weight");
203  weight = 1.0;
204  }
205 
206  // Original histograms
207  auto jarl = std::make_unique<TH1D>("jarl", "jarl", 1000, -0.05, 0.05);
208  jarl->SetDirectory(nullptr);
209  auto jarl_th23 = std::make_unique<TH2D>("jarl_th23", "jarl_th23", 500, -0.05, 0.05, 500, 0.3, 0.7);
210  jarl_th23->SetDirectory(nullptr);
211  auto jarl_dcp = std::make_unique<TH2D>("jarl_dcp", "jarl_dcp", 500, -0.05, 0.05, 500, -1. * TMath::Pi(), TMath::Pi());
212  jarl_dcp->SetDirectory(nullptr);
213 
214  jarl->SetTitle("Jarlskog Invariant;J #equiv s_{13}c_{13}^{2}s_{12}c_{12}s_{23}c_{23}sin#delta;Posterior probability");
215  jarl_th23->SetTitle("Jarlskog Invariant;J #equiv s_{13}c_{13}^{2}s_{12}c_{12}s_{23}c_{23}sin#delta;Posterior probability");
216 
217  // Clones
218  auto jarl_IH = M3::Clone(jarl.get(), "jarl_IH");
219  auto jarl_NH = M3::Clone(jarl.get(), "jarl_NH");
220 
221  auto jarl_th23_IH = M3::Clone(jarl_th23.get(), "jarl_th23_IH");
222  auto jarl_th23_NH = M3::Clone(jarl_th23.get(), "jarl_th23_NH");
223 
224  auto jarl_dcp_IH = M3::Clone(jarl_dcp.get(), "jarl_dcp_IH");
225  auto jarl_dcp_NH = M3::Clone(jarl_dcp.get(), "jarl_dcp_NH");
226 
227  auto jarl_flatsindcp = M3::Clone(jarl.get(), "jarl_flatsindcp");
228  auto jarl_IH_flatsindcp = M3::Clone(jarl.get(), "jarl_IH_flatsindcp");
229  auto jarl_NH_flatsindcp = M3::Clone(jarl.get(), "jarl_NH_flatsindcp");
230 
231  auto jarl_th23_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_flatsindcp");
232  auto jarl_th23_IH_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_IH_flatsindcp");
233  auto jarl_th23_NH_flatsindcp = M3::Clone(jarl_th23.get(), "jarl_th23_NH_flatsindcp");
234 
235  auto jarl_prior = M3::Clone(jarl.get(), "jarl_prior");
236  auto jarl_prior_flatsindcp = M3::Clone(jarl.get(), "jarl_prior_flatsindcp");
237  std::unique_ptr<TH1D> jarl_wRC_prior, jarl_wRC_prior_flatsindcp, jarl_wRC_prior_t2kth23;
238  // Only use this if chain has reweigh weight [mostly coming from Reactor Constrains]
239  if(DoReweight){
240  jarl_wRC_prior = M3::Clone(jarl.get(), "jarl_wRC_prior");
241  jarl_wRC_prior_flatsindcp = M3::Clone(jarl.get(), "jarl_wRC_prior_flatsindcp");
242  jarl_wRC_prior_t2kth23 = M3::Clone(jarl.get(), "jarl_wRC_prior_flatsindcp");
243  }
244 
245  // to apply a prior that is flat in sin(dcp) intead of dcp
246  auto prior3 = std::make_unique<TF1>("prior3", "TMath::Abs(TMath::Cos(x))");
247 
248  // T2K prior is flat (and uncorrelated) in dcp, sin^2(th13), sin^2(th23)
249  auto randGen = std::make_unique<TRandom3>(0);
250  const Long64_t countwidth = nEntries/5;
251 
252  for(int i = 0; i < nEntries; ++i) {
253  if (i % countwidth == 0) {
256  } else {
257  Chain->GetEntry(i);
258  }
259 
260  if(step < BurnInCut) continue; // burn-in cut
261 
262  const double j = CalcJarlskog(s2th13, s2th23, s2th12, dcp);
263  const double prior_weight = prior3->Eval(dcp);
264 
265  jarl->Fill(j, weight);
266  jarl_th23->Fill(j, s2th23, weight);
267  jarl_dcp->Fill(j, dcp, weight);
268 
269  jarl_flatsindcp->Fill(j, prior_weight*weight);
270  jarl_th23_flatsindcp->Fill(j, s2th23, prior_weight*weight);
271 
272  const double prior_s2th13 = SamplePriorForParam(Sin2Theta13Index, randGen, {0.,1.});
273  const double prior_s2th23 = SamplePriorForParam(Sin2Theta23Index, randGen, {0.,1.});
274  const double prior_s2th12 = SamplePriorForParam(Sin2Theta12Index, randGen, {0.,1.});
275  const double prior_dcp = SamplePriorForParam(DeltaCPIndex, randGen, {-1.*TMath::Pi(),TMath::Pi()});
276  // KS: This is hardcoded but we always assume flat in delta CP so probably fine
277  const double prior_sindcp = randGen->Uniform(-1., 1.);
278 
279  const double prior_s13 = std::sqrt(prior_s2th13);
280  const double prior_s23 = std::sqrt(prior_s2th23);
281  const double prior_s12 = std::sqrt(prior_s2th12);
282  const double prior_sdcp = std::sin(prior_dcp);
283  const double prior_c13 = std::sqrt(1.-prior_s2th13);
284  const double prior_c12 = std::sqrt(1.-prior_s2th12);
285  const double prior_c23 = std::sqrt(1.-prior_s2th23);
286  const double prior_j = prior_s13*prior_c13*prior_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sdcp;
287  const double prior_flatsindcp_j = prior_s13*prior_c13*prior_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sindcp;
288 
289  jarl_prior->Fill(prior_j);
290  jarl_prior_flatsindcp->Fill(prior_flatsindcp_j);
291 
292  if(DoReweight) {
293  const double prior_wRC_s2th13 = randGen->Gaus(Sin13_NewPrior.first, Sin13_NewPrior.second);
294  const double prior_wRC_s13 = std::sqrt(prior_wRC_s2th13);
295  const double prior_wRC_c13 = std::sqrt(1.-prior_wRC_s2th13);
296  const double prior_wRC_j = prior_wRC_s13*prior_wRC_c13*prior_wRC_c13*prior_s12*prior_c12*prior_s23*prior_c23*prior_sdcp;
297  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;
298  const double s23 = std::sqrt(s2th23);
299  const double c23 = std::sqrt(1.-s2th23);
300 
301  jarl_wRC_prior->Fill(prior_wRC_j);
302  jarl_wRC_prior_flatsindcp->Fill(prior_wRC_flatsindcp_j);
303  jarl_wRC_prior_t2kth23->Fill(prior_wRC_s13*prior_wRC_c13*prior_wRC_c13*prior_s12*prior_c12*s23*c23*prior_sdcp);
304  }
305 
306  if(dm2 > 0.) {
307  jarl_NH->Fill(j, weight);
308  jarl_th23_NH->Fill(j, s2th23, weight);
309  jarl_dcp_NH->Fill(j, dcp, weight);
310  jarl_NH_flatsindcp->Fill(j, prior_weight*weight);
311  jarl_th23_NH_flatsindcp->Fill(j, s2th23, prior_weight*weight);
312  }
313  else if(dm2 < 0.) {
314  jarl_IH->Fill(j, weight);
315  jarl_th23_IH->Fill(j, s2th23, weight);
316  jarl_dcp_IH->Fill(j, dcp, weight);
317  jarl_IH_flatsindcp->Fill(j, prior_weight*weight);
318  jarl_th23_IH_flatsindcp->Fill(j, s2th23, prior_weight*weight);
319  }
320  }
321 
322  jarl->Write("jarlskog_both");
323  jarl_NH->Write("jarlskog_NH");
324  jarl_IH->Write("jarlskog_IH");
325  jarl_th23->Write("jarlskog_th23_both");
326  jarl_th23_NH->Write("jarlskog_th23_NH");
327  jarl_th23_IH->Write("jarlskog_th23_IH");
328 
329  jarl_dcp->Write("jarlskog_dcp_both");
330  jarl_dcp_NH->Write("jarlskog_dcp_NH");
331  jarl_dcp_IH->Write("jarlskog_dcp_IH");
332 
333 
334  jarl_flatsindcp->Write("jarlskog_both_flatsindcp");
335  jarl_NH_flatsindcp->Write("jarlskog_NH_flatsindcp");
336  jarl_IH_flatsindcp->Write("jarlskog_IH_flatsindcp");
337  jarl_th23_flatsindcp->Write("jarlskog_th23_both_flatsindcp");
338  jarl_th23_NH_flatsindcp->Write("jarlskog_th23_NH_flatsindcp");
339  jarl_th23_IH_flatsindcp->Write("jarlskog_th23_IH_flatsindcp");
340 
341  jarl_prior->Write("jarl_prior");
342  jarl_prior_flatsindcp->Write("jarl_prior_flatsindcp");
343  if(DoReweight) {
344  jarl_wRC_prior->Write("jarl_wRC_prior");
345  jarl_wRC_prior_flatsindcp->Write("jarl_wRC_prior_flatsindcp");
346  jarl_wRC_prior_t2kth23->Write("jarl_wRC_prior_t2kth23");
347  }
348 
349  MakeJarlskogPlot(jarl, jarl_flatsindcp,
350  jarl_NH, jarl_NH_flatsindcp,
351  jarl_IH, jarl_IH_flatsindcp);
352 
353  // Perform Savage Dickey analysis
354  if(DoReweight) {
355  SavageDickeyPlot(jarl, jarl_wRC_prior, "Jarlskog flat #delta_{CP}", 0);
356  SavageDickeyPlot(jarl_flatsindcp, jarl_wRC_prior_flatsindcp, "Jarlskog flat sin#delta_{CP}", 0);
357  } else {
358  SavageDickeyPlot(jarl, jarl_prior, "Jarlskog flat #delta_{CP}", 0);
359  SavageDickeyPlot(jarl_flatsindcp, jarl_prior_flatsindcp, "Jarlskog flat sin#delta_{CP}", 0);
360  }
361 
362  JarlskogDir->Close();
363  delete JarlskogDir;
364 
365  Chain->SetBranchStatus("*", true);
366  OutputFile->cd();
367 }
YAML::Node TMacroToYAML(const TMacro &macro)
KS: Convert a ROOT TMacro object to a YAML node.
Definition: YamlHelper.h:147
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Definition: YamlHelper.h:55
TFile * OutputFile
The output file.
std::string MCMCFile
Name of MCMC file.
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.
double SamplePriorForParam(const int paramIndex, const std::unique_ptr< TRandom3 > &randGen, const std::vector< double > &FlatBounds) const
Draw Prior value.
double CalcJarlskog(const double s2th13, const double s2th23, const double s2th12, const double dcp) const
Calculate Jarlskog Invariant using oscillation parameters.
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.
constexpr static const double _BAD_DOUBLE_
Default value used for double initialisation.
Definition: Core.h:46
void PrintProgressBar(const Long64_t Done, const Long64_t All)
KS: Simply print progress bar.
Definition: Monitor.cpp:213
void EstimateDataTransferRate(TChain *chain, const Long64_t entry)
KS: Check what CPU you are using.
Definition: Monitor.cpp:196

◆ SamplePriorForParam()

double OscProcessor::SamplePriorForParam ( const int  paramIndex,
const std::unique_ptr< TRandom3 > &  randGen,
const std::vector< double > &  FlatBounds 
) const
protected

Draw Prior value.

Definition at line 113 of file OscProcessor.cpp.

113  {
114 // ***************
115  TString Title = "";
116  double Prior = 1.0, PriorError = 1.0;
117  bool FlatPrior = false;
118 
119  // Get info for this parameter
120  GetNthParameter(paramIndex, Prior, PriorError, Title);
121 
122  ParameterEnum ParType = ParamType[paramIndex];
123  int ParamTemp = paramIndex - ParamTypeStartPos[ParType];
124  FlatPrior = ParamFlat[ParType][ParamTemp];
125 
126  if (FlatPrior) {
127  return randGen->Uniform(FlatBounds[0], FlatBounds[1]);
128  } else {
129  // Gaussian prior centered at Prior with width PriorError
130  return randGen->Gaus(Prior, PriorError);
131  }
132 }
void GetNthParameter(const int param, double &Prior, double &PriorError, TString &Title) const
Get properties of parameter by passing it number.

Member Data Documentation

◆ DeltaCPIndex

int OscProcessor::DeltaCPIndex
protected

Index of \(\delta_{\mathrm{CP}}\) in the parameter list.

Definition at line 75 of file OscProcessor.h.

◆ DeltaCPName

std::string OscProcessor::DeltaCPName
protected

Name of the parameter representing \(\delta_{\mathrm{CP}}\) (the CP-violating phase).

Definition at line 64 of file OscProcessor.h.

◆ DeltaM2_23Index

int OscProcessor::DeltaM2_23Index
protected

Index of \(\Delta m^2_{32}\) in the parameter list.

Definition at line 77 of file OscProcessor.h.

◆ DeltaM2_23Name

std::string OscProcessor::DeltaM2_23Name
protected

Name of the parameter representing \(\Delta m^2_{32}\) (mass-squared difference).

Definition at line 66 of file OscProcessor.h.

◆ OscEnabled

bool OscProcessor::OscEnabled
protected

Will plot Jarlskog Invariant using information in the chain.

Definition at line 55 of file OscProcessor.h.

◆ PlotJarlskog

bool OscProcessor::PlotJarlskog
protected

Will plot Jarlskog Invariant using information in the chain.

Definition at line 52 of file OscProcessor.h.

◆ Sin2Theta12Index

int OscProcessor::Sin2Theta12Index
protected

Index of \(\sin^2\theta_{12}\) in the parameter list.

Definition at line 71 of file OscProcessor.h.

◆ Sin2Theta12Name

std::string OscProcessor::Sin2Theta12Name
protected

Name of the parameter representing \(\sin^2\theta_{12}\).

Definition at line 60 of file OscProcessor.h.

◆ Sin2Theta13Index

int OscProcessor::Sin2Theta13Index
protected

Index of \(\sin^2\theta_{13}\) in the parameter list.

Definition at line 69 of file OscProcessor.h.

◆ Sin2Theta13Name

std::string OscProcessor::Sin2Theta13Name
protected

Name of the parameter representing \(\sin^2\theta_{13}\).

Definition at line 58 of file OscProcessor.h.

◆ Sin2Theta23Index

int OscProcessor::Sin2Theta23Index
protected

Index of \(\sin^2\theta_{23}\) in the parameter list.

Definition at line 73 of file OscProcessor.h.

◆ Sin2Theta23Name

std::string OscProcessor::Sin2Theta23Name
protected

Name of the parameter representing \(\sin^2\theta_{23}\).

Definition at line 62 of file OscProcessor.h.


The documentation for this class was generated from the following files: