MaCh3  2.2.3
Reference Guide
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
BinnedSplineHandler Class Referenceabstract

Bin-by-bin class calculating response for spline parameters. More...

#include <Splines/BinnedSplineHandler.h>

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

Public Member Functions

 BinnedSplineHandler (ParameterHandlerGeneric *xsec_, MaCh3Modes *Modes_)
 Constructor. More...
 
virtual ~BinnedSplineHandler ()
 Destructor. More...
 
void Evaluate () override
 CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; probably the best one here! Same thing but pass parameter spline segments instead of variations. More...
 
void AddSample (const std::string &SampleName, const std::string &SampleTitle, const std::vector< std::string > &OscChanFileNames, const std::vector< std::string > &SplineVarNames)
 add oscillation channel to spline monolith More...
 
void TransferToMonolith ()
 flatten multidimensional spline array into proper monolith More...
 
void cleanUpMemory ()
 Remove setup variables not needed for spline evaluations. More...
 
virtual void FillSampleArray (std::string SampleTitle, std::vector< std::string > OscChanFileNames)
 Loads and processes splines from ROOT files for a given sample. More...
 
std::vector< std::vector< int > > StripDuplicatedModes (const std::vector< std::vector< int > > &InputVector)
 Check if there are any repeated modes. This is used to reduce the number of modes in case many interaction modes get averaged into one spline. More...
 
std::vector< std::vector< int > > GetEventSplines (const std::string &SampleTitle, int iOscChan, int EventMode, double Var1Val, double Var2Val, double Var3Val)
 Return the splines which affect a given event. More...
 
std::vector< TAxis * > FindSplineBinning (const std::string &FileName, const std::string &SampleTitle)
 Grab histograms with spline binning. More...
 
int CountNumberOfLoadedSplines (bool NonFlat=false, int Verbosity=0)
 
std::string getDimLabel (const int BinningOpt, const unsigned int Axis) const
 
int GetSampleIndex (const std::string &SampleTitle) const
 Get index of sample based on name. More...
 
bool isValidSplineIndex (const std::string &SampleTitle, int iSyst, int iOscChan, int iMode, int iVar1, int iVar2, int iVar3)
 Ensure we have spline for a given bin. More...
 
void BuildSampleIndexingArray (const std::string &SampleTitle)
 
void PrepForReweight ()
 
void getSplineCoeff_SepMany (int splineindex, M3::float_t *&xArray, M3::float_t *&manyArray)
 
void PrintBinning (TAxis *Axis) const
 
void PrintSampleDetails (const std::string &SampleTitle) const
 Print info like Sample ID of spline params etc. More...
 
void PrintArrayDetails (const std::string &SampleTitle) const
 
const M3::float_tretPointer (const int sample, const int oscchan, const int syst, const int mode, const int var1bin, const int var2bin, const int var3bin) const
 get pointer to spline weight based on bin variables More...
 
void PrepareSplineFile (std::string FileName) override
 KS: Prepare spline file that can be used for fast loading. More...
 
void LoadSplineFile (std::string FileName) override
 KS: Load preprocessed spline file. More...
 
- Public Member Functions inherited from SplineBase
 SplineBase ()
 Constructor. More...
 
virtual ~SplineBase ()
 Destructor. More...
 
virtual std::string GetName () const
 Get class name. More...
 
short int GetNParams () const
 Get number of spline parameters. More...
 

Protected Types

enum  TokenOrdering {
  kSystToken , kModeToken , kVar1BinToken , kVar2BinToken ,
  kVar3BinToken , kNTokens
}
 

Protected Member Functions

void CalcSplineWeights () override
 CPU based code which eval weight for each spline. More...
 
void ModifyWeights () override
 Calc total event weight, not used by Bin-by-bin splines. More...
 
virtual std::vector< std::string > GetTokensFromSplineName (std::string FullSplineName)=0
 
- Protected Member Functions inherited from SplineBase
void FindSplineSegment ()
 CW:Code used in step by step reweighting, Find Spline Segment for each param. More...
 
void PrepareFastSplineInfoDir (std::unique_ptr< TFile > &SplineFile) const
 KS: Prepare Fast Spline Info within SplineFile. More...
 
void LoadFastSplineInfoDir (std::unique_ptr< TFile > &SplineFile)
 KS: Load preprocessed FastSplineInfo. More...
 
void getTF1Coeff (TF1_red *&spl, int &nPoints, float *&coeffs)
 CW: Gets the polynomial coefficients for TF1. More...
 

Protected Attributes

ParameterHandlerGenericxsec
 Pointer to covariance from which we get information about spline params. More...
 
std::vector< std::string > SampleNames
 
std::vector< std::string > SampleTitles
 
std::vector< int > Dimensions
 
std::vector< std::vector< std::string > > DimensionLabels
 
std::vector< int > nSplineParams
 
std::vector< int > nOscChans
 
std::vector< std::vector< std::vector< TAxis * > > > SplineBinning
 
std::vector< std::vector< std::string > > SplineFileParPrefixNames
 
std::vector< std::vector< std::vector< int > > > SplineModeVecs
 
std::vector< std::vector< int > > GlobalSystIndex
 This holds the global spline index and is used to grab the current parameter value to evaluate splines at. Each internal vector will be of size of the number of spline systematics which affect that sample. More...
 
std::vector< std::vector< SplineInterpolation > > SplineInterpolationTypes
 spline interpolation types for each sample. These vectors are from a call to GetSplineInterpolationFromSampleID() More...
 
std::vector< std::string > UniqueSystNames
 name of each spline parameter More...
 
std::vector< int > UniqueSystIndices
 Global index of each spline param, it allows us to match spline ordering with global. More...
 
std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< int > > > > > > > indexvec
 Variables related to determined which modes have splines and which piggy-back of other modes. More...
 
std::vector< int > coeffindexvec
 
std::vector< int > uniquecoeffindices
 Unique coefficient indices. More...
 
std::vector< TSpline3_red * > splinevec_Monolith
 holds each spline object before stripping into coefficient monolith More...
 
int MonolithSize
 
int MonolithIndex
 
int CoeffIndex
 
bool * isflatarray
 Need to keep track of which splines are flat and which aren't. More...
 
M3::float_txcoeff_arr
 x coefficients for each spline More...
 
M3::float_tmanycoeff_arr
 ybcd coefficients for each spline More...
 
std::vector< M3::float_tweightvec_Monolith
 Stores weight from spline evaluation for each single spline. More...
 
std::vector< int > uniquesplinevec_Monolith
 Maps single spline object with single parameter. More...
 
MaCh3ModesModes
 pointer to MaCh3 Mode from which we get spline suffix More...
 
- Protected Attributes inherited from SplineBase
std::vector< FastSplineInfoSplineInfoArray
 
short int * SplineSegments
 
float * ParamValues
 Store parameter values they are not in FastSplineInfo as in case of GPU we need to copy paste it to GPU. More...
 
short int nParams
 Number of parameters that have splines. More...
 

Private Member Functions

void InvestigateMissingSplines () const
 This function will find missing splines in file. More...
 
void LoadSettingsDir (std::unique_ptr< TFile > &SplineFile)
 KS: Load preprocessed Settings. More...
 
void LoadMonolithDir (std::unique_ptr< TFile > &SplineFile)
 KS: Load preprocessed Monolith. More...
 
void LoadIndexDir (std::unique_ptr< TFile > &SplineFile)
 KS: Load preprocessed Index. More...
 
void PrepareSettingsDir (std::unique_ptr< TFile > &SplineFile) const
 KS: Prepare Settings Info within SplineFile. More...
 
void PrepareMonolithDir (std::unique_ptr< TFile > &SplineFile) const
 KS: Prepare Monolith Info within SplineFile. More...
 
void PrepareIndexDir (std::unique_ptr< TFile > &SplineFile) const
 KS: Prepare Index Info within SplineFile. More...
 
void PrepareOtherInfoDir (std::unique_ptr< TFile > &SplineFile) const
 KS: Prepare Other Info within SplineFile. More...
 

Detailed Description

Bin-by-bin class calculating response for spline parameters.

See also
For more details, visit the Wiki.
Author
Dan Barrow
Ed Atkin
Henry Wallace

Definition at line 12 of file BinnedSplineHandler.h.

Member Enumeration Documentation

◆ TokenOrdering

Enumerator
kSystToken 
kModeToken 
kVar1BinToken 
kVar2BinToken 
kVar3BinToken 
kNTokens 

Definition at line 138 of file BinnedSplineHandler.h.

Constructor & Destructor Documentation

◆ BinnedSplineHandler()

_MaCh3_Safe_Include_Start_ _MaCh3_Safe_Include_End_ BinnedSplineHandler::BinnedSplineHandler ( ParameterHandlerGeneric xsec_,
MaCh3Modes Modes_ 
)

Constructor.

Todo:
ETA - do all of these functions and members actually need to be public?

Definition at line 14 of file BinnedSplineHandler.cpp.

