MaCh3 DUNE 1.0.0
Reference Guide
Loading...
Searching...
No Matches
samplePDFDUNEBeamFD.h
Go to the documentation of this file.
1#ifndef _samplePDFDUNEBeamFD_h_
2#define _samplePDFDUNEBeamFD_h_
3
4#include <iostream>
5#include <TTree.h>
6#include <TH1D.h>
7#include <TH2D.h>
8#include <TMath.h>
9#include <TFile.h>
10#include <TGraph2DErrors.h>
11#include <vector>
12#include <omp.h>
13#include <list>
14
15#include "splines/splinesDUNE.h"
16#include "covariance/covarianceXsec.h"
17#include "covariance/covarianceOsc.h"
18#include "samplePDF/samplePDFFDBase.h"
19
20#include "StructsDUNE.h"
22class samplePDFDUNEBeamFD : virtual public samplePDFFDBase
23{
24public:
25
26
31 samplePDFDUNEBeamFD(std::string mc_version, covarianceXsec* xsec_cov, covarianceOsc* osc_cov);
32
35
37
38 //More robust getters to make plots in different variables, mode, osc channel, systematic weighting and with bin range
39
47 TH1D* get1DVarHist(KinematicTypes Var1, int fModeToFill=-1, int fSampleToFill=-1, int WeightStyle=0, TAxis* Axis=0);
48
55 TH1D* get1DVarHist(KinematicTypes Var1, std::vector< std::vector<double> > Selection, int WeightStyle=0, TAxis* Axis=0);
56
57 protected:
59 void Init();
60
64 int setupExperimentMC(int iSample);
65
68 void setupFDMC(int iSample);
69
72
75 void SetupSplines();
76
82 double ReturnKinematicParameter (double KinematicVariable, int iSample, int iEvent);
83
89 double ReturnKinematicParameter(std::string KinematicParameter, int iSample, int iEvent);
90
96 const double* GetPointerToKinematicParameter(std::string KinematicParameter, int iSample, int iEvent);
97
103 const double* GetPointerToKinematicParameter(double KinematicVariable, int iSample, int iEvent);
104
105 std::vector<double> ReturnKinematicParameterBinning(std::string KinematicParameter);
106 inline int ReturnKinematicParameterFromString(std::string KinematicParameterStr);
107 inline std::string ReturnStringFromKinematicParameter(int KinematicParameterStr);
108
109 //DB functions which could be initialised to do something which is non-trivial
113 double CalcXsecWeightFunc(int iSample, int iEvent) {return 1.;}
114
118 void applyShifts(int iSample, int iEvent);
119
120 // dunemc
122 std::vector<struct dunemc_base> dunemcSamples;
123
125 double pot;
126
129
131 TTree *_data;
132
133 TString _nutype;
134
135 int _mode;
136
137 //Reco Variables
138 double _erec;
139 double _erec_nue;
140 double _erec_had;
142 double _erec_lep;
144
145 double _eRecoP;
146 double _eRecoPip;
147 double _eRecoPim;
148 double _eRecoPi0;
149 double _eRecoN;
150
151 double _cvnnumu;
152 double _cvnnue;
153 double _vtx_x;
154 double _vtx_y;
155 double _vtx_z;
156
157 //Truth Variables
158 double _ev;
159 double _LepE;
160 double _eP;
161 double _ePip;
162 double _ePim;
163 double _ePi0;
164 double _eN;
165 double _NuMomX;
166 double _NuMomY;
167 double _NuMomZ;
168 double _LepMomX;
169 double _LepMomY;
170 double _LepMomZ;
175 int _isCC;
178 int _run;
180 double _LepTheta;
181 double _Q2;
182
184
185 //Positions of FD Detector systematics
207
209 std::vector<const double*> FDDetectorSystPointers;
210
213};
214
215#endif
Base class for handling FD Beam samples.
std::vector< const double * > FDDetectorSystPointers
FD Detector Systematics.
void SetupSplines()
Setup our spline file, this calls InitialseSplineObject() under the hood.
void applyShifts(int iSample, int iEvent)
Apply kinematic shifts.
TFile * _sampleFile
File containing sample objects.
samplePDFDUNEBeamFD(std::string mc_version, covarianceXsec *xsec_cov, covarianceOsc *osc_cov)
samplePDFDUNE FD beam Constructor
double CalcXsecWeightFunc(int iSample, int iEvent)
NOT IMPLEMENTED: Dunder method to calculate xsec weights.
std::vector< struct dunemc_base > dunemcSamples
DUNE MC sampels.
std::string ReturnStringFromKinematicParameter(int KinematicParameterStr)
TTree * _data
TTree containing sample Data.
TH1D * get1DVarHist(KinematicTypes Var1, int fModeToFill=-1, int fSampleToFill=-1, int WeightStyle=0, TAxis *Axis=0)
Getter to make plots in different variables, mode, osc channel, systematic weighting and with bin ran...
const double * GetPointerToKinematicParameter(std::string KinematicParameter, int iSample, int iEvent)
Returns pointer to kinemtatic parameter for event in Structs DUNE.
int ReturnKinematicParameterFromString(std::string KinematicParameterStr)
void SetupWeightPointers()
Sets up pointers weights for each event (oscillation/xsec/etc.)
double pot
Value of POT used for sample.
int setupExperimentMC(int iSample)
Function to setup MC from file.
int nFDDetectorSystPointers
Number of FD Detector Systematics.
void setupFDMC(int iSample)
Tells FD base which variables to point to/be set to.
double ReturnKinematicParameter(double KinematicVariable, int iSample, int iEvent)
Returns pointer to kinemtatic parameter for event in Structs DUNE.
void Init()
Initialises object.
std::vector< double > ReturnKinematicParameterBinning(std::string KinematicParameter)