MaCh3  2.4.2
Reference Guide
BinnedSplineHandler.cpp
Go to the documentation of this file.
1 #include "BinnedSplineHandler.h"
2 #include <memory>
3 
4 #pragma GCC diagnostic ignored "-Wuseless-cast"
5 #pragma GCC diagnostic ignored "-Wfloat-conversion"
6 
8 #include "TROOT.h"
9 #include "TKey.h"
10 #include "TH3F.h"
12 
13 //****************************************
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 }
33 
34 //****************************************
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 }
42 //****************************************
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 }
65 
66 //****************************************
67 void BinnedSplineHandler::AddSample(const std::string& SampleName,
68  const std::string& SampleTitle,
69  const std::vector<std::string>& OscChanFileNames,
70  const std::vector<std::string>& SplineVarNames)
71 //Adds samples to the large array
72 //****************************************
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 }
121 
122 //****************************************
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 }
196 
197 //****************************************
199 //****************************************
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 }
293 
294 // *****************************************
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 }
307 
308 //****************************************
310 //****************************************
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 }
346 
347 //****************************************
348 //Creates an array to be filled with monolith indexes for each sample (allows for indexing between 7D binning and 1D Vector)
349 //Only need 1 indexing array everything else interfaces with this to get binning properties
350 void BinnedSplineHandler::BuildSampleIndexingArray(const std::string& SampleTitle)
351 //****************************************
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 }
380 
381 //****************************************
382 std::vector<TAxis *> BinnedSplineHandler::FindSplineBinning(const std::string& FileName, const std::string& SampleTitle)
383 //****************************************
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 }
483 
484 //****************************************
485 int BinnedSplineHandler::CountNumberOfLoadedSplines(bool NonFlat, int Verbosity)
486 //****************************************
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 }
543 
544 //****************************************
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 }
683 
684 //****************************************
685 // Rather work with spline coefficients in the splines, let's copy ND and use coefficient arrays
686 void BinnedSplineHandler::getSplineCoeff_SepMany(int splineindex, M3::float_t* &xArray, M3::float_t* &manyArray){
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 }
719 
720 //****************************************
721 //Equally though could just use KinematicVariable to map back
722 std::string BinnedSplineHandler::getDimLabel(const int iSample, const unsigned int Axis) const
723 //****************************************
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 }
732 
733 //****************************************
734 //Returns sample index in
735 int BinnedSplineHandler::GetSampleIndex(const std::string& SampleTitle) const{
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 }
745 
746 //****************************************
747 void BinnedSplineHandler::PrintSampleDetails(const std::string& SampleTitle) const
748 //****************************************
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 }
757 
758 //****************************************
759 void BinnedSplineHandler::PrintArrayDetails(const std::string& SampleTitle) const
760 //****************************************
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 }
783 
784 //****************************************
785 bool BinnedSplineHandler::isValidSplineIndex(const std::string& SampleTitle, int iOscChan, int iSyst, int iMode, int iVar1, int iVar2, int iVar3)
786 //****************************************
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 }
822 
823 //****************************************
824 void BinnedSplineHandler::PrintBinning(TAxis *Axis) const
825 //****************************************
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 }
834 
835 //****************************************
836 std::vector< std::vector<int> > BinnedSplineHandler::GetEventSplines(const std::string& SampleTitle, int iOscChan, int EventMode, double Var1Val, double Var2Val, double Var3Val)
837 //****************************************
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 }
889 
890 // checks if there are multiple modes with the same SplineSuffix
891 // (for example if CCRES and CCCoherent are treated as one spline mode)
892 std::vector< std::vector<int> > BinnedSplineHandler::StripDuplicatedModes(const std::vector< std::vector<int> >& InputVector) {
893  //ETA - this is of size nPars from the xsec model
894  size_t InputVectorSize = InputVector.size();
895  std::vector< std::vector<int> > ReturnVec(InputVectorSize);
896 
897  //ETA - loop over all systematics
898  for (size_t iSyst=0;iSyst<InputVectorSize;iSyst++) {
899  std::vector<int> TmpVec;
900  std::vector<std::string> TestVec;
901 
902  //Loop over the modes that we've listed in xsec cov
903  for (unsigned int iMode = 0 ; iMode < InputVector[iSyst].size() ; iMode++) {
904  int Mode = InputVector[iSyst][iMode];
905  std::string ModeName = Modes->GetSplineSuffixFromMaCh3Mode(Mode);
906 
907  bool IncludeMode = true;
908  for (auto TestString : TestVec) {
909  if (ModeName == TestString) {
910  IncludeMode = false;
911  break;
912  }
913  }
914 
915  if (IncludeMode) {
916  TmpVec.push_back(Mode);
917  TestVec.push_back(ModeName);
918  }
919  }
920 
921  ReturnVec[iSyst] = TmpVec;
922  }
923  return ReturnVec;
924 }
925 
926 void BinnedSplineHandler::FillSampleArray(std::string SampleTitle, std::vector<std::string> OscChanFileNames)
927 {
928  int iSample = GetSampleIndex(SampleTitle);
929  int nOscChannels = nOscChans[iSample];
930 
931  for (int iOscChan = 0; iOscChan < nOscChannels; iOscChan++) {
932  MACH3LOG_INFO("Processing: {}", OscChanFileNames[iOscChan]);
933  TSpline3* mySpline = nullptr;
934  TSpline3_red* Spline = nullptr;
935  TString Syst, Mode;
936  int nKnots, SystNum, ModeNum, Var1Bin, Var2Bin, Var3Bin = M3::_BAD_INT_;
937  double x,y = M3::_BAD_DOUBLE_;
938  bool isFlat = true;
939 
940  std::set<std::string> SplineFileNames;
941 
942  auto File = std::unique_ptr<TFile>(TFile::Open(OscChanFileNames[iOscChan].c_str()));
943 
944  if (!File || File->IsZombie()) {
945  MACH3LOG_ERROR("File {} not found", OscChanFileNames[iOscChan]);
946  throw MaCh3Exception(__FILE__, __LINE__);
947  }
948 
949  //This is the MC specific part of the code
950  //i.e. we always assume that the splines are just store in single TDirectory and they're all in there as single objects
951  for (auto k : *File->GetListOfKeys()) {
952  auto Key = static_cast<TKey*>(k);
953  TClass *Class = gROOT->GetClass(Key->GetClassName(), false);
954  if(!Class->InheritsFrom("TSpline3")) {
955  continue;
956  }
957 
958  std::string FullSplineName = std::string(Key->GetName());
959 
960  if (SplineFileNames.count(FullSplineName) > 0) {
961  MACH3LOG_CRITICAL("Skipping spline - Found a spline whose name has already been encountered before: {}", FullSplineName);
962  continue;
963  }
964  SplineFileNames.insert(FullSplineName);
965 
966  std::vector<std::string> Tokens = GetTokensFromSplineName(FullSplineName);
967 
968  if (Tokens.size() != kNTokens) {
969  MACH3LOG_ERROR("Invalid tokens from spline name - Expected {} tokens. Check implementation in GetTokensFromSplineName()", static_cast<int>(kNTokens));
970  throw MaCh3Exception(__FILE__, __LINE__);
971  }
972 
973  Syst = Tokens[kSystToken];
974  Mode = Tokens[kModeToken];
975  Var1Bin = std::stoi(Tokens[kVar1BinToken]);
976  Var2Bin = std::stoi(Tokens[kVar2BinToken]);
977  Var3Bin = std::stoi(Tokens[kVar3BinToken]);
978 
979  SystNum = -1;
980  for (unsigned iSyst = 0; iSyst < SplineFileParPrefixNames[iSample].size(); iSyst++) {
981  if (Syst == SplineFileParPrefixNames[iSample][iSyst]) {
982  SystNum = iSyst;
983  break;
984  }
985  }
986 
987  // If the syst doesn't match any of the spline names then skip it
988  if (SystNum == -1){
989  MACH3LOG_DEBUG("Couldn't match!!");
990  MACH3LOG_DEBUG("Couldn't Match any systematic name in ParameterHandler with spline name: {}" , FullSplineName);
991  continue;
992  }
993 
994  ModeNum = -1;
995  for (unsigned int iMode = 0; iMode < SplineModeVecs[iSample][SystNum].size(); iMode++) {
996  if (Mode == Modes->GetSplineSuffixFromMaCh3Mode(SplineModeVecs[iSample][SystNum][iMode])) {
997  ModeNum = iMode;
998  break;
999  }
1000  }
1001 
1002  if (ModeNum == -1) {
1003  //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
1004  // Therefore include as debug warning and continue instead
1005  MACH3LOG_DEBUG("Couldn't find mode for {} in {}. Problem Spline is : {} ", Mode, Syst, FullSplineName);
1006  continue;
1007  }
1008 
1009  mySpline = Key->ReadObject<TSpline3>();
1010 
1011  if (isValidSplineIndex(SampleTitle, iOscChan, SystNum, ModeNum, Var1Bin, Var2Bin, Var3Bin)) { // loop over all the spline knots and check their value
1012  MACH3LOG_TRACE("Pushed back monolith for spline {}", FullSplineName);
1013  // if the value is 1 then set the flat bool to false
1014  nKnots = mySpline->GetNp();
1015  isFlat = true;
1016  for (int iKnot = 0; iKnot < nKnots; iKnot++) {
1017  mySpline->GetKnot(iKnot, x, y);
1018  if (y < 0.99999 || y > 1.00001)
1019  {
1020  isFlat = false;
1021  break;
1022  }
1023  }
1024 
1025  //Rather than keeping a mega vector of splines then converting, this should just keep everything nice in memory!
1026  indexvec[iSample][iOscChan][SystNum][ModeNum][Var1Bin][Var2Bin][Var3Bin]=MonolithIndex;
1027  coeffindexvec.push_back(CoeffIndex);
1028  // Should save memory rather saving [x_i_0 ,... x_i_maxknots] for every spline!
1029  if (isFlat) {
1030  splinevec_Monolith.push_back(nullptr);
1031  delete mySpline;
1032  } else {
1033  ApplyKnotWeightCapTSpline3(mySpline, SystNum, xsec);
1034  Spline = new TSpline3_red(mySpline, SplineInterpolationTypes[iSample][SystNum]);
1035  if(mySpline) delete mySpline;
1036 
1037  splinevec_Monolith.push_back(Spline);
1038  uniquecoeffindices.push_back(MonolithIndex); //So we can get the unique coefficients and skip flat splines later on!
1039  CoeffIndex+=nKnots;
1040  }
1041  //Incrementing MonolithIndex to keep track of number of valid spline indices
1042  MonolithIndex+=1;
1043  } else {
1044  //Potentially you are not a valid spline index
1045  delete mySpline;
1046  }
1047  }//End of loop over all TKeys in file
1048 
1049  //A bit of clean up
1050  File->Delete("*");
1051  File->Close();
1052  } //End of oscillation channel loop
1053 }
1054 
1055 // *****************************************
1056 // Load SplineMonolith from ROOT file
1057 void BinnedSplineHandler::LoadSplineFile(std::string FileName) {
1058 // *****************************************
1059  M3::AddPath(FileName);
1060 
1061  // Check for spaces in the filename
1062  size_t pos = FileName.find(' ');
1063  if (pos != std::string::npos) {
1064  MACH3LOG_WARN("Filename ({}) contains spaces. Replacing spaces with underscores.", FileName);
1065  while ((pos = FileName.find(' ')) != std::string::npos) {
1066  FileName[pos] = '_';
1067  }
1068  }
1069  auto SplineFile = std::make_unique<TFile>(FileName.c_str(), "OPEN");
1070 
1071  TMacro *ConfigCov = SplineFile->Get<TMacro>("ParameterHandler");
1072  // Config which was in MCMC from which we are starting
1073  YAML::Node CovSettings = TMacroToYAML(*ConfigCov);
1074  // Config from currently used cov object
1075  YAML::Node ConfigCurrent = xsec->GetConfig();
1076 
1077  if (!compareYAMLNodes(CovSettings, ConfigCurrent))
1078  {
1079  MACH3LOG_ERROR("Loading precomputed spline file, however encountered different YAML config, please regenerate input");
1080  throw MaCh3Exception(__FILE__ , __LINE__ );
1081  }
1082 
1083  LoadSettingsDir(SplineFile);
1084  LoadMonolithDir(SplineFile);
1085  LoadIndexDir(SplineFile);
1086  LoadFastSplineInfoDir(SplineFile);
1087 
1088  for (int iSpline = 0; iSpline < nParams; iSpline++) {
1089  SplineInfoArray[iSpline].splineParsPointer = xsec->RetPointer(UniqueSystIndices[iSpline]);
1090  }
1091  SplineFile->Close();
1092 }
1093 
1094 // *****************************************
1095 // KS: Prepare Fast Spline Info within SplineFile
1096 void BinnedSplineHandler::LoadSettingsDir(std::unique_ptr<TFile>& SplineFile) {
1097 // *****************************************
1098  TTree *Settings = SplineFile->Get<TTree>("Settings");
1099  int CoeffIndex_temp, MonolithSize_temp;
1100  short int nParams_temp;
1101  Settings->SetBranchAddress("CoeffIndex", &CoeffIndex_temp);
1102  Settings->SetBranchAddress("MonolithSize", &MonolithSize_temp);
1103  Settings->SetBranchAddress("nParams", &nParams_temp);
1104  int indexvec_sizes[7];
1105  for (int i = 0; i < 7; ++i) {
1106  Settings->SetBranchAddress(("indexvec_size" + std::to_string(i+1)).c_str(), &indexvec_sizes[i]);
1107  }
1108 
1109  int SplineBinning_size1, SplineBinning_size2, SplineBinning_size3;
1110  Settings->SetBranchAddress("SplineBinning_size1", &SplineBinning_size1);
1111  Settings->SetBranchAddress("SplineBinning_size2", &SplineBinning_size2);
1112  Settings->SetBranchAddress("SplineBinning_size3", &SplineBinning_size3);
1113  int SplineModeVecs_size1, SplineModeVecs_size2, SplineModeVecs_size3;
1114  Settings->SetBranchAddress("SplineModeVecs_size1", &SplineModeVecs_size1);
1115  Settings->SetBranchAddress("SplineModeVecs_size2", &SplineModeVecs_size2);
1116  Settings->SetBranchAddress("SplineModeVecs_size3", &SplineModeVecs_size3);
1117  std::vector<std::string>* SampleNames_temp = nullptr;
1118  Settings->SetBranchAddress("SampleNames", &SampleNames_temp);
1119  std::vector<std::string>* SampleTitles_temp = nullptr;
1120  Settings->SetBranchAddress("SampleTitles", &SampleTitles_temp);
1121  Settings->GetEntry(0);
1122 
1123  CoeffIndex = CoeffIndex_temp;
1124  MonolithSize = MonolithSize_temp;
1125  SampleNames = *SampleNames_temp;
1126  SampleTitles = *SampleTitles_temp;
1127 
1128  nParams = nParams_temp;
1129 
1130  SplineSegments = new short int[nParams]();
1131  ParamValues = new float[nParams]();
1132 
1133  // Resize indexvec according to saved dimensions
1134  indexvec.resize(indexvec_sizes[0]);
1135  for (int i = 0; i < indexvec_sizes[0]; ++i) {
1136  indexvec[i].resize(indexvec_sizes[1]);
1137  for (int j = 0; j < indexvec_sizes[1]; ++j) {
1138  indexvec[i][j].resize(indexvec_sizes[2]);
1139  for (int k = 0; k < indexvec_sizes[2]; ++k) {
1140  indexvec[i][j][k].resize(indexvec_sizes[3]);
1141  for (int l = 0; l < indexvec_sizes[3]; ++l) {
1142  indexvec[i][j][k][l].resize(indexvec_sizes[4]);
1143  for (int m = 0; m < indexvec_sizes[4]; ++m) {
1144  indexvec[i][j][k][l][m].resize(indexvec_sizes[5]);
1145  for (int n = 0; n < indexvec_sizes[5]; ++n) {
1146  indexvec[i][j][k][l][m][n].resize(indexvec_sizes[6]);
1147  }
1148  }
1149  }
1150  }
1151  }
1152  }
1153 
1154  auto Resize3D = [](auto& vec, int d1, int d2, int d3) {
1155  vec.resize(d1);
1156  for (int i = 0; i < d1; ++i) {
1157  vec[i].resize(d2);
1158  for (int j = 0; j < d2; ++j) {
1159  vec[i][j].resize(d3);
1160  }
1161  }
1162  };
1163 
1164  Resize3D(SplineBinning, SplineBinning_size1, SplineBinning_size2, SplineBinning_size3);
1165  Resize3D(SplineModeVecs, SplineModeVecs_size1, SplineModeVecs_size2, SplineModeVecs_size3);
1166 }
1167 
1168 // *****************************************
1169 // KS: Prepare Fast Spline Info within SplineFile
1170 void BinnedSplineHandler::LoadMonolithDir(std::unique_ptr<TFile>& SplineFile) {
1171 // *****************************************
1172  TTree *MonolithTree = SplineFile->Get<TTree>("MonolithTree");
1173 
1175  MonolithTree->SetBranchAddress("manycoeff", manycoeff_arr);
1176  isflatarray = new bool[MonolithSize];
1178  MonolithTree->SetBranchAddress("isflatarray", isflatarray);
1179 
1180  // Load vectors
1181  std::vector<int>* coeffindexvec_temp = nullptr;
1182  MonolithTree->SetBranchAddress("coeffindexvec", &coeffindexvec_temp);
1183  std::vector<int>* uniquecoeffindices_temp = nullptr;
1184  MonolithTree->SetBranchAddress("uniquecoeffindices", &uniquecoeffindices_temp);
1185  std::vector<int>* uniquesplinevec_Monolith_temp = nullptr;
1186  MonolithTree->SetBranchAddress("uniquesplinevec_Monolith", &uniquesplinevec_Monolith_temp);
1187  std::vector<int>* UniqueSystIndices_temp = nullptr;
1188  MonolithTree->SetBranchAddress("UniqueSystIndices", &UniqueSystIndices_temp);
1189 
1190  // Allocate and load xcoeff_arr
1192  MonolithTree->SetBranchAddress("xcoeff", xcoeff_arr);
1193 
1194  MonolithTree->GetEntry(0);
1195 
1196  coeffindexvec = *coeffindexvec_temp;
1197  uniquecoeffindices = *uniquecoeffindices_temp;
1198  uniquesplinevec_Monolith = *uniquesplinevec_Monolith_temp;
1199  UniqueSystIndices = *UniqueSystIndices_temp;
1200 }
1201 
1202 // *****************************************
1203 // KS: Prepare Fast Spline Info within SplineFile
1204 void BinnedSplineHandler::LoadIndexDir(std::unique_ptr<TFile>& SplineFile) {
1205 // *****************************************
1206  TTree *IndexTree = SplineFile->Get<TTree>("IndexVec");
1207 
1208  std::vector<int> Dim(7);
1209  int value;
1210  for (int d = 0; d < 7; ++d) {
1211  IndexTree->SetBranchAddress(Form("dim%d", d+1), &Dim[d]);
1212  }
1213  IndexTree->SetBranchAddress("value", &value);
1214 
1215  // Fill indexvec with data from IndexTree
1216  for (Long64_t i = 0; i < IndexTree->GetEntries(); ++i) {
1217  IndexTree->GetEntry(i);
1218  indexvec[Dim[0]][Dim[1]][Dim[2]][Dim[3]][Dim[4]][Dim[5]][Dim[6]] = value;
1219  }
1220 
1221  // Load SplineBinning data
1222  TTree *SplineBinningTree = SplineFile->Get<TTree>("SplineBinningTree");
1223  std::vector<int> indices(3);
1224  SplineBinningTree->SetBranchAddress("i", &indices[0]);
1225  SplineBinningTree->SetBranchAddress("j", &indices[1]);
1226  SplineBinningTree->SetBranchAddress("k", &indices[2]);
1227  TAxis* axis = nullptr;
1228  SplineBinningTree->SetBranchAddress("axis", &axis);
1229 
1230  // Reconstruct TAxis objects
1231  for (Long64_t entry = 0; entry < SplineBinningTree->GetEntries(); ++entry) {
1232  SplineBinningTree->GetEntry(entry);
1233  int i = indices[0];
1234  int j = indices[1];
1235  int k = indices[2];
1236  SplineBinning[i][j][k] = static_cast<TAxis*>(axis->Clone());
1237  }
1238 
1239  std::vector<int> indices_mode(3);
1240  int mode_value;
1241  TTree *SplineModeTree = SplineFile->Get<TTree>("SplineModeTree");
1242  SplineModeTree->SetBranchAddress("i", &indices_mode[0]);
1243  SplineModeTree->SetBranchAddress("j", &indices_mode[1]);
1244  SplineModeTree->SetBranchAddress("k", &indices_mode[2]);
1245  SplineModeTree->SetBranchAddress("value", &mode_value);
1246 
1247  // Fill SplineModeVecs with values from the tree
1248  for (Long64_t entry = 0; entry < SplineModeTree->GetEntries(); ++entry) {
1249  SplineModeTree->GetEntry(entry);
1250  int i = indices_mode[0];
1251  int j = indices_mode[1];
1252  int k = indices_mode[2];
1253  SplineModeVecs[i][j][k] = mode_value;
1254  }
1255 }
1256 
1257 // *****************************************
1258 // Save SplineMonolith into ROOT file
1259 void BinnedSplineHandler::PrepareSplineFile(std::string FileName) {
1260 // *****************************************
1261  M3::AddPath(FileName);
1262  // Check for spaces in the filename
1263  size_t pos = FileName.find(' ');
1264  if (pos != std::string::npos) {
1265  MACH3LOG_WARN("Filename ({}) contains spaces. Replacing spaces with underscores.", FileName);
1266  while ((pos = FileName.find(' ')) != std::string::npos) {
1267  FileName[pos] = '_';
1268  }
1269  }
1270  // Save ROOT File
1271  auto SplineFile = std::make_unique<TFile>(FileName.c_str(), "recreate");
1272  YAML::Node ConfigCurrent = xsec->GetConfig();
1273  TMacro ConfigSave = YAMLtoTMacro(ConfigCurrent, "ParameterHandler");
1274  ConfigSave.Write();
1275 
1276  PrepareSettingsDir(SplineFile);
1277  PrepareMonolithDir(SplineFile);
1278  PrepareIndexDir(SplineFile);
1279  PrepareOtherInfoDir(SplineFile);
1280  PrepareFastSplineInfoDir(SplineFile);
1281 
1282  SplineFile->Close();
1283 }
1284 
1285 // *****************************************
1286 void BinnedSplineHandler::PrepareSettingsDir(std::unique_ptr<TFile>& SplineFile) const {
1287 // *****************************************
1288  TTree *Settings = new TTree("Settings", "Settings");
1289  int CoeffIndex_temp = CoeffIndex;
1290  int MonolithSize_temp = MonolithSize;
1291  short int nParams_temp = nParams;
1292 
1293  Settings->Branch("CoeffIndex", &CoeffIndex_temp, "CoeffIndex/I");
1294  Settings->Branch("MonolithSize", &MonolithSize_temp, "MonolithSize/I");
1295  Settings->Branch("nParams", &nParams_temp, "nParams/S");
1296 
1297  int indexvec_sizes[7];
1298  indexvec_sizes[0] = static_cast<int>(indexvec.size());
1299  indexvec_sizes[1] = (indexvec_sizes[0] > 0) ? static_cast<int>(indexvec[0].size()) : 0;
1300  indexvec_sizes[2] = (indexvec_sizes[1] > 0) ? static_cast<int>(indexvec[0][0].size()) : 0;
1301  indexvec_sizes[3] = (indexvec_sizes[2] > 0) ? static_cast<int>(indexvec[0][0][0].size()) : 0;
1302  indexvec_sizes[4] = (indexvec_sizes[3] > 0) ? static_cast<int>(indexvec[0][0][0][0].size()) : 0;
1303  indexvec_sizes[5] = (indexvec_sizes[4] > 0) ? static_cast<int>(indexvec[0][0][0][0][0].size()) : 0;
1304  indexvec_sizes[6] = (indexvec_sizes[5] > 0) ? static_cast<int>(indexvec[0][0][0][0][0][0].size()) : 0;
1305 
1306  for (int i = 0; i < 7; ++i) {
1307  Settings->Branch(("indexvec_size" + std::to_string(i+1)).c_str(),
1308  &indexvec_sizes[i], ("indexvec_size" + std::to_string(i+1) + "/I").c_str());
1309  }
1310 
1311  int SplineBinning_size1 = static_cast<int>(SplineBinning.size());
1312  int SplineBinning_size2 = (SplineBinning_size1 > 0) ? static_cast<int>(SplineBinning[0].size()) : 0;
1313  int SplineBinning_size3 = (SplineBinning_size2 > 0) ? static_cast<int>(SplineBinning[0][0].size()) : 0;
1314 
1315  Settings->Branch("SplineBinning_size1", &SplineBinning_size1, "SplineBinning_size1/I");
1316  Settings->Branch("SplineBinning_size2", &SplineBinning_size2, "SplineBinning_size2/I");
1317  Settings->Branch("SplineBinning_size3", &SplineBinning_size3, "SplineBinning_size3/I");
1318 
1319  int SplineModeVecs_size1 = static_cast<int>(SplineModeVecs.size());
1320  int SplineModeVecs_size2 = (SplineModeVecs_size1 > 0) ? static_cast<int>(SplineModeVecs[0].size()) : 0;
1321  int SplineModeVecs_size3 = (SplineModeVecs_size2 > 0) ? static_cast<int>(SplineModeVecs[0][0].size()) : 0;
1322 
1323  Settings->Branch("SplineModeVecs_size1", &SplineModeVecs_size1, "SplineModeVecs_size1/I");
1324  Settings->Branch("SplineModeVecs_size2", &SplineModeVecs_size2, "SplineModeVecs_size2/I");
1325  Settings->Branch("SplineModeVecs_size3", &SplineModeVecs_size3, "SplineModeVecs_size3/I");
1326 
1327  std::vector<std::string> SampleNames_temp = SampleNames;
1328  Settings->Branch("SampleNames", &SampleNames_temp);
1329  std::vector<std::string> SampleTitles_temp = SampleTitles;
1330  Settings->Branch("SampleTitles", &SampleTitles_temp);
1331 
1332  Settings->Fill();
1333  SplineFile->cd();
1334  Settings->Write();
1335  delete Settings;
1336 }
1337 
1338 // *****************************************
1339 void BinnedSplineHandler::PrepareMonolithDir(std::unique_ptr<TFile>& SplineFile) const {
1340 // *****************************************
1341  TTree *MonolithTree = new TTree("MonolithTree", "MonolithTree");
1342  MonolithTree->Branch("manycoeff", manycoeff_arr, Form("manycoeff[%d]/%s", CoeffIndex * _nCoeff_, M3::float_t_str_repr));
1343  MonolithTree->Branch("isflatarray", isflatarray, Form("isflatarray[%d]/O", MonolithSize));
1344 
1345  std::vector<int> coeffindexvec_temp = coeffindexvec;
1346  MonolithTree->Branch("coeffindexvec", &coeffindexvec_temp);
1347  std::vector<int> uniquecoeffindices_temp = uniquecoeffindices;
1348  MonolithTree->Branch("uniquecoeffindices", &uniquecoeffindices_temp);
1349  std::vector<int> uniquesplinevec_Monolith_temp = uniquesplinevec_Monolith;
1350  MonolithTree->Branch("uniquesplinevec_Monolith", &uniquesplinevec_Monolith_temp);
1351  std::vector<int> UniqueSystIndices_temp = UniqueSystIndices;
1352  MonolithTree->Branch("UniqueSystIndices", &UniqueSystIndices_temp);
1353  MonolithTree->Branch("xcoeff", xcoeff_arr, Form("xcoeff[%d]/%s", CoeffIndex, M3::float_t_str_repr));
1354 
1355  MonolithTree->Fill();
1356  SplineFile->cd();
1357  MonolithTree->Write();
1358  delete MonolithTree;
1359 }
1360 
1361 // *****************************************
1362 void BinnedSplineHandler::PrepareIndexDir(std::unique_ptr<TFile>& SplineFile) const {
1363 // *****************************************
1364  // Create a TTree to store the data
1365  TTree *IndexTree = new TTree("IndexVec", "IndexVec");
1366 
1367  // Vector holding the 7 dims
1368  std::vector<int> Dim(7);
1369  int value;
1370 
1371  // Create branches for each dimension
1372  for (int d = 0; d < 7; ++d) {
1373  IndexTree->Branch(Form("dim%d", d+1), &Dim[d], Form("dim%d/I", d+1));
1374  }
1375  IndexTree->Branch("value", &value, "value/I");
1376 
1377  // Fill the tree
1378  for (size_t i = 0; i < indexvec.size(); ++i) {
1379  for (size_t j = 0; j < indexvec[i].size(); ++j) {
1380  for (size_t k = 0; k < indexvec[i][j].size(); ++k) {
1381  for (size_t l = 0; l < indexvec[i][j][k].size(); ++l) {
1382  for (size_t m = 0; m < indexvec[i][j][k][l].size(); ++m) {
1383  for (size_t n = 0; n < indexvec[i][j][k][l][m].size(); ++n) {
1384  for (size_t p = 0; p < indexvec[i][j][k][l][m][n].size(); ++p) {
1385  Dim[0] = static_cast<int>(i);
1386  Dim[1] = static_cast<int>(j);
1387  Dim[2] = static_cast<int>(k);
1388  Dim[3] = static_cast<int>(l);
1389  Dim[4] = static_cast<int>(m);
1390  Dim[5] = static_cast<int>(n);
1391  Dim[6] = static_cast<int>(p);
1392  value = static_cast<int>(indexvec[i][j][k][l][m][n][p]);
1393  IndexTree->Fill();
1394  }
1395  }
1396  }
1397  }
1398  }
1399  }
1400  }
1401 
1402  SplineFile->cd();
1403  // Write the tree to the file
1404  IndexTree->Write();
1405  delete IndexTree;
1406 }
1407 
1408 // *****************************************
1409 void BinnedSplineHandler::PrepareOtherInfoDir(std::unique_ptr<TFile>& SplineFile) const {
1410 // *****************************************
1411  // Create a new tree for SplineBinning data
1412  TTree *SplineBinningTree = new TTree("SplineBinningTree", "SplineBinningTree");
1413  std::vector<int> indices(3); // To store the 3D indices
1414  TAxis* axis = nullptr;
1415  SplineBinningTree->Branch("i", &indices[0], "i/I");
1416  SplineBinningTree->Branch("j", &indices[1], "j/I");
1417  SplineBinningTree->Branch("k", &indices[2], "k/I");
1418  SplineBinningTree->Branch("axis", "TAxis", &axis);
1419 
1420  // Fill the SplineBinningTree
1421  for (size_t i = 0; i < SplineBinning.size(); ++i) {
1422  for (size_t j = 0; j < SplineBinning[i].size(); ++j) {
1423  for (size_t k = 0; k < SplineBinning[i][j].size(); ++k) {
1424  axis = SplineBinning[i][j][k];
1425  indices[0] = static_cast<int>(i);
1426  indices[1] = static_cast<int>(j);
1427  indices[2] = static_cast<int>(k);
1428  SplineBinningTree->Fill();
1429  }
1430  }
1431  }
1432  SplineFile->cd();
1433  SplineBinningTree->Write();
1434  delete SplineBinningTree;
1435 
1436  std::vector<int> indices_mode(3); // to store 3D indices
1437  int mode_value;
1438 
1439  TTree *SplineModeTree = new TTree("SplineModeTree", "SplineModeTree");
1440  // Create branches for indices and value
1441  SplineModeTree->Branch("i", &indices_mode[0], "i/I");
1442  SplineModeTree->Branch("j", &indices_mode[1], "j/I");
1443  SplineModeTree->Branch("k", &indices_mode[2], "k/I");
1444  SplineModeTree->Branch("value", &mode_value, "value/I");
1445 
1446  // Fill the tree
1447  for (size_t i = 0; i < SplineModeVecs.size(); ++i) {
1448  for (size_t j = 0; j < SplineModeVecs[i].size(); ++j) {
1449  for (size_t k = 0; k < SplineModeVecs[i][j].size(); ++k) {
1450  indices_mode[0] = static_cast<int>(i);
1451  indices_mode[1] = static_cast<int>(j);
1452  indices_mode[2] = static_cast<int>(k);
1453  mode_value = SplineModeVecs[i][j][k];
1454  SplineModeTree->Fill();
1455  }
1456  }
1457  }
1458  // Write the tree to the file
1459  SplineFile->cd();
1460  SplineModeTree->Write();
1461  delete SplineModeTree;
1462 }
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:140
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:38
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
#define MACH3LOG_TRACE
Definition: MaCh3Logger.h:33
void CleanContainer(T &)
Base case: do nothing for non-pointer types.
void CleanVector(T &)
Base case: do nothing for non-vector types.
@ kSpline
For splined parameters (1D)
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
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...
Definition: SplineStructs.h:98
bool isFlat(TSpline3_red *&spl)
CW: Helper function used in the constructor, tests to see if the spline is flat.
TMacro YAMLtoTMacro(const YAML::Node &yaml_node, const std::string &name)
Convert a YAML node to a ROOT TMacro object.
Definition: YamlHelper.h:167
YAML::Node TMacroToYAML(const TMacro &macro)
KS: Convert a ROOT TMacro object to a YAML node.
Definition: YamlHelper.h:152
bool compareYAMLNodes(const YAML::Node &node1, const YAML::Node &node2, bool Mute=false)
Compare if yaml nodes are identical.
Definition: YamlHelper.h:186
void getSplineCoeff_SepMany(int splineindex, M3::float_t *&xArray, M3::float_t *&manyArray)
std::string getDimLabel(const int BinningOpt, const unsigned int Axis) const
void PrepareSplineFile(std::string FileName) override
KS: Prepare spline file that can be used for fast loading.
bool * isflatarray
Need to keep track of which splines are flat and which aren't.
void LoadIndexDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Index.
int CountNumberOfLoadedSplines(bool NonFlat=false, int Verbosity=0)
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.
std::vector< std::vector< std::string > > SplineFileParPrefixNames
void CalcSplineWeights() override
CPU based code which eval weight for each spline.
void LoadSettingsDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Settings.
std::vector< int > coeffindexvec
void cleanUpMemory()
Remove setup variables not needed for spline evaluations.
void PrintBinning(TAxis *Axis) const
ParameterHandlerGeneric * xsec
Pointer to covariance from which we get information about spline params.
void InvestigateMissingSplines() const
This function will find missing splines in file.
void Evaluate() override
CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; proba...
std::vector< int > nSplineParams
void ModifyWeights() override
Calc total event weight, not used by Bin-by-bin splines.
std::vector< std::vector< std::vector< int > > > SplineModeVecs
std::vector< int > nOscChans
M3::float_t * xcoeff_arr
x coefficients for each spline
std::vector< std::string > SampleTitles
void LoadSplineFile(std::string FileName) override
KS: Load preprocessed spline file.
void PrintSampleDetails(const std::string &SampleTitle) const
Print info like Sample ID of spline params etc.
std::vector< std::string > UniqueSystNames
name of each spline parameter
std::vector< std::vector< SplineInterpolation > > SplineInterpolationTypes
spline interpolation types for each sample. These vectors are from a call to GetSplineInterpolationFr...
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.
std::vector< int > Dimensions
std::vector< int > uniquesplinevec_Monolith
Maps single spline object with single parameter.
void PrepareIndexDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Index Info within SplineFile.
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< M3::float_t > weightvec_Monolith
Stores weight from spline evaluation for each single spline.
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
std::vector< TAxis * > FindSplineBinning(const std::string &FileName, const std::string &SampleTitle)
Grab histograms with spline binning.
std::vector< std::vector< std::string > > DimensionLabels
void PrepareMonolithDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Monolith Info within SplineFile.
BinnedSplineHandler(ParameterHandlerGeneric *xsec_, MaCh3Modes *Modes_)
Constructor.
void LoadMonolithDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed Monolith.
std::vector< std::vector< std::vector< TAxis * > > > SplineBinning
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.
virtual std::vector< std::string > GetTokensFromSplineName(std::string FullSplineName)=0
std::vector< TSpline3_red * > splinevec_Monolith
holds each spline object before stripping into coefficient monolith
std::vector< int > uniquecoeffindices
Unique coefficient indices.
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.
MaCh3Modes * Modes
pointer to MaCh3 Mode from which we get spline suffix
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...
virtual ~BinnedSplineHandler()
Destructor.
std::vector< std::string > SampleNames
void TransferToMonolith()
flatten multidimensional spline array into proper monolith
std::vector< int > UniqueSystIndices
Global index of each spline param, it allows us to match spline ordering with global.
int GetSampleIndex(const std::string &SampleTitle) const
Get index of sample based on name.
M3::float_t * manycoeff_arr
ybcd coefficients for each spline
Custom exception class used throughout MaCh3.
KS: Class describing MaCh3 modes used in the analysis, it is being initialised from config.
Definition: MaCh3Modes.h:135
int GetNModes() const
KS: Get number of modes, keep in mind actual number is +1 greater due to unknown category.
Definition: MaCh3Modes.h:148
std::string GetSplineSuffixFromMaCh3Mode(const int Index)
DB: Get binned spline mode suffic from MaCh3 Mode.
Definition: MaCh3Modes.cpp:243
const double * RetPointer(const int iParam)
DB Pointer return to param position.
YAML::Node GetConfig() const
Getter to return a copy of the YAML node.
Class responsible for handling of systematic error parameters with different types defined in the con...
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 > GetParsNamesFromSampleName(const std::string &SampleName, const SystType Type)
DB Grab the parameter names 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.
Base class for calculating weight from spline.
Definition: SplineBase.h:25
short int nParams
Number of parameters that have splines.
Definition: SplineBase.h:81
void FindSplineSegment()
CW:Code used in step by step reweighting, Find Spline Segment for each param.
Definition: SplineBase.cpp:24
short int * SplineSegments
Definition: SplineBase.h:77
std::vector< FastSplineInfo > SplineInfoArray
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:79
void LoadFastSplineInfoDir(std::unique_ptr< TFile > &SplineFile)
KS: Load preprocessed FastSplineInfo.
Definition: SplineBase.cpp:159
void PrepareFastSplineInfoDir(std::unique_ptr< TFile > &SplineFile) const
KS: Prepare Fast Spline Info within SplineFile.
Definition: SplineBase.cpp:131
CW: Reduced TSpline3 class.
double float_t
Definition: Core.h:37
constexpr static const double _BAD_DOUBLE_
Default value used for double initialisation.
Definition: Core.h:53
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 char * float_t_str_repr
Definition: Core.h:40
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55
constexpr T fmaf_t(T x, T y, T z)
Function template for fused multiply-add.
Definition: Core.h:45
void AddPath(std::string &FilePath)
Prepends the MACH3 environment path to FilePath if it is not already present.
Definition: Monitor.cpp:382
int int_t
Definition: Core.h:38