14  : SplineBase() {
15 //****************************************
16  if (!xsec_) {
17  MACH3LOG_ERROR("Trying to create BinnedSplineHandler with uninitialized covariance object");
18  throw MaCh3Exception(__FILE__, __LINE__);
19  }
20  xsec = xsec_;
21 
22  if (!Modes_) {
23  MACH3LOG_ERROR("Trying to create BinnedSplineHandler with uninitialized MaCh3Modes object");
24  throw MaCh3Exception(__FILE__, __LINE__);
25  }
26  Modes = Modes_;
27 
28  // Keep these in class scope, important for using 1 monolith/sample!
29  MonolithIndex = 0; //Keeps track of the monolith index we're on when filling arrays (declared here so we can have multiple FillSampleArray calls)
30  CoeffIndex = 0; //Keeps track of our indexing the coefficient arrays [x, ybcd]
31  isflatarray = nullptr;
32 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
bool * isflatarray
Need to keep track of which splines are flat and which aren't.
ParameterHandlerGeneric * xsec
Pointer to covariance from which we get information about spline params.
MaCh3Modes * Modes
pointer to MaCh3 Mode from which we get spline suffix
Custom exception class for MaCh3 errors.
SplineBase()
Constructor.
Definition: SplineBase.cpp:7

◆ ~BinnedSplineHandler()

BinnedSplineHandler::~BinnedSplineHandler ( )
virtual

Destructor.

Todo:
it need some love

Definition at line 35 of file BinnedSplineHandler.cpp.

35  {
36 //****************************************
37  if(manycoeff_arr != nullptr) delete[] manycoeff_arr;
38  if(xcoeff_arr != nullptr) delete[] xcoeff_arr;
39  if(SplineSegments != nullptr) delete[] SplineSegments;
40  if(ParamValues != nullptr) delete[] ParamValues;
41 }
M3::float_t * xcoeff_arr
x coefficients for each spline
M3::float_t * manycoeff_arr
ybcd coefficients for each spline
short int * SplineSegments
Definition: SplineBase.h:74
float * ParamValues
Store parameter values they are not in FastSplineInfo as in case of GPU we need to copy paste it to G...
Definition: SplineBase.h:76

Member Function Documentation

◆ AddSample()

void BinnedSplineHandler::AddSample ( const std::string &  SampleName,
const std::string &  SampleTitle,
const std::vector< std::string > &  OscChanFileNames,
const std::vector< std::string > &  SplineVarNames 
)

add oscillation channel to spline monolith

Definition at line 67 of file BinnedSplineHandler.cpp.

73 {
74  SampleNames.push_back(SampleName);
75  SampleTitles.push_back(SampleTitle);
76  Dimensions.push_back(static_cast<int>(SplineVarNames.size()));
77  DimensionLabels.push_back(SplineVarNames);
78 
79  int nSplineParam = xsec->GetNumParamsFromSampleName(SampleName, SystType::kSpline);
80  nSplineParams.push_back(nSplineParam);
81 
82  //This holds the global index of the spline i.e. 0 -> _fNumPar
83  std::vector<int> GlobalSystIndex_Sample = xsec->GetGlobalSystIndexFromSampleName(SampleName, SystType::kSpline);
84  //Keep track of this for all the samples
85  GlobalSystIndex.push_back(GlobalSystIndex_Sample);
86 
87  std::vector<SplineInterpolation> SplineInterpolation_Sample = xsec->GetSplineInterpolationFromSampleName(SampleName);
88  // Keep track of this for all samples
89  SplineInterpolationTypes.push_back(SplineInterpolation_Sample);
90 
91  std::vector<std::string> SplineFileParPrefixNames_Sample = xsec->GetSplineParsNamesFromSampleName(SampleName);
92  SplineFileParPrefixNames.push_back(SplineFileParPrefixNames_Sample);
93 
94  MACH3LOG_INFO("Create SplineModeVecs_Sample");
95  std::vector<std::vector<int>> SplineModeVecs_Sample = StripDuplicatedModes(xsec->GetSplineModeVecFromSampleName(SampleName));
96  MACH3LOG_INFO("SplineModeVecs_Sample is of size {}", SplineModeVecs_Sample.size());
97  SplineModeVecs.push_back(SplineModeVecs_Sample);
98 
99  MACH3LOG_INFO("SplineModeVecs is of size {}", SplineModeVecs.size());
100 
101  int nOscChan = int(OscChanFileNames.size());
102  nOscChans.push_back(nOscChan);
103 
104  PrintSampleDetails(SampleTitle);
105 
106  std::vector<std::vector<TAxis *>> SampleBinning(nOscChan);
107  for (int iOscChan = 0; iOscChan < nOscChan; iOscChan++)
108  {
109  SampleBinning[iOscChan] = FindSplineBinning(OscChanFileNames[iOscChan], SampleTitle);
110  }
111  MACH3LOG_INFO("#----------------------------------------------------------------------------------------------------------------------------------#");
112  SplineBinning.push_back(SampleBinning);
113 
114  BuildSampleIndexingArray(SampleTitle);
115  PrintArrayDetails(SampleTitle);
116  MACH3LOG_INFO("#----------------------------------------------------------------------------------------------------------------------------------#");
117 
118  FillSampleArray(SampleTitle, OscChanFileNames);
119  MACH3LOG_INFO("#----------------------------------------------------------------------------------------------------------------------------------#");
120 }
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:25
@ kSpline
For splined parameters (1D)
std::vector< std::vector< std::string > > SplineFileParPrefixNames
std::vector< int > nSplineParams
std::vector< std::vector< std::vector< int > > > SplineModeVecs
std::vector< int > nOscChans
std::vector< std::string > SampleTitles
void PrintSampleDetails(const std::string &SampleTitle) const
Print info like Sample ID of spline params etc.
std::vector< std::vector< SplineInterpolation > > SplineInterpolationTypes
spline interpolation types for each sample. These vectors are from a call to GetSplineInterpolationFr...
std::vector< int > Dimensions
void BuildSampleIndexingArray(const std::string &SampleTitle)
std::vector< std::vector< int > > StripDuplicatedModes(const std::vector< std::vector< int > > &InputVector)
Check if there are any repeated modes. This is used to reduce the number of modes in case many intera...
virtual void FillSampleArray(std::string SampleTitle, std::vector< std::string > OscChanFileNames)
Loads and processes splines from ROOT files for a given sample.
std::vector< TAxis * > FindSplineBinning(const std::string &FileName, const std::string &SampleTitle)
Grab histograms with spline binning.
std::vector< std::vector< std::string > > DimensionLabels
std::vector< std::vector< std::vector< TAxis * > > > SplineBinning
void PrintArrayDetails(const std::string &SampleTitle) const
std::vector< std::vector< int > > GlobalSystIndex
This holds the global spline index and is used to grab the current parameter value to evaluate spline...
std::vector< std::string > SampleNames
const std::vector< int > GetGlobalSystIndexFromSampleName(const std::string &SampleName, const SystType Type)
DB Get spline parameters depending on given SampleName.
int GetNumParamsFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the number of parameters for the relevant SampleName.
const std::vector< std::string > GetSplineParsNamesFromSampleName(const std::string &SampleName)
DB Get spline parameters depending on given SampleName.
const std::vector< SplineInterpolation > GetSplineInterpolationFromSampleName(const std::string &SampleName)
Get the interpolation types for splines affecting a particular SampleName.
const std::vector< std::vector< int > > GetSplineModeVecFromSampleName(const std::string &SampleName)
DB Grab the Spline Modes for the relevant SampleName.

◆ BuildSampleIndexingArray()

void BinnedSplineHandler::BuildSampleIndexingArray ( const std::string &  SampleTitle)

Definition at line 350 of file BinnedSplineHandler.cpp.

352 {
353  int iSample = GetSampleIndex(SampleTitle);
354  int nSplineSysts = nSplineParams[iSample];
355  int nOscChannels = nOscChans[iSample];
356 
357  // Resize the main indexing structure
358  indexvec.emplace_back(nOscChannels);
359 
360  for (int iOscChan = 0; iOscChan < nOscChannels; ++iOscChan)
361  {
362  indexvec.back()[iOscChan].resize(nSplineSysts);
363  for (int iSyst = 0; iSyst < nSplineSysts; ++iSyst)
364  {
365  int nModesInSyst = static_cast<int>(SplineModeVecs[iSample][iSyst].size());
366  indexvec.back()[iOscChan][iSyst].resize(nModesInSyst);
367 
368  for (int iMode = 0; iMode < nModesInSyst; ++iMode)
369  {
370  const int nBins1 = SplineBinning[iSample][iOscChan][0]->GetNbins();
371  const int nBins2 = SplineBinning[iSample][iOscChan][1]->GetNbins();
372  const int nBins3 = SplineBinning[iSample][iOscChan][2]->GetNbins();
373 
374  indexvec.back()[iOscChan][iSyst][iMode]
375  .resize(nBins1,std::vector<std::vector<int>>(nBins2, std::vector<int>(nBins3, 0)));
376  }
377  } // end of iSyst loop
378  } // end of iOscChan loop
379 }
std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< int > > > > > > > indexvec
Variables related to determined which modes have splines and which piggy-back of other modes.
int GetSampleIndex(const std::string &SampleTitle) const
Get index of sample based on name.

◆ CalcSplineWeights()

void BinnedSplineHandler::CalcSplineWeights ( )
overrideprotectedvirtual

CPU based code which eval weight for each spline.

Todo:
KS: Once could use "ParamValues" but this will result in tiny bit different results due to floating point precision

Implements SplineBase.

Definition at line 309 of file BinnedSplineHandler.cpp.

311 {
312  #ifdef MULTITHREAD
313  #pragma omp parallel for simd
314  #endif
315  for (size_t iCoeff = 0; iCoeff < uniquecoeffindices.size(); ++iCoeff)
316  {
317  const int iSpline = uniquecoeffindices[iCoeff];
318  const short int uniqueIndex = short(uniquesplinevec_Monolith[iSpline]);
319  const short int currentsegment = short(SplineSegments[uniqueIndex]);
320 
321  const int segCoeff = coeffindexvec[iSpline]+currentsegment;
322  const int coeffOffset = segCoeff * _nCoeff_;
323  // These are what we can extract from the TSpline3
324  const M3::float_t y = manycoeff_arr[coeffOffset+kCoeffY];
325  const M3::float_t b = manycoeff_arr[coeffOffset+kCoeffB];
326  const M3::float_t c = manycoeff_arr[coeffOffset+kCoeffC];
327  const M3::float_t d = manycoeff_arr[coeffOffset+kCoeffD];
328 
329  // Get the variation for this reconfigure for the ith parameter
331  const M3::float_t xvar = (*SplineInfoArray[uniqueIndex].splineParsPointer);
332  // The Delta(x) = xvar - x
333  const M3::float_t dx = xvar - xcoeff_arr[segCoeff];
334 
335  //Speedy 1% time boost https://en.cppreference.com/w/c/numeric/math/fma (see ND code!)
336  M3::float_t weight = M3::fmaf_t(dx, M3::fmaf_t(dx, M3::fmaf_t(dx, d, c), b), y);
337  //This is the speedy version of writing dx^3+b*dx^2+c*dx+d
338 
339  //ETA - do we need this? We check later for negative weights and I wonder if this is even
340  //possible with the fmaf line above?
341  if(weight < 0){weight = 0.;} //Stops is getting negative weights
342 
343  weightvec_Monolith[iSpline] = weight;
344  }
345 }
constexpr int _nCoeff_
KS: We store coefficients {y,b,c,d} in one array one by one, this is only to define it once rather th...
Definition: SplineCommon.h:13
@ kCoeffB
Coefficient B.
Definition: SplineCommon.h:21
@ kCoeffD
Coefficient D.
Definition: SplineCommon.h:23
@ kCoeffY
Coefficient Y.
Definition: SplineCommon.h:20
@ kCoeffC
Coefficient C.
Definition: SplineCommon.h:22
std::vector< int > coeffindexvec
std::vector< int > uniquesplinevec_Monolith
Maps single spline object with single parameter.
std::vector< M3::float_t > weightvec_Monolith
Stores weight from spline evaluation for each single spline.
std::vector< int > uniquecoeffindices
Unique coefficient indices.
std::vector< FastSplineInfo > SplineInfoArray
Definition: SplineBase.h:71
double float_t
Definition: Core.h:30
constexpr T fmaf_t(T x, T y, T z)
Function template for fused multiply-add.
Definition: Core.h:38

◆ cleanUpMemory()

void BinnedSplineHandler::cleanUpMemory ( )

Remove setup variables not needed for spline evaluations.

Definition at line 43 of file BinnedSplineHandler.cpp.

43  {
44 //****************************************
45  //Call once everything's been allocated in SampleHandlerFDBase, cleans up junk from memory!
46  //Not a huge saving but it's better than leaving everything up to the compiler
47  MACH3LOG_INFO("Cleaning up spline memory");
63  if(isflatarray) delete [] isflatarray;
64 }
void CleanVector(std::vector< T > &vec)
Generic cleanup function.
void CleanContainer(std::vector< T * > &container)
Generic cleanup function.
std::vector< std::string > UniqueSystNames
name of each spline parameter
std::vector< TSpline3_red * > splinevec_Monolith
holds each spline object before stripping into coefficient monolith
std::vector< int > UniqueSystIndices
Global index of each spline param, it allows us to match spline ordering with global.

◆ CountNumberOfLoadedSplines()

int BinnedSplineHandler::CountNumberOfLoadedSplines ( bool  NonFlat = false,
int  Verbosity = 0 
)

Definition at line 485 of file BinnedSplineHandler.cpp.

487 {
488  int SampleCounter_NonFlat = 0;
489  int SampleCounter_All = 0;
490  int FullCounter_NonFlat = 0;
491  int FullCounter_All = 0;
492 
493  for (unsigned int iSample = 0; iSample < indexvec.size(); iSample++)
494  { // Loop over systematics
495  SampleCounter_NonFlat = 0;
496  SampleCounter_All = 0;
497  std::string SampleTitle = SampleTitles[iSample];
498  for (unsigned int iOscChan = 0; iOscChan < indexvec[iSample].size(); iOscChan++)
499  { // Loop over oscillation channels
500  for (unsigned int iSyst = 0; iSyst < indexvec[iSample][iOscChan].size(); iSyst++)
501  { // Loop over systematics
502  for (unsigned int iMode = 0; iMode < indexvec[iSample][iOscChan][iSyst].size(); iMode++)
503  { // Loop over modes
504  for (unsigned int iVar1 = 0; iVar1 < indexvec[iSample][iOscChan][iSyst][iMode].size(); iVar1++)
505  { // Loop over first dimension
506  for (unsigned int iVar2 = 0; iVar2 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1].size(); iVar2++)
507  { // Loop over second dimension
508  for (unsigned int iVar3 = 0; iVar3 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2].size(); iVar3++)
509  { // Loop over third dimension
510  if (isValidSplineIndex(SampleTitle, iOscChan, iSyst, iMode, iVar1, iVar2, iVar3))
511  {
512  int splineindex = indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2][iVar3];
513  if (splinevec_Monolith[splineindex])
514  {
515  SampleCounter_NonFlat += 1;
516  }
517  SampleCounter_All += 1;
518  }
519  }
520  }
521  }
522  }
523  }
524  }
525  MACH3LOG_DEBUG("{:<10} has {:<10} splines, of which {:<10} are not flat", SampleTitles[iSample], SampleCounter_All, SampleCounter_NonFlat);
526 
527  FullCounter_NonFlat += SampleCounter_NonFlat;
528  FullCounter_All += SampleCounter_All;
529  }
530 
531  if (Verbosity > 0)
532  {
533  MACH3LOG_INFO("Total number of splines loaded: {}", FullCounter_All);
534  MACH3LOG_INFO("Total number of non-flat splines loaded: {}", FullCounter_NonFlat);
535  }
536 
537  if (NonFlat) {
538  return FullCounter_NonFlat;
539  } else {
540  return FullCounter_All;
541  }
542 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:24
bool isValidSplineIndex(const std::string &SampleTitle, int iSyst, int iOscChan, int iMode, int iVar1, int iVar2, int iVar3)
Ensure we have spline for a given bin.

