MaCh3  2.2.3
Reference Guide
ParameterStructs.h
Go to the documentation of this file.
1 #pragma once
2 
3 // C++ includes
4 #include <set>
5 #include <list>
6 #include <unordered_map>
7 
8 // MaCh3 includes
10 #include "Manager/MaCh3Logger.h"
11 #include "Manager/Core.h"
12 
14 // ROOT include
15 #include "TSpline.h"
16 #include "TObjString.h"
17 #include "TFile.h"
18 #include "TF1.h"
19 #include "TH2Poly.h"
20 #include "TH1.h"
22 
30 
31 // *******************
33 template< typename T, size_t N >
34 std::vector<T> MakeVector( const T (&data)[N] ) {
35 // *******************
36  return std::vector<T>(data, data+N);
37 }
38 
39 // *******************
41 template <typename T>
42 void CleanVector(std::vector<T>& vec) {
43 // *******************
44  vec.clear();
45  vec.shrink_to_fit();
46 }
47 
48 // *******************
50 template <typename T>
51 void CleanVector(std::vector<std::vector<T>>& vec) {
52 // *******************
53  for (auto& innerVec : vec) {
54  innerVec.clear();
55  innerVec.shrink_to_fit();
56  }
57  vec.clear();
58  vec.shrink_to_fit();
59 }
60 
61 // *******************
63 template <typename T>
64 void CleanVector(std::vector<std::vector<std::vector<T>>>& vec) {
65 // *******************
66  for (auto& inner2DVec : vec) {
67  for (auto& innerVec : inner2DVec) {
68  innerVec.clear();
69  innerVec.shrink_to_fit();
70  }
71  inner2DVec.clear();
72  inner2DVec.shrink_to_fit();
73  }
74  vec.clear();
75  vec.shrink_to_fit();
76 }
77 
78 // *******************
81 template <typename T>
82 void CleanVector(std::vector<std::vector<std::vector<std::vector<std::vector<std::vector<std::vector<T>>>>>>> &vec) {
83 // *******************
84  for (auto& v6 : vec) {
85  for (auto& v5 : v6) {
86  for (auto& v4 : v5) {
87  for (auto& v3 : v4) {
88  for (auto& v2 : v3) {
89  for (auto& v1 : v2) {
90  v1.clear();
91  v1.shrink_to_fit();
92  }
93  v2.clear();
94  v2.shrink_to_fit();
95  }
96  v3.clear();
97  v3.shrink_to_fit();
98  }
99  v4.clear();
100  v4.shrink_to_fit();
101  }
102  v5.clear();
103  v5.shrink_to_fit();
104  }
105  v6.clear();
106  v6.shrink_to_fit();
107  }
108  vec.clear();
109  vec.shrink_to_fit();
110 }
111 
112 // *******************
114 template <typename T>
115 void CleanContainer(std::vector<T*>& container) {
116 // *******************
117  for (T* ptr : container) {
118  if(ptr) delete ptr;
119  ptr = nullptr;
120  }
121  container.clear();
122  container.shrink_to_fit();
123 }
124 
125 // *******************
127 template <typename T>
128 void CleanContainer(std::vector<std::vector<T*>>& container) {
129 // *******************
130  for (auto& inner : container) {
131  CleanContainer(inner);
132  }
133  container.clear();
134  container.shrink_to_fit();
135 }
136 
137 // *******************
139 template <typename T>
140 void CleanContainer(std::vector< std::vector< std::vector<T*> > >& container) {
141 // *******************
142  for (auto& container2D : container) {
143  for (auto& container1D : container2D) {
144  for (T* ptr : container1D) {
145  if(ptr) delete ptr;
146  }
147  container1D.clear();
148  container1D.shrink_to_fit();
149  }
150  container2D.clear();
151  container2D.shrink_to_fit();
152  }
153  container.clear();
154  container.shrink_to_fit();
155 }
156 
157 // *******************
160 // *******************
162  std::string name;
163 
166 };
167 
168 // *******************
173  std::vector<int> modes;
175  std::vector<int> horncurrents;
177  std::vector<int> pdgs;
179  std::vector<int> preoscpdgs;
181  std::vector<int> targets;
186  std::vector< std::vector< std::vector<double> > > Selection;
187 
192  std::vector< std::string > KinematicVarStr;
193 };
194 
195 // HH - a shorthand type for funcpar functions
196 using FuncParFuncType = std::function<void (const double*, std::size_t)>;
197 // *******************
201 // *******************
203  std::vector<int> modes;
205  std::vector<int> horncurrents;
207  std::vector<int> pdgs;
209  std::vector<int> preoscpdgs;
211  std::vector<int> targets;
216  std::vector< std::vector< std::vector<double> > > Selection;
217 
222  std::vector< std::string > KinematicVarStr;
223 
225  const double* valuePtr;
226 
229 };
230 
236 };
237 
247 };
248 
249 // **************************************************
252 inline std::string GetTF1(const SplineInterpolation i) {
253  // **************************************************
254  std::string Func = "";
255  switch(i) {
257  Func = "([1]+[0]*x)";
258  break;
265  MACH3LOG_ERROR("Interpolation type {} not supported for TF1!", static_cast<int>(i));
266  throw MaCh3Exception(__FILE__, __LINE__);
267  default:
268  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
269  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
270  throw MaCh3Exception(__FILE__ , __LINE__ );
271  }
272  return Func;
273 }
274 
275 // **************************************************
279 // **************************************************
281  switch(i) {
288  break;
290  Type = RespFuncType::kTF1_red;
291  break;
293  MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
294  throw MaCh3Exception(__FILE__, __LINE__);
295  default:
296  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
297  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
298  throw MaCh3Exception(__FILE__, __LINE__);
299  }
300  return Type;
301 }
302 
303 // **************************************************
306 // **************************************************
307  std::string name = "";
308  switch(i) {
309  // TSpline3 (third order spline in ROOT)
311  name = "TSpline3";
312  break;
314  name = "Linear";
315  break;
317  name = "Monotonic";
318  break;
319  // (Experimental) Akima_Spline (crd order spline which is allowed to be discontinuous in 2nd deriv)
321  name = "Akima";
322  break;
324  name = "KochanekBartels";
325  break;
327  name = "LinearFunc";
328  break;
330  MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
331  throw MaCh3Exception(__FILE__, __LINE__);
332  default:
333  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
334  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
335  throw MaCh3Exception(__FILE__ , __LINE__ );
336  }
337  return name;
338 }
339 
342 enum SystType {
347  kSystTypes
348 };
349 
350 // *******************
353 // *******************
356 
358  std::vector<int> _fSplineModes;
359 
364 };
365 
366 // **************************************************
368 inline std::string SystType_ToString(const SystType i) {
369 // **************************************************
370  std::string name = "";
371  switch(i) {
372  case SystType::kNorm:
373  name = "Norm";
374  break;
375  case SystType::kSpline:
376  name = "Spline";
377  break;
378  case SystType::kFunc:
379  name = "Functional";
380  break;
381  case SystType::kOsc:
382  name = "Oscillation";
383  break;
385  MACH3LOG_ERROR("kSystTypes is not a valid SystType!");
386  throw MaCh3Exception(__FILE__, __LINE__);
387  default:
388  MACH3LOG_ERROR("UNKNOWN SYST TYPE SPECIFIED!");
389  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
390  throw MaCh3Exception(__FILE__ , __LINE__ );
391  }
392  return name;
393 }
394 
395 // *******************
399 // *******************
400 
401 };
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:109
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:120
KS: Based on this https://github.com/gabime/spdlog/blob/a2b4262090fd3f005c2315dcb5be2f0f1774a005/incl...
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:27
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.
SplineInterpolation
Make an enum of the spline interpolation type.
@ kTSpline3
Default TSpline3 interpolation.
@ kMonotonic
EM: DOES NOT make the entire spline monotonic, only the segments.
@ kSplineInterpolations
This only enumerates.
@ kKochanekBartels
KS: Kochanek-Bartels spline: allows local control of tension, continuity, and bias.
@ kLinear
Linear interpolation between knots.
@ kLinearFunc
Liner interpolation using TF1 not spline.
@ kAkima
EM: Akima spline iis allowed to be discontinuous in 2nd derivative and coefficients in any segment.
RespFuncType
Make an enum of the spline interpolation type.
@ kTF1_red
Uses TF1_red for interpolation.
@ kTSpline3_red
Uses TSpline3_red for interpolation.
@ kRespFuncTypes
This only enumerates.
void CleanVector(std::vector< T > &vec)
Generic cleanup function.
std::string GetTF1(const SplineInterpolation i)
Get function for TF1_red.
std::function< void(const double *, std::size_t)> FuncParFuncType
std::vector< T > MakeVector(const T(&data)[N])
Template to make vector out of an array of any length.
std::string SystType_ToString(const SystType i)
Convert a Syst type type to a string.
void CleanContainer(std::vector< T * > &container)
Generic cleanup function.
RespFuncType SplineInterpolation_ToRespFuncType(const SplineInterpolation i)
Convert a RespFuncType type to a SplineInterpolation.
SystType
@ kNorm
For normalisation parameters.
@ kSpline
For splined parameters (1D)
@ kSystTypes
This only enumerates.
@ kOsc
For oscillation parameters.
@ kFunc
For functional parameters.
Custom exception class for MaCh3 errors.
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:48
HH - Functional parameters Carrier for whether you want to apply a systematic to an event or not.
std::vector< int > modes
Mode which parameter applies to.
std::vector< std::string > KinematicVarStr
const double * valuePtr
Parameter value pointer.
bool hasKinBounds
Does this parameter have kinematic bounds.
std::vector< int > targets
Targets which parameter applies to.
std::vector< int > horncurrents
Horn currents which parameter applies to.
std::vector< int > pdgs
PDG which parameter applies to.
FuncParFuncType * funcPtr
Function pointer.
std::vector< std::vector< std::vector< double > > > Selection
std::vector< int > preoscpdgs
Preosc PDG which parameter applies to.
ETA - Normalisations for cross-section parameters Carrier for whether you want to apply a systematic ...
std::vector< int > preoscpdgs
Preosc PDG which parameter applies to.
bool hasKinBounds
Does this parameter have kinematic bounds.
std::vector< std::vector< std::vector< double > > > Selection
std::vector< int > modes
Mode which parameter applies to.
std::vector< int > horncurrents
Horn currents which parameter applies to.
std::vector< int > pdgs
PDG which parameter applies to.
std::vector< std::string > KinematicVarStr
std::vector< int > targets
Targets which parameter applies to.
KS: Struct holding info about oscillation Systematics.
KS: Struct holding info about Spline Systematics.
SplineInterpolation _SplineInterpolationType
Spline interpolation vector.
double _SplineKnotUpBound
EM: Cap spline knot higher value.
double _SplineKnotLowBound
EM: Cap spline knot lower value.
std::vector< int > _fSplineModes
Modes to which spline applies (valid only for binned splines)
Base class storing info for parameters types, helping unify codebase.
int index
Parameter number of this normalisation in current systematic model.
std::string name
Name of parameters.