MaCh3 2.2.1
Reference Guide
Loading...
Searching...
No Matches
Classes | Typedefs | Enumerations | Functions
ParameterStructs.h File Reference
#include <set>
#include <list>
#include <unordered_map>
#include "Manager/MaCh3Exception.h"
#include "Manager/MaCh3Logger.h"
#include "Manager/Core.h"
#include "TSpline.h"
#include "TObjString.h"
#include "TFile.h"
#include "TF1.h"
#include "TH2Poly.h"
#include "TH1.h"
Include dependency graph for ParameterStructs.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  TypeParameterBase
 Base class storing info for parameters types, helping unify codebase. More...
 
struct  NormParameter
 ETA - Normalisations for cross-section parameters Carrier for whether you want to apply a systematic to an event or not. More...
 
struct  FunctionalParameter
 HH - Functional parameters Carrier for whether you want to apply a systematic to an event or not. More...
 
struct  SplineParameter
 KS: Struct holding info about Spline Systematics. More...
 
struct  OscillationParameter
 KS: Struct holding info about oscillation Systematics. More...
 

Typedefs

using FuncParFuncType = std::function< void(const double *, std::size_t)>
 

Enumerations

enum  RespFuncType { kTSpline3_red , kTF1_red , kRespFuncTypes }
 Make an enum of the spline interpolation type. More...
 
enum  SplineInterpolation {
  kTSpline3 , kLinear , kMonotonic , kAkima ,
  kLinearFunc , kSplineInterpolations
}
 Make an enum of the spline interpolation type. More...
 
enum  SystType {
  kNorm , kSpline , kFunc , kOsc ,
  kSystTypes
}
 

Functions

template<typename T , size_t N>
std::vector< T > MakeVector (const T(&data)[N])
 Template to make vector out of an array of any length.
 
template<typename T >
void CleanVector (std::vector< T > &vec)
 Generic cleanup function.
 
template<typename T >
void CleanVector (std::vector< std::vector< T > > &vec)
 Generic cleanup function.
 
template<typename T >
void CleanVector (std::vector< std::vector< std::vector< T > > > &vec)
 Generic cleanup function.
 
template<typename T >
void CleanVector (std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< T > > > > > > > &vec)
 Generic cleanup function.
 
template<typename T >
void CleanContainer (std::vector< T * > &container)
 Generic cleanup function.
 
template<typename T >
void CleanContainer (std::vector< std::vector< std::vector< T * > > > &container)
 Generic cleanup function.
 
std::string GetTF1 (const SplineInterpolation i)
 Get function for TF1_red.
 
RespFuncType SplineInterpolation_ToRespFuncType (const SplineInterpolation i)
 Convert a RespFuncType type to a SplineInterpolation.
 
std::string SplineInterpolation_ToString (const SplineInterpolation i)
 Convert a LLH type to a string.
 
std::string SystType_ToString (const SystType i)
 Convert a Syst type type to a string.
 

Detailed Description

Author
Asher Kaboth
Clarence Wret
Patrick Dunne
Dan Barrow
Ed Atkin
Kamil Skwarczynski

Definition in file ParameterStructs.h.

Typedef Documentation

◆ FuncParFuncType

using FuncParFuncType = std::function<void (const double*, std::size_t)>

Definition at line 184 of file ParameterStructs.h.

Enumeration Type Documentation

◆ RespFuncType

Make an enum of the spline interpolation type.

Enumerator
kTSpline3_red 

Uses TSpline3_red for interpolation.

kTF1_red 

Uses TF1_red for interpolation.

kRespFuncTypes 

This only enumerates.

Definition at line 220 of file ParameterStructs.h.

220 {
222 kTF1_red,
224};
@ kTF1_red
Uses TF1_red for interpolation.
@ kTSpline3_red
Uses TSpline3_red for interpolation.
@ kRespFuncTypes
This only enumerates.

◆ SplineInterpolation

Make an enum of the spline interpolation type.

Enumerator
kTSpline3 

Default TSpline3 interpolation.

kLinear 

Linear interpolation between knots.

kMonotonic 

EM: DOES NOT make the entire spline monotonic, only the segments.

kAkima 

EM: Akima spline iis allowed to be discontinuous in 2nd derivative and coefficients in any segment.

kLinearFunc 

Liner interpolation using TF1 not spline.

kSplineInterpolations 

This only enumerates.

Definition at line 227 of file ParameterStructs.h.

227 {
228 kTSpline3,
229 kLinear,
230 kMonotonic,
231 kAkima,
234};
@ kTSpline3
Default TSpline3 interpolation.
@ kMonotonic
EM: DOES NOT make the entire spline monotonic, only the segments.
@ kSplineInterpolations
This only enumerates.
@ 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.

◆ SystType

enum SystType

Make an enum of systematic type recognised by covariance class

Todo:
KS: Consider using enum class, it is generally recommended as safer. It will require many static_cast
Enumerator
kNorm 

For normalisation parameters.

kSpline 

For splined parameters (1D)

kFunc 

For functional parameters.

kOsc 

For oscillation parameters.

kSystTypes 

This only enumerates.

Definition at line 324 of file ParameterStructs.h.

324 {
325 kNorm,
326 kSpline,
327 kFunc,
328 kOsc,
330};
@ kNorm
For normalisation parameters.
@ kSpline
For splined parameters (1D)
@ kSystTypes
This only enumerates.
@ kOsc
For oscillation parameters.
@ kFunc
For functional parameters.

Function Documentation

◆ CleanContainer() [1/2]

template<typename T >
void CleanContainer ( std::vector< std::vector< std::vector< T * > > > &  container)

Generic cleanup function.

Definition at line 128 of file ParameterStructs.h.