◆ Evaluate()

void BinnedSplineHandler::Evaluate ( )
overridevirtual

CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; probably the best one here! Same thing but pass parameter spline segments instead of variations.

Implements SplineBase.

Definition at line 295 of file BinnedSplineHandler.cpp.

295  {
296 // *****************************************
297  // There's a parameter mapping that goes from spline parameter to a global parameter index
298  // Find the spline segments
300 
301  //KS: Huge MP loop over all valid splines
303 
304  //KS: Huge MP loop over all events calculating total weight
305  ModifyWeights();
306 }
void CalcSplineWeights() override
CPU based code which eval weight for each spline.
void ModifyWeights() override
Calc total event weight, not used by Bin-by-bin splines.
void FindSplineSegment()
CW:Code used in step by step reweighting, Find Spline Segment for each param.
Definition: SplineBase.cpp:24

◆ FillSampleArray()

void BinnedSplineHandler::FillSampleArray ( std::string  SampleTitle,
std::vector< std::string >  OscChanFileNames 
)
virtual

Loads and processes splines from ROOT files for a given sample.

Note
DB Add virtual so it can be overridden in experiment specific (if needed)

Definition at line 927 of file BinnedSplineHandler.cpp.

928 {
929  int iSample = GetSampleIndex(SampleTitle);
930  int nOscChannels = nOscChans[iSample];
931 
932  for (int iOscChan = 0; iOscChan < nOscChannels; iOscChan++) {
933  MACH3LOG_INFO("Processing: {}", OscChanFileNames[iOscChan]);
934  TSpline3* mySpline = nullptr;
935  TSpline3_red* Spline = nullptr;
936  TString Syst, Mode;
937  int nKnots, SystNum, ModeNum, Var1Bin, Var2Bin, Var3Bin = M3::_BAD_INT_;
938  double x,y = M3::_BAD_DOUBLE_;
939  bool isFlat = true;
940 
941  std::set<std::string> SplineFileNames;
942 
943  auto File = std::unique_ptr<TFile>(TFile::Open(OscChanFileNames[iOscChan].c_str()));
944 
945  if (!File || File->IsZombie()) {
946  MACH3LOG_ERROR("File {} not found", OscChanFileNames[iOscChan]);
947  throw MaCh3Exception(__FILE__, __LINE__);
948  }
949 
950  //This is the MC specific part of the code
951  //i.e. we always assume that the splines are just store in single TDirectory and they're all in there as single objects
952  for (auto k : *File->GetListOfKeys()) {
953  auto Key = static_cast<TKey*>(k);
954  TClass *Class = gROOT->GetClass(Key->GetClassName(), false);
955  if(!Class->InheritsFrom("TSpline3")) {
956  continue;
957  }
958 
959  std::string FullSplineName = std::string(Key->GetName());
960 
961  if (SplineFileNames.count(FullSplineName) > 0) {
962  MACH3LOG_CRITICAL("Skipping spline - Found a spline whose name has already been encountered before: {}", FullSplineName);
963  continue;
964  }
965  SplineFileNames.insert(FullSplineName);
966 
967  std::vector<std::string> Tokens = GetTokensFromSplineName(FullSplineName);
968 
969  if (Tokens.size() != kNTokens) {
970  MACH3LOG_ERROR("Invalid tokens from spline name - Expected {} tokens. Check implementation in GetTokensFromSplineName()", static_cast<int>(kNTokens));
971  throw MaCh3Exception(__FILE__, __LINE__);
972  }
973 
974  Syst = Tokens[kSystToken];
975  Mode = Tokens[kModeToken];
976  Var1Bin = std::stoi(Tokens[kVar1BinToken]);
977  Var2Bin = std::stoi(Tokens[kVar2BinToken]);
978  Var3Bin = std::stoi(Tokens[kVar3BinToken]);
979 
980  SystNum = -1;
981  for (unsigned iSyst = 0; iSyst < SplineFileParPrefixNames[iSample].size(); iSyst++) {
982  if (Syst == SplineFileParPrefixNames[iSample][iSyst]) {
983  SystNum = iSyst;
984  break;
985  }
986  }
987 
988  // If the syst doesn't match any of the spline names then skip it
989  if (SystNum == -1){
990  MACH3LOG_DEBUG("Couldn't match!!");
991  MACH3LOG_DEBUG("Couldn't Match any systematic name in ParameterHandler with spline name: {}" , FullSplineName);
992  continue;
993  }
994 
995  ModeNum = -1;
996  for (unsigned int iMode = 0; iMode < SplineModeVecs[iSample][SystNum].size(); iMode++) {
997  if (Mode == Modes->GetSplineSuffixFromMaCh3Mode(SplineModeVecs[iSample][SystNum][iMode])) {
998  ModeNum = iMode;
999  break;
1000  }
1001  }
1002 
1003  if (ModeNum == -1) {
1004  //DB - If you have splines in the root file that you don't want to use (e.g. removing a mode from a syst), this will cause a throw
1005  // Therefore include as debug warning and continue instead
1006  MACH3LOG_DEBUG("Couldn't find mode for {} in {}. Problem Spline is : {} ", Mode, Syst, FullSplineName);
1007  continue;
1008  }
1009 
1010  mySpline = Key->ReadObject<TSpline3>();
1011 
1012  if (isValidSplineIndex(SampleTitle, iOscChan, SystNum, ModeNum, Var1Bin, Var2Bin, Var3Bin)) { // loop over all the spline knots and check their value
1013  MACH3LOG_DEBUG("Pushed back monolith for spline {}", FullSplineName);
1014  // if the value is 1 then set the flat bool to false
1015  nKnots = mySpline->GetNp();
1016  isFlat = true;
1017  for (int iKnot = 0; iKnot < nKnots; iKnot++) {
1018  mySpline->GetKnot(iKnot, x, y);
1019  if (y < 0.99999 || y > 1.00001)
1020  {
1021  isFlat = false;
1022  break;
1023  }
1024  }
1025 
1026  //Rather than keeping a mega vector of splines then converting, this should just keep everything nice in memory!
1027  indexvec[iSample][iOscChan][SystNum][ModeNum][Var1Bin][Var2Bin][Var3Bin]=MonolithIndex;
1028  coeffindexvec.push_back(CoeffIndex);
1029  // Should save memory rather saving [x_i_0 ,... x_i_maxknots] for every spline!
1030  if (isFlat) {
1031  splinevec_Monolith.push_back(nullptr);
1032  delete mySpline;
1033  } else {
1034  ApplyKnotWeightCapTSpline3(mySpline, SystNum, xsec);
1035  Spline = new TSpline3_red(mySpline, SplineInterpolationTypes[iSample][SystNum]);
1036  if(mySpline) delete mySpline;
1037 
1038  splinevec_Monolith.push_back(Spline);
1039  uniquecoeffindices.push_back(MonolithIndex); //So we can get the unique coefficients and skip flat splines later on!
1040  CoeffIndex+=nKnots;
1041  }
1042  //Incrementing MonolithIndex to keep track of number of valid spline indices
1043  MonolithIndex+=1;
1044  } else {
1045  //Potentially you are not a valid spline index
1046  delete mySpline;
1047  }
1048  }//End of loop over all TKeys in file
1049 
1050  //A bit of clean up
1051  File->Delete("*");
1052  File->Close();
1053  } //End of oscillation channel loop
1054 }
std::vector< bool > isFlat
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:28
void ApplyKnotWeightCapTSpline3(TSpline3 *&Spline, const int splineParsIndex, ParameterHandlerGeneric *XsecCov)
EM: Apply capping to knot weight for specified spline parameter. param graph needs to have been set i...
virtual std::vector< std::string > GetTokensFromSplineName(std::string FullSplineName)=0
std::string GetSplineSuffixFromMaCh3Mode(const int Index)
DB: Get binned spline mode suffic from MaCh3 Mode.
Definition: MaCh3Modes.cpp:243
CW: Reduced TSpline3 class.
constexpr static const double _BAD_DOUBLE_
Default value used for double initialisation.
Definition: Core.h:46
TFile * Open(const std::string &Name, const std::string &Type, const std::string &File, const int Line)
Opens a ROOT file with the given name and mode.
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:48

◆ FindSplineBinning()

std::vector< TAxis * > BinnedSplineHandler::FindSplineBinning ( const std::string &  FileName,
const std::string &  SampleTitle 
)

Grab histograms with spline binning.

Definition at line 382 of file BinnedSplineHandler.cpp.

