MaCh3  2.5.1
Reference Guide
Functions
FitterBase.cpp File Reference
#include "FitterBase.h"
#include "TRandom.h"
#include "TStopwatch.h"
#include "TTree.h"
#include "TGraphAsymmErrors.h"
Include dependency graph for FitterBase.cpp:

Go to the source code of this file.

Functions

void WriteHistograms (TH1 *hist, const std::string &baseName)
 Helper to write histograms. More...
 
void WriteHistogramsByMode (SampleHandlerInterface *sample, const std::string &suffix, const bool by_mode, const bool by_channel, const std::vector< TDirectory * > &SampleDir)
 Generic histogram writer - should make main code more palatable. More...
 

Function Documentation

◆ WriteHistograms()

void WriteHistograms ( TH1 *  hist,
const std::string &  baseName 
)

Helper to write histograms.

Definition at line 1317 of file FitterBase.cpp.

1317  {
1318 // *************************
1319  if (!hist) return;
1320  hist->SetTitle(baseName.c_str());
1321  // Get the class name of the histogram
1322  TString className = hist->ClassName();
1323 
1324  // Set the appropriate axis title based on the histogram type
1325  if (className.Contains("TH1")) {
1326  hist->GetYaxis()->SetTitle("Events");
1327  } else if (className.Contains("TH2")) {
1328  hist->GetZaxis()->SetTitle("Events");
1329  }
1330  hist->Write(baseName.c_str());
1331 }

◆ WriteHistogramsByMode()

void WriteHistogramsByMode ( SampleHandlerInterface sample,
const std::string &  suffix,
const bool  by_mode,
const bool  by_channel,
const std::vector< TDirectory * > &  SampleDir 
)

Generic histogram writer - should make main code more palatable.

Definition at line 1335 of file FitterBase.cpp.

1339  {
1340 // *************************
1341  MaCh3Modes *modes = sample->GetMaCh3Modes();
1342  for (int iSample = 0; iSample < sample->GetNSamples(); ++iSample) {
1343  SampleDir[iSample]->cd();
1344  const std::string sampleName = sample->GetSampleTitle(iSample);
1345  for(int iDim1 = 0; iDim1 < sample->GetNDim(iSample); iDim1++) {
1346  std::string ProjectionName = sample->GetKinVarName(iSample, iDim1);
1347  std::string ProjectionSuffix = "_1DProj" + std::to_string(iDim1);
1348 
1349  // Probably a better way of handling this logic
1350  if (by_mode) {
1351  for (int iMode = 0; iMode < modes->GetNModes(); ++iMode) {
1352  auto modeHist = sample->Get1DVarHistByModeAndChannel(iSample, ProjectionName, iMode);
1353  WriteHistograms(modeHist.get(), sampleName + "_" + modes->GetMaCh3ModeName(iMode) + ProjectionSuffix + suffix);
1354  }
1355  }
1356 
1357  if (by_channel) {
1358  for (int iChan = 0; iChan < sample->GetNOscChannels(iSample); ++iChan) {
1359  auto chanHist = sample->Get1DVarHistByModeAndChannel(iSample, ProjectionName, -1, iChan); // -1 skips over mode plotting
1360  WriteHistograms(chanHist.get(), sampleName + "_" + sample->GetFlavourName(iSample, iChan) + ProjectionSuffix + suffix);
1361  }
1362  }
1363 
1364  if (by_mode && by_channel) {
1365  for (int iMode = 0; iMode < modes->GetNModes(); ++iMode) {
1366  for (int iChan = 0; iChan < sample->GetNOscChannels(iSample); ++iChan) {
1367  auto hist = sample->Get1DVarHistByModeAndChannel(iSample, ProjectionName, iMode, iChan);
1368  WriteHistograms(hist.get(), sampleName + "_" + modes->GetMaCh3ModeName(iMode) + "_" + sample->GetFlavourName(iSample, iChan) + ProjectionSuffix + suffix);
1369  }
1370  }
1371  }
1372 
1373  if (!by_mode && !by_channel) {
1374  auto hist = sample->Get1DVarHist(iSample, ProjectionName);
1375  WriteHistograms(hist.get(), sampleName + ProjectionSuffix + suffix);
1376  // Only for 2D and Beyond
1377  for (int iDim2 = iDim1 + 1; iDim2 < sample->GetNDim(iSample); ++iDim2) {
1378  // Get the names for the two dimensions
1379  std::string XVarName = sample->GetKinVarName(iSample, iDim1);
1380  std::string YVarName = sample->GetKinVarName(iSample, iDim2);
1381 
1382  // Get the 2D histogram for this pair
1383  auto hist2D = sample->Get2DVarHist(iSample, XVarName, YVarName);
1384 
1385  // Write the histogram
1386  std::string suffix2D = "_2DProj_" + std::to_string(iDim1) + "_vs_" + std::to_string(iDim2) + suffix;
1387  WriteHistograms(hist2D.get(), sampleName + suffix2D);
1388  }
1389  }
1390  }
1391  }
1392 }
void WriteHistograms(TH1 *hist, const std::string &baseName)
Helper to write histograms.
KS: Class describing MaCh3 modes used in the analysis, it is being initialised from config.
Definition: MaCh3Modes.h:142
int GetNModes() const
KS: Get number of modes, keep in mind actual number is +1 greater due to unknown category.
Definition: MaCh3Modes.h:155
std::string GetMaCh3ModeName(const int Index) const
KS: Get normal name of mode, if mode not known you will get UNKNOWN_BAD.
Definition: MaCh3Modes.cpp:156
virtual std::unique_ptr< TH1 > Get1DVarHistByModeAndChannel(const int iSample, const std::string &ProjectionVar_Str, const int kModeToFill=-1, const int kChannelToFill=-1, const int WeightStyle=0)=0
Build a 1D histogram for a given variable, optionally filtered by mode and channel.
virtual std::unique_ptr< TH1 > Get1DVarHist(const int iSample, const std::string &ProjectionVar, const std::vector< KinematicCut > &EventSelectionVec={}, int WeightStyle=0, const std::vector< KinematicCut > &SubEventSelectionVec={})=0
Return 1D projection of MC into given 1D variable (doesn't have to be variable used in the fit)
virtual std::string GetFlavourName(const int iSample, const int iChannel) const =0
Get the flavour name for a given sample and oscillation channel.
MaCh3Modes * GetMaCh3Modes() const
Return pointer to MaCh3 modes.
virtual int GetNOscChannels(const int iSample) const =0
Get number of oscillation channels for a single sample.
virtual M3::int_t GetNSamples()
returns total number of samples
virtual std::string GetSampleTitle(const int iSample) const =0
Get fancy title for specified samples.
virtual std::string GetKinVarName(const int iSample, const int Dimension) const =0
Return Kinematic Variable name for specified sample and dimension for example "Reconstructed_Neutrino...
virtual std::unique_ptr< TH2 > Get2DVarHist(const int iSample, const std::string &ProjectionVarX, const std::string &ProjectionVarY, const std::vector< KinematicCut > &EventSelectionVec={}, const int WeightStyle=0, const std::vector< KinematicCut > &SubEventSelectionVec={})=0
Build a 2D projection of MC events into specified variables.
virtual int GetNDim(const int Sample) const =0
DB Get what dimensionality binning for given sample has.