128 {
129// *******************
130 for (auto& container2D : container) {
131 for (auto& container1D : container2D) {
132 for (T* ptr : container1D) {
133 if(ptr) delete ptr;
134 }
135 container1D.clear();
136 container1D.shrink_to_fit();
137 }
138 container2D.clear();
139 container2D.shrink_to_fit();
140 }
141 container.clear();
142 container.shrink_to_fit();
143}

◆ CleanContainer() [2/2]

template<typename T >
void CleanContainer ( std::vector< T * > &  container)

Generic cleanup function.

Definition at line 115 of file ParameterStructs.h.

115 {
116// *******************
117 for (T* ptr : container) {
118 if(ptr) delete ptr;
119 ptr = nullptr;
120 }
121 container.clear();
122 container.shrink_to_fit();
123}

◆ CleanVector() [1/4]

template<typename T >
void CleanVector ( std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< std::vector< T > > > > > > > &  vec)

Generic cleanup function.

Todo:
Use recursive to make it more scalable in future

Definition at line 82 of file ParameterStructs.h.

82 {
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}

◆ CleanVector() [2/4]

template<typename T >
void CleanVector ( std::vector< std::vector< std::vector< T > > > &  vec)

Generic cleanup function.

Definition at line 64 of file ParameterStructs.h.

64 {
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}

◆ CleanVector() [3/4]

template<typename T >
void CleanVector ( std::vector< std::vector< T > > &  vec)

Generic cleanup function.

Definition at line 51 of file ParameterStructs.h.

51 {
52// *******************
53 for (auto& innerVec : vec) {
54 innerVec.clear();
55 innerVec.shrink_to_fit();
56 }
57 vec.clear();
58 vec.shrink_to_fit();
59}

◆ CleanVector() [4/4]

template<typename T >
void CleanVector ( std::vector< T > &  vec)

Generic cleanup function.

Definition at line 42 of file ParameterStructs.h.

42 {
43// *******************
44 vec.clear();
45 vec.shrink_to_fit();
46}

◆ GetTF1()

std::string GetTF1 ( const SplineInterpolation  i)
inline

Get function for TF1_red.

Parameters
iInterpolation type

Definition at line 239 of file ParameterStructs.h.

239 {
240 // **************************************************
241 std::string Func = "";
242 switch(i) {
244 Func = "([1]+[0]*x)";
245 break;
251 MACH3LOG_ERROR("Interpolation type {} not supported for TF1!", static_cast<int>(i));
252 throw MaCh3Exception(__FILE__, __LINE__);
253 default:
254 MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
255 MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
256 throw MaCh3Exception(__FILE__ , __LINE__ );
257 }
258 return Func;
259}
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:25
Custom exception class for MaCh3 errors.

◆ MakeVector()

template<typename T , size_t N>
std::vector< T > MakeVector ( const T(&)  data[N])

Template to make vector out of an array of any length.

Definition at line 34 of file ParameterStructs.h.

34 {
35// *******************
36 return std::vector<T>(data, data+N);
37}

◆ SplineInterpolation_ToRespFuncType()

RespFuncType SplineInterpolation_ToRespFuncType ( const SplineInterpolation  i)
inline

Convert a RespFuncType type to a SplineInterpolation.

Parameters
iInterpolation type

Definition at line 264 of file ParameterStructs.h.

264 {
265// **************************************************
267 switch(i) {
273 break;
276 break;
278 MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
279 throw MaCh3Exception(__FILE__, __LINE__);
280 default:
281 MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
282 MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
283 throw MaCh3Exception(__FILE__, __LINE__);
284 }
285 return Type;
286}
RespFuncType
Make an enum of the spline interpolation type.

◆ SplineInterpolation_ToString()

std::string SplineInterpolation_ToString ( const SplineInterpolation  i)
inline

Convert a LLH type to a string.

Definition at line 290 of file ParameterStructs.h.

290 {
291// **************************************************
292 std::string name = "";
293 switch(i) {
294 // TSpline3 (third order spline in ROOT)
296 name = "TSpline3";
297 break;
299 name = "Linear";
300 break;
302 name = "Monotonic";
303 break;
304 // (Experimental) Akima_Spline (crd order spline which is allowed to be discontinuous in 2nd deriv)
306 name = "Akima";
307 break;
309 name = "LinearFunc";
310 break;
312 MACH3LOG_ERROR("kSplineInterpolations is not a valid interpolation type!");
313 throw MaCh3Exception(__FILE__, __LINE__);
314 default:
315 MACH3LOG_ERROR("UNKNOWN SPLINE INTERPOLATION SPECIFIED!");
316 MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
317 throw MaCh3Exception(__FILE__ , __LINE__ );
318 }
319 return name;
320}

◆ SystType_ToString()

std::string SystType_ToString ( const SystType  i)
inline

Convert a Syst type type to a string.

Definition at line 350 of file ParameterStructs.h.

350 {
351// **************************************************
352 std::string name = "";
353 switch(i) {
354 case SystType::kNorm:
355 name = "Norm";
356 break;
358 name = "Spline";
359 break;
360 case SystType::kFunc:
361 name = "Functional";
362 break;
363 case SystType::kOsc:
364 name = "Oscillation";
365 break;
367 MACH3LOG_ERROR("kSystTypes is not a valid SystType!");
368 throw MaCh3Exception(__FILE__, __LINE__);
369 default:
370 MACH3LOG_ERROR("UNKNOWN SYST TYPE SPECIFIED!");
371 MACH3LOG_ERROR("You gave {}", static_cast<int>(i));
372 throw MaCh3Exception(__FILE__ , __LINE__ );
373 }
374 return name;
375}