384 {
385  int iSample = GetSampleIndex(SampleTitle);
386 
387  //Try declaring these outside of TFile so they aren't owned by File
388  constexpr int nDummyBins = 1;
389  constexpr double DummyEdges[2] = {-1e15, 1e15};
390  TAxis* DummyAxis = new TAxis(nDummyBins, DummyEdges);
391  TH2F* Hist2D = nullptr;
392  TH3F* Hist3D = nullptr;
393 
394  auto File = std::unique_ptr<TFile>(TFile::Open(FileName.c_str(), "READ"));
395  if (!File || File->IsZombie())
396  {
397  MACH3LOG_ERROR("File {} not found", FileName);
398  MACH3LOG_ERROR("This is caused by something here! {} : {}", __FILE__, __LINE__);
399  throw MaCh3Exception(__FILE__ , __LINE__ );
400  }
401 
402  MACH3LOG_INFO("Finding binning for:");
403  MACH3LOG_INFO("{}", FileName);
404 
405  std::string TemplateName = "dev_tmp_0_0";
406  TObject *Obj = File->Get(TemplateName.c_str());
407  //If you can't find dev_tmp_0_0 then this will cause a problem
408  if (!Obj)
409  {
410  TemplateName = "dev_tmp.0.0";
411  Obj = File->Get(TemplateName.c_str());
412  if (!Obj)
413  {
414  MACH3LOG_ERROR("Could not find dev_tmp_0_0 in spline file. Spline binning cannot be set!");
415  MACH3LOG_ERROR("FileName: {}", FileName);
416  throw MaCh3Exception(__FILE__ , __LINE__ );
417  }
418  }
419 
420  //Now check if dev_tmp_0_0 is a TH2 i.e. specifying the dimensions of the splines is 2D
421  bool isHist2D = Obj->IsA() == TH2F::Class();
422  //For T2K annoyingly all objects are TH3Fs
423  bool isHist3D = Obj->IsA() == TH3F::Class();
424  if (!isHist2D && !isHist3D)
425  {
426  MACH3LOG_ERROR("Object doesn't inherit from either TH2D and TH3D - Odd A");
427  throw MaCh3Exception(__FILE__ , __LINE__ );
428  }
429 
430  if (isHist2D)
431  {
432  if (Dimensions[iSample] != 2)
433  {
434  MACH3LOG_ERROR("Trying to load a 2D spline template when nDim={}", Dimensions[iSample]);
435  throw MaCh3Exception(__FILE__, __LINE__);
436  }
437  Hist2D = File->Get<TH2F>(TemplateName.c_str());
438  }
439 
440  if (isHist3D)
441  {
442  Hist3D = File->Get<TH3F>((TemplateName.c_str()));
443 
444  if (Dimensions[iSample] != 3 && Hist3D->GetZaxis()->GetNbins() != 1)
445  {
446  MACH3LOG_ERROR("Trying to load a 3D spline template when nDim={}", Dimensions[iSample]);
447  throw MaCh3Exception(__FILE__ , __LINE__ );
448  }
449  Hist3D = File->Get<TH3F>(TemplateName.c_str());
450  }
451 
452  std::vector<TAxis*> ReturnVec;
453  // KS: Resize to reduce impact of push back and memory fragmentation
454  ReturnVec.resize(3);
455  if (Dimensions[iSample] == 2) {
456  if (isHist2D) {
457  ReturnVec[0] = static_cast<TAxis*>(Hist2D->GetXaxis()->Clone());
458  ReturnVec[1] = static_cast<TAxis*>(Hist2D->GetYaxis()->Clone());
459  ReturnVec[2] = static_cast<TAxis*>(DummyAxis->Clone());
460  } else if (isHist3D) {
461  ReturnVec[0] = static_cast<TAxis*>(Hist3D->GetXaxis()->Clone());
462  ReturnVec[1] = static_cast<TAxis*>(Hist3D->GetYaxis()->Clone());
463  ReturnVec[2] = static_cast<TAxis*>(DummyAxis->Clone());
464  }
465  } else if (Dimensions[iSample] == 3) {
466  ReturnVec[0] = static_cast<TAxis*>(Hist3D->GetXaxis()->Clone());
467  ReturnVec[1] = static_cast<TAxis*>(Hist3D->GetYaxis()->Clone());
468  ReturnVec[2] = static_cast<TAxis*>(Hist3D->GetZaxis()->Clone());
469  } else {
470  MACH3LOG_ERROR("Number of dimensions not valid! Given: {}", Dimensions[iSample]);
471  throw MaCh3Exception(__FILE__, __LINE__);
472  }
473 
474  for (unsigned int iAxis = 0; iAxis < ReturnVec.size(); ++iAxis) {
475  PrintBinning(ReturnVec[iAxis]);
476  }
477 
478  MACH3LOG_INFO("Left PrintBinning now tidying up");
479  delete DummyAxis;
480 
481  return ReturnVec;
482 }
void PrintBinning(TAxis *Axis) const

◆ getDimLabel()

std::string BinnedSplineHandler::getDimLabel ( const int  BinningOpt,
const unsigned int  Axis 
) const

Definition at line 722 of file BinnedSplineHandler.cpp.

724 {
725  if(Axis > DimensionLabels[iSample].size()){
726  MACH3LOG_ERROR("The spline Axis you are trying to get the label of is larger than the number of dimensions");
727  MACH3LOG_ERROR("You are trying to get axis {} but have only got {}", Axis, Dimensions[iSample]);
728  throw MaCh3Exception(__FILE__ , __LINE__ );
729  }
730  return DimensionLabels.at(iSample).at(Axis);
731 }
int size

◆ GetEventSplines()

std::vector< std::vector< int > > BinnedSplineHandler::GetEventSplines ( const std::string &  SampleTitle,
int  iOscChan,
int  EventMode,
double  Var1Val,
double  Var2Val,
double  Var3Val 
)

Return the splines which affect a given event.

Definition at line 836 of file BinnedSplineHandler.cpp.

838 {
839  std::vector<std::vector<int>> ReturnVec;
840  int SampleIndex = GetSampleIndex(SampleTitle);
841  int nSplineSysts = static_cast<int>(indexvec[SampleIndex][iOscChan].size());
842 
843  int Mode = -1;
844  std::string SuffixForEventMode = Modes->GetSplineSuffixFromMaCh3Mode(EventMode);
845  for (int iMode=0;iMode<Modes->GetNModes();iMode++) {
846  if (SuffixForEventMode == Modes->GetSplineSuffixFromMaCh3Mode(iMode)) {
847  Mode = iMode;
848  break;
849  }
850  }
851  if (Mode == -1) {
852  return ReturnVec;
853  }
854 
855  int Var1Bin = SplineBinning[SampleIndex][iOscChan][0]->FindBin(Var1Val)-1;
856  if (Var1Bin < 0 || Var1Bin >= SplineBinning[SampleIndex][iOscChan][0]->GetNbins()) {
857  return ReturnVec;
858  }
859 
860  int Var2Bin = SplineBinning[SampleIndex][iOscChan][1]->FindBin(Var2Val)-1;
861  if (Var2Bin < 0 || Var2Bin >= SplineBinning[SampleIndex][iOscChan][1]->GetNbins()) {
862  return ReturnVec;
863  }
864 
865  int Var3Bin = SplineBinning[SampleIndex][iOscChan][2]->FindBin(Var3Val)-1;
866  if (Var3Bin < 0 || Var3Bin >= SplineBinning[SampleIndex][iOscChan][2]->GetNbins()){
867  return ReturnVec;
868  }
869 
870  for(int iSyst=0; iSyst<nSplineSysts; iSyst++){
871  std::vector<int> spline_modes = SplineModeVecs[SampleIndex][iSyst];
872  int nSampleModes = static_cast<int>(spline_modes.size());
873 
874  //ETA - look here at the length of spline_modes and what you're actually comparing against
875  for(int iMode = 0; iMode<nSampleModes ; iMode++){
876  //Only consider if the event mode (Mode) matches ones of the spline modes
877  if (Mode == spline_modes[iMode]) {
878  int splineID=indexvec[SampleIndex][iOscChan][iSyst][iMode][Var1Bin][Var2Bin][Var3Bin];
879  //Also check that the spline isn't flat
880  if(!isflatarray[splineID]){
881  ReturnVec.push_back({SampleIndex, iOscChan, iSyst, iMode, Var1Bin, Var2Bin, Var3Bin});
882  }
883  }
884  }
885  }
886 
887  return ReturnVec;
888 }
int GetNModes() const
KS: Get number of modes, keep in mind actual number is +1 greater due to unknown category.
Definition: MaCh3Modes.h:54

◆ GetSampleIndex()

int BinnedSplineHandler::GetSampleIndex ( const std::string &  SampleTitle) const

Get index of sample based on name.

Parameters
SampleTitleThe title of the sample to search for.

Definition at line 735 of file BinnedSplineHandler.cpp.

735  {
736 //****************************************
737  for (size_t iSample = 0; iSample < SampleTitles.size(); ++iSample) {
738  if (SampleTitle == SampleTitles[iSample]) {
739  return static_cast<int>(iSample);
740  }
741  }
742  MACH3LOG_ERROR("Sample name not found: {}", SampleTitle);
743  throw MaCh3Exception(__FILE__, __LINE__);
744 }

◆ getSplineCoeff_SepMany()

void BinnedSplineHandler::getSplineCoeff_SepMany ( int  splineindex,
M3::float_t *&  xArray,
M3::float_t *&  manyArray 
)

Definition at line 686 of file BinnedSplineHandler.cpp.

686  {
687 //****************************************
688  //No point evaluating a flat spline
689  int nPoints = splinevec_Monolith[splineindex]->GetNp();
690 
691  for (int i = 0; i < nPoints; i++) {
692  xArray[i] = 1.0;
693  for (int j = 0; j < _nCoeff_; j++) {
694  manyArray[i*_nCoeff_+j] = 1.0;
695  }
696  }
697 
698  for(int i=0; i<nPoints; i++) {
699  M3::float_t x = M3::float_t(-999.99);
700  M3::float_t y = M3::float_t(-999.99);
701  M3::float_t b = M3::float_t(-999.99);
702  M3::float_t c = M3::float_t(-999.99);
703  M3::float_t d = M3::float_t(-999.99);
704  splinevec_Monolith[splineindex]->GetCoeff(i, x, y, b, c, d);
705 
706  // Store the coefficients for each knot contiguously in memory
707  // 4 because manyArray stores y,b,c,d
708  xArray[i] = x;
709  manyArray[i * _nCoeff_ + kCoeffY] = y;
710  manyArray[i * _nCoeff_ + kCoeffB] = b;
711  manyArray[i * _nCoeff_ + kCoeffC] = c;
712  manyArray[i * _nCoeff_ + kCoeffD] = d;
713  }
714 
715  //We now clean up the splines!
716  delete splinevec_Monolith[splineindex];
717  splinevec_Monolith[splineindex] = nullptr;
718 }

◆ GetTokensFromSplineName()

virtual std::vector<std::string> BinnedSplineHandler::GetTokensFromSplineName ( std::string  FullSplineName)
protectedpure virtual

◆ InvestigateMissingSplines()

void BinnedSplineHandler::InvestigateMissingSplines ( ) const
private

This function will find missing splines in file.

Definition at line 123 of file BinnedSplineHandler.cpp.

123  {
124 //****************************************
125  // Map: iSample → iSyst → modeSuffix → {totalSplines, zeroCount}
126  std::map<unsigned int, std::map<unsigned int, std::map<std::string, std::pair<unsigned int, unsigned int>>>> systZeroCounts;
127 
128  for (unsigned int iSample = 0; iSample < indexvec.size(); iSample++) {
129  std::string SampleName = SampleNames[iSample];
130 
131  // Get list of systematic names for this sample
132  std::vector<std::string> SplineFileParPrefixNames_Sample =
134 
135  for (unsigned int iOscChan = 0; iOscChan < indexvec[iSample].size(); iOscChan++) {
136  for (unsigned int iSyst = 0; iSyst < indexvec[iSample][iOscChan].size(); iSyst++) {
137  unsigned int zeroCount = 0;
138  unsigned int totalSplines = 0;
139 
140  for (unsigned int iMode = 0; iMode < indexvec[iSample][iOscChan][iSyst].size(); iMode++) {
141  // Get the mode suffix string
142  std::string modeSuffix =
143  Modes->GetSplineSuffixFromMaCh3Mode(SplineModeVecs[iSample][iSyst][iMode]);
144 
145  for (unsigned int iVar1 = 0; iVar1 < indexvec[iSample][iOscChan][iSyst][iMode].size(); iVar1++) {
146  for (unsigned int iVar2 = 0; iVar2 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1].size(); iVar2++) {
147  for (unsigned int iVar3 = 0; iVar3 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2].size(); iVar3++) {
148  totalSplines++;
149  if (indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2][iVar3] == 0) {
150  zeroCount++;
151  if(zeroCount > 1){
152  systZeroCounts[iSample][iSyst][modeSuffix] = {totalSplines, zeroCount};
153  }
155  "Sample '{}' | OscChan {} | Syst '{}' | Mode '{}' | Var1 {} | Var2 {} | Var3 {} => Value: {}",
156  SampleTitles[iSample],
157  iOscChan,
158  SplineFileParPrefixNames_Sample[iSyst],
159  modeSuffix,
160  iVar1,
161  iVar2,
162  iVar3,
163  indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2][iVar3]
164  );
165  }
166  }
167  }
168  }
169  }
170  }
171  }
172  }
173 
174  // KS: Let's print this atrocious mess...
175  for (const auto& samplePair : systZeroCounts) {
176  unsigned int iSample = samplePair.first;
177  std::vector<std::string> SplineFileParPrefixNames_Sample = xsec->GetParsNamesFromSampleName(SampleNames[iSample], kSpline);
178  for (const auto& systPair : samplePair.second) {
179  unsigned int iSyst = systPair.first;
180  const auto& systName = SplineFileParPrefixNames_Sample[iSyst];
181  for (const auto& modePair : systPair.second) {
182  const auto& modeSuffix = modePair.first;
183  const auto& counts = modePair.second;
185  "Sample '{}': Systematic '{}' has missing splines in mode '{}'. Expected Splines: {}, Missing Splines: {}",
186  SampleTitles[iSample],
187  systName,
188  modeSuffix,
189  counts.first,
190  counts.second
191  );
192  }
193  }
194  }
195 }
const std::vector< std::string > GetParsNamesFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the parameter names for the relevant SampleName.

