MaCh3 2.2.1
Reference Guide
Loading...
Searching...
No Matches
SplineBase.h
Go to the documentation of this file.
1#pragma once
2
3// C++ includes
4#include <cstdlib>
5#include <iomanip>
6
7// MaCh3 includes
10
12// ROOT include
13#include "TFile.h"
14#include "TH1F.h"
15#include "TKey.h"
16#include "TString.h"
17#include "TIterator.h"
18#include "TStopwatch.h"
19#include "TCanvas.h"
20#include "TStyle.h"
21#include "TTree.h"
23
26 public:
28 SplineBase();
30 virtual ~SplineBase();
31
33 virtual void Evaluate() = 0;
34
36 virtual inline std::string GetName()const {return "SplineBase";};
37
39 short int GetNParams()const {return nParams;};
40
41 protected:
43 void FindSplineSegment();
45 virtual void CalcSplineWeights() = 0;
47 virtual void ModifyWeights() = 0;
48
53 void getTF1Coeff(TF1_red* &spl, int &nPoints, float *&coeffs);
54
57 std::vector<FastSplineInfo> SplineInfoArray;
60 short int *SplineSegments;
64 short int nParams;
65};
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:106
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:117
Contains definitions for spline coefficients and structure used in both CPU and GPU code.
Contains structures and helper functions for handling spline representations of systematic parameters...
Base class for calculating weight from spline.
Definition: SplineBase.h:25
virtual ~SplineBase()
Destructor.
Definition: SplineBase.cpp:16
virtual void CalcSplineWeights()=0
CPU based code which eval weight for each spline.
short int nParams
Number of parameters that have splines.
Definition: SplineBase.h:64
virtual std::string GetName() const
Get class name.
Definition: SplineBase.h:36
virtual void ModifyWeights()=0
Calc total event weight.
void FindSplineSegment()
CW:Code used in step by step reweighting, Find Spline Segment for each param.
Definition: SplineBase.cpp:24
virtual void Evaluate()=0
CW: This Eval should be used when using two separate x,{y,a,b,c,d} arrays to store the weights; proba...
short int * SplineSegments
Definition: SplineBase.h:60
std::vector< FastSplineInfo > SplineInfoArray
Definition: SplineBase.h:57
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:62
SplineBase()
Constructor.
Definition: SplineBase.cpp:7
void getTF1Coeff(TF1_red *&spl, int &nPoints, float *&coeffs)
CW: Gets the polynomial coefficients for TF1.
Definition: SplineBase.cpp:106
short int GetNParams() const
Get number of spline parameters.
Definition: SplineBase.h:39
CW: A reduced TF1 class only. Only saves parameters for each TF1 and how many parameters each paramet...