MaCh3  2.5.0
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 // *******************
101 // *******************
103  std::string name;
104 
106  bool hasKinBounds = false;
109  std::vector< std::vector< std::vector<double> > > Selection;
114  std::vector< std::string > KinematicVarStr;
115 
118 };
119 
120 // *******************
126  std::vector<int> modes;
128  std::vector<int> pdgs;
130  std::vector<int> preoscpdgs;
132  std::vector<int> targets;
133 };
134 
136 using FuncParFuncType = std::function<void (const M3::float_t*, std::size_t)>;
137 // *******************
142 // *******************
144  std::vector<int> modes;
145 
147  const M3::float_t* valuePtr = nullptr;
148 
151 };
152 
158 };
159 
169 };
170 
171 // **************************************************
174 inline std::string GetTF1(const SplineInterpolation i) {
175  // **************************************************
176  std::string Func = "";
177  switch(i) {
179  Func = "([1]+[0]*x)";
180  break;
187  MACH3LOG_ERROR("Interpolation type {} not supported for TF1!", static_cast<int>(i));
188  throw MaCh3Exception(__FILE__, __LINE__);
189  default:
190  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
191  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
192  throw MaCh3Exception(__FILE__ , __LINE__ );
193  }
194  return Func;
195 }
196 
197 // **************************************************
201 // **************************************************
203  switch(i) {
210  break;
212  Type = RespFuncType::kTF1_red;
213  break;
215  MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
216  throw MaCh3Exception(__FILE__, __LINE__);
217  default:
218  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
219  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
220  throw MaCh3Exception(__FILE__, __LINE__);
221  }
222  return Type;
223 }
224 
225 // **************************************************
228 // **************************************************
229  std::string name = "";
230  switch(i) {
231  // TSpline3 (third order spline in ROOT)
233  name = "TSpline3";
234  break;
236  name = "Linear";
237  break;
239  name = "Monotonic";
240  break;
241  // (Experimental) Akima_Spline (crd order spline which is allowed to be discontinuous in 2nd deriv)
243  name = "Akima";
244  break;
246  name = "KochanekBartels";
247  break;
249  name = "LinearFunc";
250  break;
252  MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
253  throw MaCh3Exception(__FILE__, __LINE__);
254  default:
255  MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
256  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
257  throw MaCh3Exception(__FILE__ , __LINE__ );
258  }
259  return name;
260 }
261 
264 enum SystType {
269  kSystTypes
270 };
271 
272 // *******************
275 // *******************
278 
280  std::vector<int> _fSplineModes;
281 
283  std::string _fSplineNames;
284 
289 };
290 
291 // **************************************************
293 inline std::string SystType_ToString(const SystType i) {
294 // **************************************************
295  std::string name = "";
296  switch(i) {
297  case SystType::kNorm:
298  name = "Norm";
299  break;
300  case SystType::kSpline:
301  name = "Spline";
302  break;
303  case SystType::kFunc:
304  name = "Functional";
305  break;
306  case SystType::kOsc:
307  name = "Oscillation";
308  break;
310  MACH3LOG_ERROR("kSystTypes is not a valid SystType!");
311  throw MaCh3Exception(__FILE__, __LINE__);
312  default:
313  MACH3LOG_ERROR("UNKNOWN SYST TYPE SPECIFIED!");
314  MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
315  throw MaCh3Exception(__FILE__ , __LINE__ );
316  }
317  return name;
318 }
319 
320 // *******************
324 // *******************
325 };
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:141
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::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.
std::function< void(const M3::float_t *, std::size_t)> FuncParFuncType
HH - a shorthand type for funcpar functions.
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.
double float_t
Definition: Core.h:37
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.
const M3::float_t * valuePtr
Parameter value pointer.
FuncParFuncType * funcPtr
Function pointer.
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.
std::vector< int > modes
Mode which parameter applies to.
std::vector< int > pdgs
PDG which parameter applies to.
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)
std::string _fSplineNames
Name of spline in TTree (TBranch),.
Base class storing info for parameters types, helping unify codebase.
int index
Parameter number of this normalisation in current systematic model.
bool hasKinBounds
Does this parameter have kinematic bounds.
std::string name
Name of parameters.
std::vector< std::string > KinematicVarStr
std::vector< std::vector< std::vector< double > > > Selection