◆ isValidSplineIndex()

bool BinnedSplineHandler::isValidSplineIndex ( const std::string &  SampleTitle,
int  iSyst,
int  iOscChan,
int  iMode,
int  iVar1,
int  iVar2,
int  iVar3 
)

Ensure we have spline for a given bin.

Definition at line 785 of file BinnedSplineHandler.cpp.

787 {
788  int iSample = GetSampleIndex(SampleTitle);
789  bool isValid = true;
790 
791  // Lambda to check if an index is valid for a specific dimension
792  auto checkIndex = [&isValid](int index, size_t size, const std::string& name) {
793  if (index < 0 || index >= int(size)) {
794  MACH3LOG_ERROR("{} index is invalid! 0 <= Index < {} ", name, size);
795  isValid = false;
796  }
797  };
798 
799  checkIndex(iSample, indexvec.size(), "Sample");
800  if (isValid) checkIndex(iOscChan, indexvec[iSample].size(), "OscChan");
801  if (isValid) checkIndex(iSyst, indexvec[iSample][iOscChan].size(), "Syst");
802  if (isValid) checkIndex(iMode, indexvec[iSample][iOscChan][iSyst].size(), "Mode");
803  if (isValid) checkIndex(iVar1, indexvec[iSample][iOscChan][iSyst][iMode].size(), "Var1");
804  if (isValid) checkIndex(iVar2, indexvec[iSample][iOscChan][iSyst][iMode][iVar1].size(), "Var2");
805  if (isValid) checkIndex(iVar3, indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2].size(), "Var3");
806 
807  if (!isValid)
808  {
809  MACH3LOG_ERROR("Given iSample: {}", iSample);
810  MACH3LOG_ERROR("Given iOscChan: {}", iOscChan);
811  MACH3LOG_ERROR("Given iSyst: {}", iSyst);
812  MACH3LOG_ERROR("Given iMode: {}", iMode);
813  MACH3LOG_ERROR("Given iVar1: {}", iVar1);
814  MACH3LOG_ERROR("Given iVar2: {}", iVar2);
815  MACH3LOG_ERROR("Given iVar3: {}", iVar3);
816  MACH3LOG_ERROR("Come visit me at : {} : {}", __FILE__, __LINE__);
817  throw MaCh3Exception(__FILE__, __LINE__);
818  }
819 
820  return true;
821 }

◆ LoadIndexDir()

void BinnedSplineHandler::LoadIndexDir ( std::unique_ptr< TFile > &  SplineFile)
private

KS: Load preprocessed Index.

Parameters
FileFile from which we load new tree

Definition at line 1192 of file BinnedSplineHandler.cpp.

1192  {
1193 // *****************************************
1194  TTree *IndexTree = SplineFile->Get<TTree>("IndexVec");
1195 
1196  std::vector<int> Dim(7);
1197  int value;
1198  for (int d = 0; d < 7; ++d) {
1199  IndexTree->SetBranchAddress(Form("dim%d", d+1), &Dim[d]);
1200  }
1201  IndexTree->SetBranchAddress("value", &value);
1202 
1203  // Fill indexvec with data from IndexTree
1204  for (Long64_t i = 0; i < IndexTree->GetEntries(); ++i) {
1205  IndexTree->GetEntry(i);
1206  indexvec[Dim[0]][Dim[1]][Dim[2]][Dim[3]][Dim[4]][Dim[5]][Dim[6]] = value;
1207  }
1208 
1209  // Load SplineBinning data
1210  TTree *SplineBinningTree = SplineFile->Get<TTree>("SplineBinningTree");
1211  std::vector<int> indices(3);
1212  SplineBinningTree->SetBranchAddress("i", &indices[0]);
1213  SplineBinningTree->SetBranchAddress("j", &indices[1]);
1214  SplineBinningTree->SetBranchAddress("k", &indices[2]);
1215  TAxis* axis = nullptr;
1216  SplineBinningTree->SetBranchAddress("axis", &axis);
1217 
1218  // Reconstruct TAxis objects
1219  for (Long64_t entry = 0; entry < SplineBinningTree->GetEntries(); ++entry) {
1220  SplineBinningTree->GetEntry(entry);
1221  int i = indices[0];
1222  int j = indices[1];
1223  int k = indices[2];
1224  SplineBinning[i][j][k] = static_cast<TAxis*>(axis->Clone());
1225  }
1226 
1227  std::vector<int> indices_mode(3);
1228  int mode_value;
1229  TTree *SplineModeTree = SplineFile->Get<TTree>("SplineModeTree");
1230  SplineModeTree->SetBranchAddress("i", &indices_mode[0]);
1231  SplineModeTree->SetBranchAddress("j", &indices_mode[1]);
1232  SplineModeTree->SetBranchAddress("k", &indices_mode[2]);
1233  SplineModeTree->SetBranchAddress("value", &mode_value);
1234 
1235  // Fill SplineModeVecs with values from the tree
1236  for (Long64_t entry = 0; entry < SplineModeTree->GetEntries(); ++entry) {
1237  SplineModeTree->GetEntry(entry);
1238  int i = indices_mode[0];
1239  int j = indices_mode[1];
1240  int k = indices_mode[2];
1241  SplineModeVecs[i][j][k] = mode_value;
1242  }
1243 }

◆ LoadMonolithDir()

void BinnedSplineHandler::LoadMonolithDir ( std::unique_ptr< TFile > &  SplineFile)
private

KS: Load preprocessed Monolith.

Parameters
FileFile from which we load new tree

Definition at line 1158 of file BinnedSplineHandler.cpp.

1158  {
1159 // *****************************************
1160  TTree *MonolithTree = SplineFile->Get<TTree>("MonolithTree");
1161 
1163  MonolithTree->SetBranchAddress("manycoeff", manycoeff_arr);
1164  isflatarray = new bool[MonolithSize];
1166  MonolithTree->SetBranchAddress("isflatarray", isflatarray);
1167 
1168  // Load vectors
1169  std::vector<int>* coeffindexvec_temp = nullptr;
1170  MonolithTree->SetBranchAddress("coeffindexvec", &coeffindexvec_temp);
1171  std::vector<int>* uniquecoeffindices_temp = nullptr;
1172  MonolithTree->SetBranchAddress("uniquecoeffindices", &uniquecoeffindices_temp);
1173  std::vector<int>* uniquesplinevec_Monolith_temp = nullptr;
1174  MonolithTree->SetBranchAddress("uniquesplinevec_Monolith", &uniquesplinevec_Monolith_temp);
1175  std::vector<int>* UniqueSystIndices_temp = nullptr;
1176  MonolithTree->SetBranchAddress("UniqueSystIndices", &UniqueSystIndices_temp);
1177 
1178  // Allocate and load xcoeff_arr
1180  MonolithTree->SetBranchAddress("xcoeff", xcoeff_arr);
1181 
1182  MonolithTree->GetEntry(0);
1183 
1184  coeffindexvec = *coeffindexvec_temp;
1185  uniquecoeffindices = *uniquecoeffindices_temp;
1186  uniquesplinevec_Monolith = *uniquesplinevec_Monolith_temp;
1187  UniqueSystIndices = *UniqueSystIndices_temp;
1188 }

◆ LoadSettingsDir()

void BinnedSplineHandler::LoadSettingsDir ( std::unique_ptr< TFile > &  SplineFile)
private

KS: Load preprocessed Settings.

Parameters
FileFile from which we load new tree

Definition at line 1084 of file BinnedSplineHandler.cpp.

1084  {
1085 // *****************************************
1086  TTree *Settings = SplineFile->Get<TTree>("Settings");
1087  int CoeffIndex_temp, MonolithSize_temp;
1088  short int nParams_temp;
1089  Settings->SetBranchAddress("CoeffIndex", &CoeffIndex_temp);
1090  Settings->SetBranchAddress("MonolithSize", &MonolithSize_temp);
1091  Settings->SetBranchAddress("nParams", &nParams_temp);
1092  int indexvec_sizes[7];
1093  for (int i = 0; i < 7; ++i) {
1094  Settings->SetBranchAddress(("indexvec_size" + std::to_string(i+1)).c_str(), &indexvec_sizes[i]);
1095  }
1096 
1097  int SplineBinning_size1, SplineBinning_size2, SplineBinning_size3;
1098  Settings->SetBranchAddress("SplineBinning_size1", &SplineBinning_size1);
1099  Settings->SetBranchAddress("SplineBinning_size2", &SplineBinning_size2);
1100  Settings->SetBranchAddress("SplineBinning_size3", &SplineBinning_size3);
1101  int SplineModeVecs_size1, SplineModeVecs_size2, SplineModeVecs_size3;
1102  Settings->SetBranchAddress("SplineModeVecs_size1", &SplineModeVecs_size1);
1103  Settings->SetBranchAddress("SplineModeVecs_size2", &SplineModeVecs_size2);
1104  Settings->SetBranchAddress("SplineModeVecs_size3", &SplineModeVecs_size3);
1105  std::vector<std::string>* SampleNames_temp = nullptr;
1106  Settings->SetBranchAddress("SampleNames", &SampleNames_temp);
1107  std::vector<std::string>* SampleTitles_temp = nullptr;
1108  Settings->SetBranchAddress("SampleTitles", &SampleTitles_temp);
1109  Settings->GetEntry(0);
1110 
1111  CoeffIndex = CoeffIndex_temp;
1112  MonolithSize = MonolithSize_temp;
1113  SampleNames = *SampleNames_temp;
1114  SampleTitles = *SampleTitles_temp;
1115 
1116  nParams = nParams_temp;
1117 
1118  SplineSegments = new short int[nParams]();
1119  ParamValues = new float[nParams]();
1120 
1121  // Resize indexvec according to saved dimensions
1122  indexvec.resize(indexvec_sizes[0]);
1123  for (int i = 0; i < indexvec_sizes[0]; ++i) {
1124  indexvec[i].resize(indexvec_sizes[1]);
1125  for (int j = 0; j < indexvec_sizes[1]; ++j) {
1126  indexvec[i][j].resize(indexvec_sizes[2]);
1127  for (int k = 0; k < indexvec_sizes[2]; ++k) {
1128  indexvec[i][j][k].resize(indexvec_sizes[3]);
1129  for (int l = 0; l < indexvec_sizes[3]; ++l) {
1130  indexvec[i][j][k][l].resize(indexvec_sizes[4]);
1131  for (int m = 0; m < indexvec_sizes[4]; ++m) {
1132  indexvec[i][j][k][l][m].resize(indexvec_sizes[5]);
1133  for (int n = 0; n < indexvec_sizes[5]; ++n) {
1134  indexvec[i][j][k][l][m][n].resize(indexvec_sizes[6]);
1135  }
1136  }
1137  }
1138  }
1139  }
1140  }
1141 
1142  auto Resize3D = [](auto& vec, int d1, int d2, int d3) {
1143  vec.resize(d1);
1144  for (int i = 0; i < d1; ++i) {
1145  vec[i].resize(d2);
1146  for (int j = 0; j < d2; ++j) {
1147  vec[i][j].resize(d3);
1148  }
1149  }
1150  };
1151 
1152  Resize3D(SplineBinning, SplineBinning_size1, SplineBinning_size2, SplineBinning_size3);
1153  Resize3D(SplineModeVecs, SplineModeVecs_size1, SplineModeVecs_size2, SplineModeVecs_size3);
1154 }
short int nParams
Number of parameters that have splines.
Definition: SplineBase.h:78

◆ LoadSplineFile()

void BinnedSplineHandler::LoadSplineFile ( std::string  FileName)
overridevirtual

KS: Load preprocessed spline file.

Parameters
FileNamePath to ROOT file with predefined reduced Spline Monolith

