MaCh3  2.4.2
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 
31 
32 // *******************
34 template< typename T, size_t N >
35 std::vector<T> MakeVector( const T (&data)[N] ) {
36 // *******************
37  return std::vector<T>(data, data+N);
38 }
39 
41 template <typename T>
42 void CleanVector(T&) {}
43 
44 // *******************
51 template <typename T>
52 void CleanVector(std::vector<T>& v) {
53 // *******************
54  for (auto& x : v)
55  CleanVector(x);
56 
57  v.clear();
58  v.shrink_to_fit();
59 }
60 
61 // *******************
63 template <typename T>
64 void CleanContainer(T&) {}
65 // *******************
66 
67 // *******************
71 template <typename T>
72 void CleanContainer(T*& ptr) {
73  // *******************
74  delete ptr;
75  ptr = nullptr;
76 }
77 
78 // *******************
85 template <typename T>
86 void CleanContainer(std::vector<T>& container) {
87  // *******************
88  for (auto& element : container)
89  CleanContainer(element);
90 
91  container.clear();
92  container.shrink_to_fit();
93 }
94 
95 // *******************
99 // *******************
101  std::string name;
102 
105 };
106 
107 // *******************
113  std::vector<int> modes;
115  std::vector<int> pdgs;
117  std::vector<int> preoscpdgs;
119  std::vector<int> targets;
124  std::vector< std::vector< std::vector<double> > > Selection;
125 
130  std::vector< std::string > KinematicVarStr;
131 };
132 
134 using FuncParFuncType = std::function<void (const double*, std::size_t)>;
135 // *******************
140 // *******************
142  std::vector<int> modes;
144  std::vector<int> pdgs;
146  std::vector<int> preoscpdgs;
148  std::vector<int> targets;
153  std::vector< std::vector< std::vector<double> > > Selection;
154 
159  std::vector< std::string > KinematicVarStr;
160 
162  const double* valuePtr;
163 
166 };
167 
173 };
174 
184 };
185 
186 // **************************************************
189 inline std::string GetTF1(const SplineInterpolation i) {
190  // **************************************************
191  std::string Func = "";
192  switch(i) {
194  Func = "([1]+[0]*x)";
195  break;
202  MACH3LOG_ERROR("Interpolation type {} not supported for TF1!", static_cast<int>(i));
203  throw MaCh3Exception(__FILE__, __LINE__);
204  default:
205  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
206  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
207  throw MaCh3Exception(__FILE__ , __LINE__ );
208  }
209  return Func;
210 }
211 
212 // **************************************************
216 // **************************************************
218  switch(i) {
225  break;
227  Type = RespFuncType::kTF1_red;
228  break;
230  MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
231  throw MaCh3Exception(__FILE__, __LINE__);
232  default:
233  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
234  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
235  throw MaCh3Exception(__FILE__, __LINE__);
236  }
237  return Type;
238 }
239 
240 // **************************************************
243 // **************************************************
244  std::string name = "";
245  switch(i) {
246  // TSpline3 (third order spline in ROOT)
248  name = "TSpline3";
249  break;
251  name = "Linear";
252  break;
254  name = "Monotonic";
255  break;
256  // (Experimental) Akima_Spline (crd order spline which is allowed to be discontinuous in 2nd deriv)
258  name = "Akima";
259  break;
261  name = "KochanekBartels";
262  break;
264  name = "LinearFunc";
265  break;
267  MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
268  throw MaCh3Exception(__FILE__, __LINE__);
269  default:
270  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
271  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
272  throw MaCh3Exception(__FILE__ , __LINE__ );
273  }
274  return name;
275 }
276 
279 enum SystType {
284  kSystTypes
285 };
286 
287 // *******************
290 // *******************
293 
295  std::vector<int> _fSplineModes;
296 
301 };
302 
303 // **************************************************
305 inline std::string SystType_ToString(const SystType i) {
306 // **************************************************
307  std::string name = "";
308  switch(i) {
309  case SystType::kNorm:
310  name = "Norm";
311  break;
312  case SystType::kSpline:
313  name = "Spline";
314  break;
315  case SystType::kFunc:
316  name = "Functional";
317  break;
318  case SystType::kOsc:
319  name = "Oscillation";
320  break;
322  MACH3LOG_ERROR("kSystTypes is not a valid SystType!");
323  throw MaCh3Exception(__FILE__, __LINE__);
324  default:
325  MACH3LOG_ERROR("UNKNOWN SYST TYPE SPECIFIED!");
326  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
327  throw MaCh3Exception(__FILE__ , __LINE__ );
328  }
329  return name;
330 }
331 
332 // *******************
336 // *******************
337 };
KS: Core MaCh3 definitions and compile-time configuration utilities.
#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
Defines the custom exception class used throughout MaCh3.
MaCh3 Logging utilities built on top of SPDLOG.
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
std::string SplineInterpolation_ToString(const SplineInterpolation i)
Convert a LLH type to a string.
void CleanContainer(T &)
Base case: do nothing for non-pointer types.
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.
std::string GetTF1(const SplineInterpolation i)
Get function for TF1_red.
std::function< void(const double *, std::size_t)> FuncParFuncType
HH - a shorthand type for funcpar functions.
std::vector< T > MakeVector(const T(&data)[N])
Template to make vector out of an array of any length.
void CleanVector(T &)
Base case: do nothing for non-vector types.
std::string SystType_ToString(const SystType i)
Convert a Syst type type to a string.
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 used throughout MaCh3.
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55
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 > 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 > 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.