Implements SplineBase.

Definition at line 1058 of file BinnedSplineHandler.cpp.

1058  {
1059 // *****************************************
1060  M3::AddPath(FileName);
1061 
1062  // Check for spaces in the filename
1063  size_t pos = FileName.find(' ');
1064  if (pos != std::string::npos) {
1065  MACH3LOG_WARN("Filename ({}) contains spaces. Replacing spaces with underscores.", FileName);
1066  while ((pos = FileName.find(' ')) != std::string::npos) {
1067  FileName[pos] = '_';
1068  }
1069  }
1070  auto SplineFile = std::make_unique<TFile>(FileName.c_str(), "OPEN");
1071  LoadSettingsDir(SplineFile);
1072  LoadMonolithDir(SplineFile);
1073  LoadIndexDir(SplineFile);
1074  LoadFastSplineInfoDir(SplineFile);
1075 
1076  for (int iSpline = 0; iSpline < nParams; iSpline++) {
1077  SplineInfoArray[iSpline].splineParsPointer = xsec->RetPointer(UniqueSystIndices[iSpline]);
1078  }
1079  SplineFile->Close();
1080 }
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:26
void LoadIndexDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Index.
void LoadSettingsDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Settings.
void LoadMonolithDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Monolith.
const double * RetPointer(const int iParam)
DB Pointer return to param position.
void LoadFastSplineInfoDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed FastSplineInfo.
Definition: SplineBase.cpp:159
void AddPath(std::string &FilePath)
Prepends the MACH3 environment path to FilePath if it is not already present.
Definition: Monitor.cpp:367

◆ ModifyWeights()

void BinnedSplineHandler::ModifyWeights ( )
inlineoverrideprotectedvirtual

Calc total event weight, not used by Bin-by-bin splines.

Implements SplineBase.

Definition at line 80 of file BinnedSplineHandler.h.

80 {return;};

◆ PrepareIndexDir()

void BinnedSplineHandler::PrepareIndexDir ( std::unique_ptr< TFile > &  SplineFile) const
private

KS: Prepare Index Info within SplineFile.

Definition at line 1346 of file BinnedSplineHandler.cpp.

1346  {
1347 // *****************************************
1348  // Create a TTree to store the data
1349  TTree *IndexTree = new TTree("IndexVec", "IndexVec");
1350 
1351  // Vector holding the 7 dims
1352  std::vector<int> Dim(7);
1353  int value;
1354 
1355  // Create branches for each dimension
1356  for (int d = 0; d < 7; ++d) {
1357  IndexTree->Branch(Form("dim%d", d+1), &Dim[d], Form("dim%d/I", d+1));
1358  }
1359  IndexTree->Branch("value", &value, "value/I");
1360 
1361  // Fill the tree
1362  for (size_t i = 0; i < indexvec.size(); ++i) {
1363  for (size_t j = 0; j < indexvec[i].size(); ++j) {
1364  for (size_t k = 0; k < indexvec[i][j].size(); ++k) {
1365  for (size_t l = 0; l < indexvec[i][j][k].size(); ++l) {
1366  for (size_t m = 0; m < indexvec[i][j][k][l].size(); ++m) {
1367  for (size_t n = 0; n < indexvec[i][j][k][l][m].size(); ++n) {
1368  for (size_t p = 0; p < indexvec[i][j][k][l][m][n].size(); ++p) {
1369  Dim[0] = static_cast<int>(i);
1370  Dim[1] = static_cast<int>(j);
1371  Dim[2] = static_cast<int>(k);
1372  Dim[3] = static_cast<int>(l);
1373  Dim[4] = static_cast<int>(m);
1374  Dim[5] = static_cast<int>(n);
1375  Dim[6] = static_cast<int>(p);
1376  value = static_cast<int>(indexvec[i][j][k][l][m][n][p]);
1377  IndexTree->Fill();
1378  }
1379  }
1380  }
1381  }
1382  }
1383  }
1384  }
1385 
1386  SplineFile->cd();
1387  // Write the tree to the file
1388  IndexTree->Write();
1389  delete IndexTree;
1390 }

◆ PrepareMonolithDir()

void BinnedSplineHandler::PrepareMonolithDir ( std::unique_ptr< TFile > &  SplineFile) const
private

KS: Prepare Monolith Info within SplineFile.

Definition at line 1323 of file BinnedSplineHandler.cpp.

1323  {
1324 // *****************************************
1325  TTree *MonolithTree = new TTree("MonolithTree", "MonolithTree");
1326  MonolithTree->Branch("manycoeff", manycoeff_arr, Form("manycoeff[%d]/%s", CoeffIndex * _nCoeff_, M3::float_t_str_repr));
1327  MonolithTree->Branch("isflatarray", isflatarray, Form("isflatarray[%d]/O", MonolithSize));
1328 
1329  std::vector<int> coeffindexvec_temp = coeffindexvec;
1330  MonolithTree->Branch("coeffindexvec", &coeffindexvec_temp);
1331  std::vector<int> uniquecoeffindices_temp = uniquecoeffindices;
1332  MonolithTree->Branch("uniquecoeffindices", &uniquecoeffindices_temp);
1333  std::vector<int> uniquesplinevec_Monolith_temp = uniquesplinevec_Monolith;
1334  MonolithTree->Branch("uniquesplinevec_Monolith", &uniquesplinevec_Monolith_temp);
1335  std::vector<int> UniqueSystIndices_temp = UniqueSystIndices;
1336  MonolithTree->Branch("UniqueSystIndices", &UniqueSystIndices_temp);
1337  MonolithTree->Branch("xcoeff", xcoeff_arr, Form("xcoeff[%d]/%s", CoeffIndex, M3::float_t_str_repr));
1338 
1339  MonolithTree->Fill();
1340  SplineFile->cd();
1341  MonolithTree->Write();
1342  delete MonolithTree;
1343 }
constexpr static const char * float_t_str_repr
Definition: Core.h:33

◆ PrepareOtherInfoDir()

void BinnedSplineHandler::PrepareOtherInfoDir ( std::unique_ptr< TFile > &  SplineFile) const
private

KS: Prepare Other Info within SplineFile.

Definition at line 1393 of file BinnedSplineHandler.cpp.

1393  {
1394 // *****************************************
1395  // Create a new tree for SplineBinning data
1396  TTree *SplineBinningTree = new TTree("SplineBinningTree", "SplineBinningTree");
1397  std::vector<int> indices(3); // To store the 3D indices
1398  TAxis* axis = nullptr;
1399  SplineBinningTree->Branch("i", &indices[0], "i/I");
1400  SplineBinningTree->Branch("j", &indices[1], "j/I");
1401  SplineBinningTree->Branch("k", &indices[2], "k/I");
1402  SplineBinningTree->Branch("axis", "TAxis", &axis);
1403 
1404  // Fill the SplineBinningTree
1405  for (size_t i = 0; i < SplineBinning.size(); ++i) {
1406  for (size_t j = 0; j < SplineBinning[i].size(); ++j) {
1407  for (size_t k = 0; k < SplineBinning[i][j].size(); ++k) {
1408  axis = SplineBinning[i][j][k];
1409  indices[0] = static_cast<int>(i);
1410  indices[1] = static_cast<int>(j);
1411  indices[2] = static_cast<int>(k);
1412  SplineBinningTree->Fill();
1413  }
1414  }
1415  }
1416  SplineFile->cd();
1417  SplineBinningTree->Write();
1418  delete SplineBinningTree;
1419 
1420  std::vector<int> indices_mode(3); // to store 3D indices
1421  int mode_value;
1422 
1423  TTree *SplineModeTree = new TTree("SplineModeTree", "SplineModeTree");
1424  // Create branches for indices and value
1425  SplineModeTree->Branch("i", &indices_mode[0], "i/I");
1426  SplineModeTree->Branch("j", &indices_mode[1], "j/I");
1427  SplineModeTree->Branch("k", &indices_mode[2], "k/I");
1428  SplineModeTree->Branch("value", &mode_value, "value/I");
1429 
1430  // Fill the tree
1431  for (size_t i = 0; i < SplineModeVecs.size(); ++i) {
1432  for (size_t j = 0; j < SplineModeVecs[i].size(); ++j) {
1433  for (size_t k = 0; k < SplineModeVecs[i][j].size(); ++k) {
1434  indices_mode[0] = static_cast<int>(i);
1435  indices_mode[1] = static_cast<int>(j);
1436  indices_mode[2] = static_cast<int>(k);
1437  mode_value = SplineModeVecs[i][j][k];
1438  SplineModeTree->Fill();
1439  }
1440  }
1441  }
1442  // Write the tree to the file
1443  SplineFile->cd();
1444  SplineModeTree->Write();
1445  delete SplineModeTree;
1446 }

◆ PrepareSettingsDir()

void BinnedSplineHandler::PrepareSettingsDir ( std::unique_ptr< TFile > &  SplineFile) const
private

KS: Prepare Settings Info within SplineFile.

Definition at line 1270 of file BinnedSplineHandler.cpp.

1270  {
1271 // *****************************************
1272  TTree *Settings = new TTree("Settings", "Settings");
1273  int CoeffIndex_temp = CoeffIndex;
1274  int MonolithSize_temp = MonolithSize;
1275  short int nParams_temp = nParams;
1276 
1277  Settings->Branch("CoeffIndex", &CoeffIndex_temp, "CoeffIndex/I");
1278  Settings->Branch("MonolithSize", &MonolithSize_temp, "MonolithSize/I");
1279  Settings->Branch("nParams", &nParams_temp, "nParams/S");
1280 
1281  int indexvec_sizes[7];
1282  indexvec_sizes[0] = static_cast<int>(indexvec.size());
1283  indexvec_sizes[1] = (indexvec_sizes[0] > 0) ? static_cast<int>(indexvec[0].size()) : 0;
1284  indexvec_sizes[2] = (indexvec_sizes[1] > 0) ? static_cast<int>(indexvec[0][0].size()) : 0;
1285  indexvec_sizes[3] = (indexvec_sizes[2] > 0) ? static_cast<int>(indexvec[0][0][0].size()) : 0;
1286  indexvec_sizes[4] = (indexvec_sizes[3] > 0) ? static_cast<int>(indexvec[0][0][0][0].size()) : 0;
1287  indexvec_sizes[5] = (indexvec_sizes[4] > 0) ? static_cast<int>(indexvec[0][0][0][0][0].size()) : 0;
1288  indexvec_sizes[6] = (indexvec_sizes[5] > 0) ? static_cast<int>(indexvec[0][0][0][0][0][0].size()) : 0;
1289 
1290  for (int i = 0; i < 7; ++i) {
1291  Settings->Branch(("indexvec_size" + std::to_string(i+1)).c_str(),
1292  &indexvec_sizes[i], ("indexvec_size" + std::to_string(i+1) + "/I").c_str());
1293  }
1294 
1295  int SplineBinning_size1 = static_cast<int>(SplineBinning.size());
1296  int SplineBinning_size2 = (SplineBinning_size1 > 0) ? static_cast<int>(SplineBinning[0].size()) : 0;
1297  int SplineBinning_size3 = (SplineBinning_size2 > 0) ? static_cast<int>(SplineBinning[0][0].size()) : 0;
1298 
1299  Settings->Branch("SplineBinning_size1", &SplineBinning_size1, "SplineBinning_size1/I");
1300  Settings->Branch("SplineBinning_size2", &SplineBinning_size2, "SplineBinning_size2/I");
1301  Settings->Branch("SplineBinning_size3", &SplineBinning_size3, "SplineBinning_size3/I");
1302 
1303  int SplineModeVecs_size1 = static_cast<int>(SplineModeVecs.size());
1304  int SplineModeVecs_size2 = (SplineModeVecs_size1 > 0) ? static_cast<int>(SplineModeVecs[0].size()) : 0;
1305  int SplineModeVecs_size3 = (SplineModeVecs_size2 > 0) ? static_cast<int>(SplineModeVecs[0][0].size()) : 0;
1306 
1307  Settings->Branch("SplineModeVecs_size1", &SplineModeVecs_size1, "SplineModeVecs_size1/I");
1308  Settings->Branch("SplineModeVecs_size2", &SplineModeVecs_size2, "SplineModeVecs_size2/I");
1309  Settings->Branch("SplineModeVecs_size3", &SplineModeVecs_size3, "SplineModeVecs_size3/I");
1310 
1311  std::vector<std::string> SampleNames_temp = SampleNames;
1312  Settings->Branch("SampleNames", &SampleNames_temp);
1313  std::vector<std::string> SampleTitles_temp = SampleTitles;
1314  Settings->Branch("SampleTitles", &SampleTitles_temp);
1315 
1316  Settings->Fill();
1317  SplineFile->cd();
1318  Settings->Write();
1319  delete Settings;
1320 }

◆ PrepareSplineFile()

void BinnedSplineHandler::PrepareSplineFile ( std::string  FileName)
overridevirtual

KS: Prepare spline file that can be used for fast loading.

Implements SplineBase.

Definition at line 1247 of file BinnedSplineHandler.cpp.

1247  {
1248 // *****************************************
1249  M3::AddPath(FileName);
1250  // Check for spaces in the filename
1251  size_t pos = FileName.find(' ');
1252  if (pos != std::string::npos) {
1253  MACH3LOG_WARN("Filename ({}) contains spaces. Replacing spaces with underscores.", FileName);
1254  while ((pos = FileName.find(' ')) != std::string::npos) {
1255  FileName[pos] = '_';
1256  }
1257  }
1258 
1259  auto SplineFile = std::make_unique<TFile>(FileName.c_str(), "recreate");
1260  PrepareSettingsDir(SplineFile);
1261  PrepareMonolithDir(SplineFile);
1262  PrepareIndexDir(SplineFile);
1263  PrepareOtherInfoDir(SplineFile);
1264  PrepareFastSplineInfoDir(SplineFile);
1265 
1266  SplineFile->Close();
1267 }
void PrepareIndexDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Index Info within SplineFile.
void PrepareMonolithDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Monolith Info within SplineFile.
void PrepareOtherInfoDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Other Info within SplineFile.
void PrepareSettingsDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Settings Info within SplineFile.
void PrepareFastSplineInfoDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Fast Spline Info within SplineFile.
Definition: SplineBase.cpp:131

◆ PrepForReweight()

void BinnedSplineHandler::PrepForReweight ( )

Definition at line 545 of file BinnedSplineHandler.cpp.

545  {
546 //****************************************
547  std::vector<TSpline3_red*> UniqueSystSplines;
548 
549  // DB Find all the Unique systs across each sample and oscillation channel
550  // This assumes that each occurence of the same systematic spline has the same knot spacing
551  // Which is a reasonable assumption for the current implementation of spline evaluations
552  for (unsigned int iSample = 0; iSample < indexvec.size(); iSample++)
553  { // Loop over systematics
554  for (unsigned int iSyst = 0; iSyst < indexvec[iSample][0].size(); iSyst++)
555  { // Loop over systematics
556  std::string SystName = SplineFileParPrefixNames[iSample][iSyst];
557  bool FoundSyst = false;
558 
559  //ETA - this always seems to be empty to begin with??
560  //so this loop never gets used?
561  for (unsigned int iFoundSyst = 0; iFoundSyst < UniqueSystNames.size(); iFoundSyst++)
562  {
563  if (SystName == UniqueSystNames[iFoundSyst])
564  {
565  FoundSyst = true;
566  }
567  }
568  if (!FoundSyst)
569  {
570  bool FoundNonFlatSpline = false;
571  //This is all basically to check that you have some resposne from a spline
572  //systematic somewhere
573  for (unsigned int iOscChan = 0; iOscChan < indexvec[iSample].size(); iOscChan++)
574  { // Loop over oscillation channels
575  for (unsigned int iMode = 0; iMode < indexvec[iSample][iOscChan][iSyst].size(); iMode++)
576  { // Loop over modes
577  for (unsigned int iVar1 = 0; iVar1 < indexvec[iSample][iOscChan][iSyst][iMode].size(); iVar1++)
578  { // Loop over first dimension
579  for (unsigned int iVar2 = 0; iVar2 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1].size(); iVar2++)
580  { // Loop over second dimension
581  for (unsigned int iVar3 = 0; iVar3 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2].size(); iVar3++)
582  { // Loop over third dimension
583  int splineindex=indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2][iVar3];
584  if (splinevec_Monolith[splineindex])
585  {
586  UniqueSystSplines.push_back(splinevec_Monolith[splineindex]);
587  UniqueSystIndices.push_back(GlobalSystIndex[iSample][iSyst]);
588  FoundNonFlatSpline = true;
589  }
590  if (FoundNonFlatSpline) { break;}
591  }//3D loop end
592  if (FoundNonFlatSpline) { break;}
593  }//2D loop end
594  if (FoundNonFlatSpline){ break; }
595  }//1D loop end
596  if (FoundNonFlatSpline){ break; }
597  }//mode loop end
598  if (FoundNonFlatSpline) { break; }
599  }//osc loop end
600  //ETA - only push back unique name if a non-flat response has been found
601  if(FoundNonFlatSpline){
602  UniqueSystNames.push_back(SystName);
603  }
604 
605  if (!FoundNonFlatSpline)
606  {
607  MACH3LOG_INFO("{} syst has no response in sample {}", SystName, iSample);
608  MACH3LOG_INFO("Whilst this isn't necessarily a problem, it seems odd");
609  continue;
610  }
611  }
612  }//Syst loop end
613  }
614 
615  nParams = static_cast<short int>(UniqueSystSplines.size());
616 
617  // DB Find the number of splines knots which assumes each instance of the syst has the same number of knots
618  SplineSegments = new short int[nParams]();
619  ParamValues = new float[nParams]();
620  SplineInfoArray.resize(nParams);
621  for (int iSpline = 0; iSpline < nParams; iSpline++)
622  {
623  SplineInfoArray[iSpline].nPts = static_cast<M3::int_t>(UniqueSystSplines[iSpline]->GetNp());
624  SplineInfoArray[iSpline].xPts.resize(SplineInfoArray[iSpline].nPts);
625  SplineInfoArray[iSpline].splineParsPointer = xsec->RetPointer(UniqueSystIndices[iSpline]);
626  for (int iKnot = 0; iKnot < SplineInfoArray[iSpline].nPts; iKnot++)
627  {
628  M3::float_t xPoint;
629  M3::float_t yPoint;
630  UniqueSystSplines[iSpline]->GetKnot(iKnot, xPoint, yPoint);
631  SplineInfoArray[iSpline].xPts[iKnot] = xPoint;
632  }
633  //ETA - let this just be set as the first segment by default
634  SplineSegments[iSpline] = 0;
635  ParamValues[iSpline] = 0.;
636  }
637 
638  MACH3LOG_INFO("nUniqueSysts: {}", nParams);
639  MACH3LOG_INFO("{:<15} | {:<20} | {:<6}", "Spline Index", "Syst Name", "nKnots");
640  for (int iUniqueSyst = 0; iUniqueSyst < nParams; iUniqueSyst++)
641  {
642  MACH3LOG_INFO("{:<15} | {:<20} | {:<6}", iUniqueSyst, UniqueSystNames[iUniqueSyst], SplineInfoArray[iUniqueSyst].nPts);
643  }
644 
645  //ETA
646  //Isn't this just doing what CountNumberOfLoadedSplines() does?
647  int nCombinations_FlatSplines = 0;
648  int nCombinations_All = 0;
649  // DB Now actually loop over splines to determine which are all null i.e. flat
650  for (unsigned int iSample = 0; iSample < indexvec.size(); iSample++)
651  { // Loop over systematics
652  for (unsigned int iOscChan = 0; iOscChan < indexvec[iSample].size(); iOscChan++)
653  { // Loop over oscillation channels
654  for (unsigned int iSyst = 0; iSyst < indexvec[iSample][iOscChan].size(); iSyst++)
655  { // Loop over systematics
656  for (unsigned int iMode = 0; iMode < indexvec[iSample][iOscChan][iSyst].size(); iMode++)
657  { // Loop over modes
658  //bool isFlat = false;
659  for (unsigned int iVar1 = 0; iVar1 < indexvec[iSample][iOscChan][iSyst][iMode].size(); iVar1++)
660  { // Loop over first dimension
661  for (unsigned int iVar2 = 0; iVar2 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1].size(); iVar2++)
662  { // Loop over second dimension
663  for (unsigned int iVar3 = 0; iVar3 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2].size(); iVar3++)
664  { // Loop over third dimension
665  int splineindex=indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2][iVar3];
666  if (splinevec_Monolith[splineindex])
667  {
668  nCombinations_All += 1;
669  } else{
670  nCombinations_FlatSplines += 1;
671  nCombinations_All += 1;
672  }
673  }
674  }
675  }
676  }
677  }
678  }
679  }
680  // We need to grab the maximum number of knots
681  MACH3LOG_INFO("Number of combinations of Sample, OscChan, Syst and Mode which have entirely flat response: {} / {}", nCombinations_FlatSplines, nCombinations_All);
682 }
int int_t
Definition: Core.h:31

◆ PrintArrayDetails()

void BinnedSplineHandler::PrintArrayDetails ( const std::string &  SampleTitle) const

Definition at line 759 of file BinnedSplineHandler.cpp.

761 {
762  int iSample = GetSampleIndex(SampleTitle);
763  int nOscChannels = int(indexvec[iSample].size());
764  MACH3LOG_INFO("Sample {} has {} oscillation channels", SampleTitle, nOscChannels);
765 
766  for (int iOscChan = 0; iOscChan < nOscChannels; iOscChan++)
767  {
768  int nSysts = int(indexvec[iSample][iOscChan].size());
769  MACH3LOG_INFO("Oscillation channel {} has {} systematics", iOscChan, nSysts);
770  }
771 
772  MACH3LOG_INFO("#----------------------------------------------------------------------------------------------------------------------------------#");
773  MACH3LOG_INFO("Printing no. of modes affected by each systematic for each oscillation channel");
774  for (unsigned int iOscChan = 0; iOscChan < indexvec[iSample].size(); iOscChan++) {
775  std::string modes = fmt::format("OscChan: {}\t", iOscChan);
776  for (unsigned int iSyst = 0; iSyst < indexvec[iSample][iOscChan].size(); iSyst++) {
777  modes += fmt::format("{} ", indexvec[iSample][iOscChan][iSyst].size());
778  }
779  MACH3LOG_INFO("{}", modes);
780  }
781  MACH3LOG_INFO("#----------------------------------------------------------------------------------------------------------------------------------#");
782 }

◆ PrintBinning()

void BinnedSplineHandler::PrintBinning ( TAxis *  Axis) const

Definition at line 824 of file BinnedSplineHandler.cpp.

826 {
827  const int NBins = Axis->GetNbins();
828  std::string text = "";
829  for (int iBin = 0; iBin <= NBins; iBin++) {
830  text += fmt::format("{} ", Axis->GetXbins()->GetAt(iBin));
831  }
832  MACH3LOG_INFO("{}", text);
833 }

◆ PrintSampleDetails()

void BinnedSplineHandler::PrintSampleDetails ( const std::string &  SampleTitle) const

Print info like Sample ID of spline params etc.

Definition at line 747 of file BinnedSplineHandler.cpp.

749 {
750  const int iSample = GetSampleIndex(SampleTitle);
751 
752  MACH3LOG_INFO("Details about sample: {:<20}", SampleTitles[iSample]);
753  MACH3LOG_INFO("\t Dimension: {:<35}", Dimensions[iSample]);
754  MACH3LOG_INFO("\t nSplineParam: {:<35}", nSplineParams[iSample]);
755  MACH3LOG_INFO("\t nOscChan: {:<35}", nOscChans[iSample]);
756 }

◆ retPointer()

const M3::float_t* BinnedSplineHandler::retPointer ( const int  sample,
const int  oscchan,
const int  syst,
const int  mode,
const int  var1bin,
const int  var2bin,
const int  var3bin 
) const
inline

get pointer to spline weight based on bin variables

Definition at line 65 of file BinnedSplineHandler.h.

66  {
67  int index = indexvec[sample][oscchan][syst][mode][var1bin][var2bin][var3bin];
68  return &weightvec_Monolith[index];
69  }

◆ StripDuplicatedModes()

std::vector< std::vector< int > > BinnedSplineHandler::StripDuplicatedModes ( const std::vector< std::vector< int > > &  InputVector)

Check if there are any repeated modes. This is used to reduce the number of modes in case many interaction modes get averaged into one spline.

Definition at line 892 of file BinnedSplineHandler.cpp.

892  {
893 
894  //ETA - this is of size nPars from the xsec model
895  size_t InputVectorSize = InputVector.size();
896  std::vector< std::vector<int> > ReturnVec(InputVectorSize);
897 
898  //ETA - loop over all systematics
899  for (size_t iSyst=0;iSyst<InputVectorSize;iSyst++) {
900  std::vector<int> TmpVec;
901  std::vector<std::string> TestVec;
902 
903  //Loop over the modes that we've listed in xsec cov
904  for (unsigned int iMode = 0 ; iMode < InputVector[iSyst].size() ; iMode++) {
905  int Mode = InputVector[iSyst][iMode];
906  std::string ModeName = Modes->GetSplineSuffixFromMaCh3Mode(Mode);
907 
908  bool IncludeMode = true;
909  for (auto TestString : TestVec) {
910  if (ModeName == TestString) {
911  IncludeMode = false;
912  break;
913  }
914  }
915 
916  if (IncludeMode) {
917  TmpVec.push_back(Mode);
918  TestVec.push_back(ModeName);
919  }
920  }
921 
922  ReturnVec[iSyst] = TmpVec;
923  }
924  return ReturnVec;
925 }

◆ TransferToMonolith()

void BinnedSplineHandler::TransferToMonolith ( )

flatten multidimensional spline array into proper monolith

Definition at line 198 of file BinnedSplineHandler.cpp.

200 {
201  PrepForReweight();
203 
206  MACH3LOG_ERROR("Something's gone wrong when we tried to get the size of your monolith");
207  MACH3LOG_ERROR("MonolithSize is {}", MonolithSize);
208  MACH3LOG_ERROR("MonolithIndex is {}", MonolithIndex);
209  throw MaCh3Exception(__FILE__ , __LINE__ );
210  }
211 
212  MACH3LOG_INFO("Now transferring splines to a monolith if size {}", MonolithSize);
213 
216  isflatarray = new bool[MonolithSize];
217 
220 
221  for (unsigned int iSample = 0; iSample < indexvec.size(); iSample++)
222  { // Loop over sample
223  for (unsigned int iOscChan = 0; iOscChan < indexvec[iSample].size(); iOscChan++)
224  { // Loop over oscillation channels
225  for (unsigned int iSyst = 0; iSyst < indexvec[iSample][iOscChan].size(); iSyst++)
226  { // Loop over systematics
227  for (unsigned int iMode = 0; iMode < indexvec[iSample][iOscChan][iSyst].size(); iMode++)
228  { // Loop over modes
229  for (unsigned int iVar1 = 0; iVar1 < indexvec[iSample][iOscChan][iSyst][iMode].size(); iVar1++)
230  { // Loop over first dimension
231  for (unsigned int iVar2 = 0; iVar2 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1].size(); iVar2++)
232  { // Loop over second dimension
233  for (unsigned int iVar3 = 0; iVar3 < indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2].size(); iVar3++)
234  {
235  int splineindex=indexvec[iSample][iOscChan][iSyst][iMode][iVar1][iVar2][iVar3];
236  weightvec_Monolith[splineindex] = 1.0;
237 
238  bool foundUniqueSpline = false;
239  for (int iUniqueSyst = 0; iUniqueSyst < nParams; iUniqueSyst++)
240  {
241  if (SplineFileParPrefixNames[iSample][iSyst] == UniqueSystNames[iUniqueSyst])
242  {
243  uniquesplinevec_Monolith[splineindex] = iUniqueSyst;
244  foundUniqueSpline = true;
245  }
246  }//unique syst loop end
247 
248  if (!foundUniqueSpline)
249  {
250  MACH3LOG_ERROR("Unique spline index not found");
251  MACH3LOG_ERROR("For Spline {}", SplineFileParPrefixNames[iSample][iSyst]);
252  MACH3LOG_ERROR("Couldn't match {} with any of the following {} systs:", SplineFileParPrefixNames[iSample][iSyst], nParams);
253  for (int iUniqueSyst = 0; iUniqueSyst < nParams; iUniqueSyst++)
254  {
255  MACH3LOG_ERROR("{},", UniqueSystNames.at(iUniqueSyst));
256  }//unique syst loop end
257  throw MaCh3Exception(__FILE__ , __LINE__ );
258  }
259 
260  int splineKnots;
261  if(splinevec_Monolith[splineindex]){
262  isflatarray[splineindex]=false;
263  splineKnots=splinevec_Monolith[splineindex]->GetNp();
264 
265  //Now to fill up our coefficient arrayss
266  M3::float_t* tmpXCoeffArr = new M3::float_t[splineKnots];
267  M3::float_t* tmpManyCoeffArr = new M3::float_t[splineKnots*_nCoeff_];
268 
269  int iCoeff=coeffindexvec[splineindex];
270  getSplineCoeff_SepMany(splineindex, tmpXCoeffArr, tmpManyCoeffArr);
271 
272  for(int i = 0; i < splineKnots; i++){
273  xcoeff_arr[iCoeff+i]=tmpXCoeffArr[i];
274 
275  for(int j = 0; j < _nCoeff_; j++){
276  manycoeff_arr[(iCoeff+i)*_nCoeff_+j]=tmpManyCoeffArr[i*_nCoeff_+j];
277  }
278  }
279  delete[] tmpXCoeffArr;
280  delete[] tmpManyCoeffArr;
281  }
282  else {
283  isflatarray[splineindex]=true;
284  }
285  }//3d loop end
286  }//2d loop end
287  }//1d loop end
288  }//mode loop end
289  }//syst2 loop end
290  }//osc loop end
291  }//syst1 loop end
292 }
void getSplineCoeff_SepMany(int splineindex, M3::float_t *&xArray, M3::float_t *&manyArray)
int CountNumberOfLoadedSplines(bool NonFlat=false, int Verbosity=0)
void InvestigateMissingSplines() const
This function will find missing splines in file.

Member Data Documentation

◆ CoeffIndex

int BinnedSplineHandler::CoeffIndex
protected

Definition at line 122 of file BinnedSplineHandler.h.

◆ coeffindexvec

std::vector<int > BinnedSplineHandler::coeffindexvec
protected

Definition at line 113 of file BinnedSplineHandler.h.

◆ DimensionLabels

std::vector<std::vector<std::string> > BinnedSplineHandler::DimensionLabels
protected

Definition at line 88 of file BinnedSplineHandler.h.

◆ Dimensions

std::vector<int> BinnedSplineHandler::Dimensions
protected

Definition at line 87 of file BinnedSplineHandler.h.

◆ GlobalSystIndex

std::vector< std::vector<int> > BinnedSplineHandler::GlobalSystIndex
protected

This holds the global spline index and is used to grab the current parameter value to evaluate splines at. Each internal vector will be of size of the number of spline systematics which affect that sample.

Definition at line 101 of file BinnedSplineHandler.h.

◆ indexvec

std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< int > > > > > > > BinnedSplineHandler::indexvec
protected

Variables related to determined which modes have splines and which piggy-back of other modes.

Definition at line 112 of file BinnedSplineHandler.h.

◆ isflatarray

bool* BinnedSplineHandler::isflatarray
protected

Need to keep track of which splines are flat and which aren't.

Definition at line 125 of file BinnedSplineHandler.h.

◆ manycoeff_arr

M3::float_t* BinnedSplineHandler::manycoeff_arr
protected

ybcd coefficients for each spline

Definition at line 129 of file BinnedSplineHandler.h.

◆ Modes

MaCh3Modes* BinnedSplineHandler::Modes
protected

pointer to MaCh3 Mode from which we get spline suffix

Definition at line 137 of file BinnedSplineHandler.h.

◆ MonolithIndex

int BinnedSplineHandler::MonolithIndex
protected

Definition at line 121 of file BinnedSplineHandler.h.

◆ MonolithSize

int BinnedSplineHandler::MonolithSize
protected

Definition at line 120 of file BinnedSplineHandler.h.

◆ nOscChans

std::vector<int> BinnedSplineHandler::nOscChans
protected

Definition at line 90 of file BinnedSplineHandler.h.

◆ nSplineParams

std::vector<int> BinnedSplineHandler::nSplineParams
protected

Definition at line 89 of file BinnedSplineHandler.h.

◆ SampleNames

std::vector<std::string> BinnedSplineHandler::SampleNames
protected

Definition at line 85 of file BinnedSplineHandler.h.

◆ SampleTitles

std::vector<std::string> BinnedSplineHandler::SampleTitles
protected

Definition at line 86 of file BinnedSplineHandler.h.

◆ SplineBinning

std::vector< std::vector< std::vector<TAxis*> > > BinnedSplineHandler::SplineBinning
protected

Definition at line 92 of file BinnedSplineHandler.h.

◆ SplineFileParPrefixNames

std::vector< std::vector<std::string> > BinnedSplineHandler::SplineFileParPrefixNames
protected

Definition at line 93 of file BinnedSplineHandler.h.

◆ SplineInterpolationTypes

std::vector< std::vector<SplineInterpolation> > BinnedSplineHandler::SplineInterpolationTypes
protected

spline interpolation types for each sample. These vectors are from a call to GetSplineInterpolationFromSampleID()

Definition at line 104 of file BinnedSplineHandler.h.

◆ SplineModeVecs

std::vector< std::vector< std::vector<int> > > BinnedSplineHandler::SplineModeVecs
protected

A vector of vectors of the spline modes that a systematic applies to This gets compared against the event mode to figure out if a syst should apply to an event or not

Definition at line 97 of file BinnedSplineHandler.h.

◆ splinevec_Monolith

std::vector< TSpline3_red* > BinnedSplineHandler::splinevec_Monolith
protected

holds each spline object before stripping into coefficient monolith

Definition at line 118 of file BinnedSplineHandler.h.

◆ uniquecoeffindices

std::vector<int> BinnedSplineHandler::uniquecoeffindices
protected

Unique coefficient indices.

Definition at line 115 of file BinnedSplineHandler.h.

◆ uniquesplinevec_Monolith

std::vector<int> BinnedSplineHandler::uniquesplinevec_Monolith
protected

Maps single spline object with single parameter.

Definition at line 134 of file BinnedSplineHandler.h.

◆ UniqueSystIndices

std::vector<int> BinnedSplineHandler::UniqueSystIndices
protected

Global index of each spline param, it allows us to match spline ordering with global.

Definition at line 109 of file BinnedSplineHandler.h.

◆ UniqueSystNames

std::vector<std::string> BinnedSplineHandler::UniqueSystNames
protected

name of each spline parameter

Definition at line 107 of file BinnedSplineHandler.h.

◆ weightvec_Monolith

std::vector<M3::float_t> BinnedSplineHandler::weightvec_Monolith
protected

Stores weight from spline evaluation for each single spline.

Definition at line 132 of file BinnedSplineHandler.h.

◆ xcoeff_arr

M3::float_t* BinnedSplineHandler::xcoeff_arr
protected

x coefficients for each spline

Definition at line 127 of file BinnedSplineHandler.h.

◆ xsec

ParameterHandlerGeneric* BinnedSplineHandler::xsec
protected

Pointer to covariance from which we get information about spline params.

Definition at line 82 of file BinnedSplineHandler.h.


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