MaCh3 DUNE 1.0.0
Reference Guide
Loading...
Searching...
No Matches
samplePDFDUNEBeamND Class Reference

Base class for handling beam ND LAR samples. More...

#include <samplePDFDUNE/samplePDFDUNEBeamND.h>

Inheritance diagram for samplePDFDUNEBeamND:
Collaboration diagram for samplePDFDUNEBeamND:

Public Types

enum  KinematicTypes { kTrueNeutrinoEnergy , kRecoQ }
 Enum to identify kinematics. More...
 

Public Member Functions

 samplePDFDUNEBeamND (std::string mc_version, covarianceXsec *xsec_cov, covarianceOsc *osc_cov)
 samplePDFDUNE ND beam Constructor
 
 ~samplePDFDUNEBeamND ()
 destructor
 

Protected Member Functions

void Init ()
 Initialises object.
 
int setupExperimentMC (int iSample)
 Function to setup MC from file.
 
void setupFDMC (int iSample)
 Tells FD base which variables to point to/be set to.
 
void SetupWeightPointers ()
 Sets up pointers weights for each event (oscillation/xsec/etc.)
 
void SetupSplines ()
 Sets up splines.
 
const double * GetPointerToKinematicParameter (KinematicTypes KinPar, int iSample, int iEvent)
 Returns pointer to kinemtatic parameter for event in Structs DUNE.
 
const double * GetPointerToKinematicParameter (double KinematicVariable, int iSample, int iEvent)
 Returns pointer to kinemtatic parameter for event in Structs DUNE.
 
const double * GetPointerToKinematicParameter (std::string KinematicParameter, int iSample, int iEvent)
 Returns pointer to kinemtatic parameter for event in Structs DUNE.
 
double ReturnKinematicParameter (double KinematicVariable, int iSample, int iEvent)
 Returns pointer to kinemtatic parameter for event in Structs DUNE.
 
double ReturnKinematicParameter (std::string KinematicParameter, int iSample, int iEvent)
 Returns pointer to kinemtatic parameter for event in Structs DUNE.
 
std::vector< double > ReturnKinematicParameterBinning (std::string KinematicParameter)
 Gets binning for a given parameter.
 
std::string ReturnStringFromKinematicParameter (int KinematicParameter)
 Gets name of kinematic parmaeter.
 
int ReturnKinematicParameterFromString (std::string KinematicParameterStr)
 Get kinematic parameter ID from string name.
 
double CalcXsecWeightFunc (int iSample, int iEvent)
 NOT IMPLEMENTED: Dunder method to calculate xsec weights.
 
void applyShifts (int iSample, int iEvent)
 NOT IMPLEMENTED: Apply kinematic shifts.
 

Protected Attributes

std::vector< struct dunemc_basedunendmcSamples
 Array filled with MC samples for each oscillation channel.
 
TFile * _sampleFile
 File containing sample objects.
 
TTree * _data
 TTree containing sample Data.
 
double pot
 Value of POT used for sample.
 
TString _nutype
 
int _mode
 
double _ev
 
double _erec
 
double _erec_lep
 
double _erec_had
 
int _reco_numu
 
int _reco_nue
 
double _eRecoP
 
double _eRecoPip
 
double _eRecoPim
 
double _eRecoPi0
 
double _eRecoN
 
double _LepNuAngle
 
double _LepE
 
double _eP
 
double _ePip
 
double _ePim
 
double _ePi0
 
double _eN
 
double _BeRPA_cvwgt
 
int _isCC
 
int _nuPDGunosc
 
int _nuPDG
 
int _run
 
int _isND
 
int _isFHC
 
double _vtx_x
 
double _vtx_y
 
double _vtx_z
 
double _LepTheta
 
double _Q2
 
int _reco_q
 
bool iselike
 
bool isND
 
bool IsRHC
 
double tot_escale_nd_pos
 
double tot_escale_sqrt_nd_pos
 
double tot_escale_invsqrt_nd_pos
 
double had_escale_nd_pos
 
double had_escale_sqrt_nd_pos
 
double had_escale_invsqrt_nd_pos
 
double mu_escale_nd_pos
 
double mu_escale_sqrt_nd_pos
 
double mu_escale_invsqrt_nd_pos
 
double n_escale_nd_pos
 
double n_escale_sqrt_nd_pos
 
double n_escale_invsqrt_nd_pos
 
double em_escale_nd_pos
 
double em_escale_sqrt_nd_pos
 
double em_escale_invsqrt_nd_pos
 
double had_res_nd_pos
 
double mu_res_nd_pos
 
double n_res_nd_pos
 
double em_res_nd_pos
 
std::vector< const double * > NDDetectorSystPointers
 ND Detector Systematics.
 
int nNDDetectorSystPointers
 Number of FD Detector Systematics.
 

Detailed Description

Base class for handling beam ND LAR samples.

Definition at line 23 of file samplePDFDUNEBeamND.h.

Member Enumeration Documentation

◆ KinematicTypes

Enum to identify kinematics.

Enumerator
kTrueNeutrinoEnergy 
kRecoQ 

Definition at line 36 of file samplePDFDUNEBeamND.h.

Constructor & Destructor Documentation

◆ samplePDFDUNEBeamND()

samplePDFDUNEBeamND::samplePDFDUNEBeamND ( std::string mc_version,
covarianceXsec * xsec_cov,
covarianceOsc * osc_cov )

samplePDFDUNE ND beam Constructor

Parameters
mc_versionConfig Name
xsec_covCross-section covariance matrix
osc_covOscillation covariance matrix

Definition at line 10 of file samplePDFDUNEBeamND.cpp.

10 : samplePDFFDBase(mc_version_, xsec_cov_, osc_cov_) {
11 Initialise();
12}

◆ ~samplePDFDUNEBeamND()

samplePDFDUNEBeamND::~samplePDFDUNEBeamND ( )

destructor

Definition at line 14 of file samplePDFDUNEBeamND.cpp.

14 {
15}

Member Function Documentation

◆ applyShifts()

void samplePDFDUNEBeamND::applyShifts ( int iSample,
int iEvent )
protected

NOT IMPLEMENTED: Apply kinematic shifts.

Parameters
iSampleSample Number
iEventEvent number

Definition at line 373 of file samplePDFDUNEBeamND.cpp.

373 {
374 // reset erec back to original value
375
376 dunendmcSamples[iSample].rw_erec_shifted[iEvent] = dunendmcSamples[iSample].rw_erec[iEvent];
377
378 //Calculate values needed
379 double sqrtErecHad = sqrt(dunendmcSamples[iSample].rw_erec_had[iEvent]);
380 double sqrtErecLep = sqrt(dunendmcSamples[iSample].rw_erec_lep[iEvent]);
381 double sqrteRecoPi0 = sqrt(dunendmcSamples[iSample].rw_eRecoPi0[iEvent]);
382 double sqrteRecoN = sqrt(dunendmcSamples[iSample].rw_eRecoN[iEvent]);
383 double sumEhad = dunendmcSamples[iSample].rw_eRecoP[iEvent] + dunendmcSamples[iSample].rw_eRecoPip[iEvent] + dunendmcSamples[iSample].rw_eRecoPim[iEvent];
384 double sqrtSumEhad = sqrt(sumEhad);
385
386 double invSqrtErecHad = 1/(sqrtErecHad+0.1);
387 double invSqrtErecLep = 1/(sqrtErecLep+0.1);
388 double invSqrteRecoPi0 = 1/(sqrteRecoPi0+0.1);
389 double invSqrteRecoN = 1/(sqrteRecoN+0.1);
390 double invSqrtSumEhad = 1/(sqrtSumEhad+0.1);
391
392 bool CCnumu {dunendmcSamples[iSample].rw_isCC[iEvent]==1 && abs(dunendmcSamples[iSample].rw_nuPDG[iEvent])==14 && dunendmcSamples[iSample].nupdgUnosc[iEvent]==2};
393 bool CCnue {dunendmcSamples[iSample].rw_isCC[iEvent]==1 && abs(dunendmcSamples[iSample].rw_nuPDG[iEvent])==12 && dunendmcSamples[iSample].nupdgUnosc[iEvent]==1};
394 bool NotCCnumu {!(dunendmcSamples[iSample].rw_isCC[iEvent]==1 && abs(dunendmcSamples[iSample].rw_nuPDG[iEvent])==14) && dunendmcSamples[iSample].nupdgUnosc[iEvent]==2};
395
396
397 TotalEScaleND(NDDetectorSystPointers[0], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_erec_had[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], NotCCnumu);
398
399 TotalEScaleSqrtND(NDDetectorSystPointers[1], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_erec_had[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], sqrtErecHad, sqrtErecLep, NotCCnumu);
400
401 TotalEScaleInvSqrtND(NDDetectorSystPointers[2], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_erec_had[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], invSqrtErecHad, invSqrtErecLep, NotCCnumu);
402
403 HadEScaleND(NDDetectorSystPointers[3], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], sumEhad);
404
405 HadEScaleSqrtND(NDDetectorSystPointers[4], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], sumEhad, sqrtSumEhad);
406
407 HadEScaleInvSqrtND(NDDetectorSystPointers[5], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], sumEhad, invSqrtSumEhad);
408
409 MuEScaleND(NDDetectorSystPointers[6], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], CCnumu);
410
411 MuEScaleSqrtND(NDDetectorSystPointers[7], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], sqrtErecLep, CCnumu);
412
413 MuEScaleInvSqrtND(NDDetectorSystPointers[8], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], invSqrtErecLep, CCnumu);
414
415 NEScaleND(NDDetectorSystPointers[9], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoN[iEvent]);
416
417 NEScaleSqrtND(NDDetectorSystPointers[10], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoN[iEvent], sqrteRecoN);
418
419 NEScaleInvSqrtND(NDDetectorSystPointers[11], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoN[iEvent], invSqrteRecoN);
420
421 EMEScaleND(NDDetectorSystPointers[12], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoPi0[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], CCnue);
422
423 EMEScaleSqrtND(NDDetectorSystPointers[13], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoPi0[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], sqrtErecLep, sqrteRecoPi0, CCnue);
424
425 EMEScaleInvSqrtND(NDDetectorSystPointers[14], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoPi0[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], invSqrtErecLep, invSqrteRecoPi0, CCnue);
426
427 HadResND(NDDetectorSystPointers[15], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoP[iEvent], dunendmcSamples[iSample].rw_eRecoPip[iEvent], dunendmcSamples[iSample].rw_eRecoPim[iEvent], dunendmcSamples[iSample].rw_eP[iEvent], dunendmcSamples[iSample].rw_ePip[iEvent], dunendmcSamples[iSample].rw_ePim[iEvent]);
428
429 MuResND(NDDetectorSystPointers[16], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], dunendmcSamples[iSample].rw_LepE[iEvent], CCnumu);
430
431 NResND(NDDetectorSystPointers[17], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoN[iEvent], dunendmcSamples[iSample].rw_eN[iEvent]);
432
433 EMResND(NDDetectorSystPointers[18], &dunendmcSamples[iSample].rw_erec_shifted[iEvent], dunendmcSamples[iSample].rw_eRecoPi0[iEvent], dunendmcSamples[iSample].rw_ePi0[iEvent], dunendmcSamples[iSample].rw_erec_lep[iEvent], dunendmcSamples[iSample].rw_LepE[iEvent], CCnue);
434
435}
void TotalEScaleInvSqrtND(const double *par, double *erec, double erecHad, double erecLep, double invSqrtErecHad, double invSqrtErecLep, bool NotCCnumu)
void EMEScaleSqrtND(const double *par, double *erec, double eRecoPi0, double erecLep, double sqrtErecLep, double sqrteRecoPi0, bool CCnue)
void MuEScaleSqrtND(const double *par, double *erec, double erecLep, double sqrtErecLep, bool CCnumu)
void EMResND(const double *par, double *erec, double eRecoPi0, double ePi0, double erecLep, double LepE, bool CCnue)
void NEScaleND(const double *par, double *erec, double eRecoN)
void NEScaleSqrtND(const double *par, double *erec, double eRecoN, double sqrteRecoN)
void EMEScaleND(const double *par, double *erec, double eRecoPi0, double erecLep, bool CCnue)
void TotalEScaleND(const double *par, double *erec, double erecHad, double erecLep, bool NotCCnumu)
void HadEScaleND(const double *par, double *erec, double sumEhad)
void NResND(const double *par, double *erec, double eRecoN, double eN)
void HadResND(const double *par, double *erec, double eRecoP, double eRecoPip, double eRecoPim, double eP, double ePip, double ePim)
void TotalEScaleSqrtND(const double *par, double *erec, double erecHad, double erecLep, double sqrtErecHad, double sqrtErecLep, bool NotCCnumu)
void EMEScaleInvSqrtND(const double *par, double *erec, double eRecoPi0, double erecLep, double invSqrtErecLep, double invSqrteRecoPi0, bool CCnue)
void HadEScaleSqrtND(const double *par, double *erec, double sumEhad, double sqrtSumEhad)
void HadEScaleInvSqrtND(const double *par, double *erec, double sumEhad, double invSqrtSumEhad)
void MuResND(const double *par, double *erec, double erecLep, double LepE, bool CCnumu)
void MuEScaleND(const double *par, double *erec, double erecLep, bool CCnumu)
void MuEScaleInvSqrtND(const double *par, double *erec, double erecLep, double invSqrtErecLep, bool CCnumu)
void NEScaleInvSqrtND(const double *par, double *erec, double eRecoN, double invSqrteRecoN)
std::vector< const double * > NDDetectorSystPointers
ND Detector Systematics.
std::vector< struct dunemc_base > dunendmcSamples
Array filled with MC samples for each oscillation channel.

◆ CalcXsecWeightFunc()

double samplePDFDUNEBeamND::CalcXsecWeightFunc ( int iSample,
int iEvent )
inlineprotected

NOT IMPLEMENTED: Dunder method to calculate xsec weights.

Parameters
iSamplesample ID
iEventEvent number

Definition at line 112 of file samplePDFDUNEBeamND.h.

112{return 1.;}

◆ GetPointerToKinematicParameter() [1/3]

const double * samplePDFDUNEBeamND::GetPointerToKinematicParameter ( double KinematicVariable,
int iSample,
int iEvent )
protected

Returns pointer to kinemtatic parameter for event in Structs DUNE.

Parameters
KinematicVariableKinematic parameter as double (gets cast -> int)
iSampleSample ID
iEventEvent ID
Returns
Pointer to KinPar for a given event

Definition at line 341 of file samplePDFDUNEBeamND.cpp.

341 {
342 KinematicTypes KinPar = (KinematicTypes) std::round(KinematicVariable);
343 return GetPointerToKinematicParameter(KinPar,iSample,iEvent);
344}
KinematicTypes
Enum to identify kinematics.
const double * GetPointerToKinematicParameter(KinematicTypes KinPar, int iSample, int iEvent)
Returns pointer to kinemtatic parameter for event in Structs DUNE.

◆ GetPointerToKinematicParameter() [2/3]

const double * samplePDFDUNEBeamND::GetPointerToKinematicParameter ( KinematicTypes KinPar,
int iSample,
int iEvent )
protected

Returns pointer to kinemtatic parameter for event in Structs DUNE.

Parameters
KinParKinematic parameter enum val
iSampleSample ID
iEventEvent ID
Returns
Pointer to KinPar for a given event

Definition at line 323 of file samplePDFDUNEBeamND.cpp.

323 {
324 double* KinematicValue;
325
326 switch(KinPar){
328 KinematicValue = &dunendmcSamples[iSample].rw_etru[iEvent];
329 break;
330 case kRecoQ:
331 KinematicValue = &dunendmcSamples[iSample].rw_reco_q[iEvent];
332 break;
333 default:
334 MACH3LOG_ERROR("Did not recognise Kinematic Parameter type...");
335 throw MaCh3Exception(__FILE__, __LINE__);
336 }
337
338 return KinematicValue;
339}

◆ GetPointerToKinematicParameter() [3/3]

const double * samplePDFDUNEBeamND::GetPointerToKinematicParameter ( std::string KinematicParameter,
int iSample,
int iEvent )
protected

Returns pointer to kinemtatic parameter for event in Structs DUNE.

Parameters
KinematicParameterKinematic parameter name as string (gets cast -> int)
iSampleSample ID
iEventEvent ID
Returns
Pointer to KinPar for a given event

Definition at line 346 of file samplePDFDUNEBeamND.cpp.

346 {
347 KinematicTypes KinPar = static_cast<KinematicTypes>(ReturnKinematicParameterFromString(KinematicParameter));
348 return GetPointerToKinematicParameter(KinPar,iSample,iEvent);
349}
int ReturnKinematicParameterFromString(std::string KinematicParameterStr)
Get kinematic parameter ID from string name.

◆ Init()

void samplePDFDUNEBeamND::Init ( )
protected

Initialises object.

Definition at line 17 of file samplePDFDUNEBeamND.cpp.

17 {
18 dunendmcSamples.resize(nSamples,dunemc_base());
19
20 IsRHC = SampleManager->raw()["SampleBools"]["isrhc"].as<bool>();
21 SampleDetID = SampleManager->raw()["DetID"].as<int>();
22 iselike = SampleManager->raw()["SampleBools"]["iselike"].as<bool>();
23
24 tot_escale_nd_pos = -999;
27 had_escale_nd_pos = -999;
30 mu_escale_nd_pos = -999;
33 n_escale_nd_pos = -999;
36 em_escale_nd_pos = -999;
39 had_res_nd_pos = -999;
40 mu_res_nd_pos = -999;
41 n_res_nd_pos = -999;
42 em_res_nd_pos = -999;
43
44 /*
45 nNDDetectorSystPointers = funcParsIndex.size();
46 NDDetectorSystPointers = std::vector<const double*>(nNDDetectorSystPointers);
47
48 int func_it = 0;
49 for (std::vector<int>::iterator it = funcParsIndex.begin(); it != funcParsIndex.end(); ++it, ++func_it) {
50 std::string name = funcParsNames.at(func_it);
51
52 if (name == "TotalEScaleND") {
53 tot_escale_nd_pos = *it;
54 NDDetectorSystPointers[func_it] = XsecCov->retPointer(tot_escale_nd_pos);
55 }
56 else if (name == "TotalEScaleSqrtND") {
57 tot_escale_sqrt_nd_pos = *it;
58 NDDetectorSystPointers[func_it] = XsecCov->retPointer(tot_escale_sqrt_nd_pos);
59 }
60 else if (name == "TotalEScaleInvSqrtND") {
61 tot_escale_invsqrt_nd_pos = *it;
62 NDDetectorSystPointers[func_it] = XsecCov->retPointer(tot_escale_invsqrt_nd_pos);
63 }
64 else if (name == "HadEScaleND") {
65 had_escale_nd_pos = *it;
66 NDDetectorSystPointers[func_it] = XsecCov->retPointer(had_escale_nd_pos);
67 }
68 else if (name == "HadEScaleSqrtND") {
69 had_escale_sqrt_nd_pos = *it;
70 NDDetectorSystPointers[func_it] = XsecCov->retPointer(had_escale_sqrt_nd_pos);
71 }
72 else if (name == "HadEScaleInvSqrtND") {
73 had_escale_invsqrt_nd_pos = *it;
74 NDDetectorSystPointers[func_it] = XsecCov->retPointer(had_escale_invsqrt_nd_pos);
75 }
76 else if (name == "MuEScaleND") {
77 mu_escale_nd_pos = *it;
78 NDDetectorSystPointers[func_it] = XsecCov->retPointer(mu_escale_nd_pos);
79 }
80 else if (name == "MuEScaleSqrtND") {
81 mu_escale_sqrt_nd_pos = *it;
82 NDDetectorSystPointers[func_it] = XsecCov->retPointer(mu_escale_sqrt_nd_pos);
83 }
84 else if (name == "MuEScaleInvSqrtND") {
85 mu_escale_invsqrt_nd_pos = *it;
86 NDDetectorSystPointers[func_it] = XsecCov->retPointer(mu_escale_invsqrt_nd_pos);
87 }
88 else if (name == "NEScaleND") {
89 n_escale_nd_pos = *it;
90 NDDetectorSystPointers[func_it] = XsecCov->retPointer(n_escale_nd_pos);
91 }
92 else if (name == "NEScaleSqrtND") {
93 n_escale_sqrt_nd_pos = *it;
94 NDDetectorSystPointers[func_it] = XsecCov->retPointer(n_escale_sqrt_nd_pos);
95 }
96 else if (name == "NEScaleInvSqrtND") {
97 n_escale_invsqrt_nd_pos = *it;
98 NDDetectorSystPointers[func_it] = XsecCov->retPointer(n_escale_invsqrt_nd_pos);
99 }
100 else if (name == "EMEScaleND") {
101 em_escale_nd_pos = *it;
102 NDDetectorSystPointers[func_it] = XsecCov->retPointer(em_escale_nd_pos);
103 }
104 else if (name == "EMEScaleSqrtND") {
105 em_escale_sqrt_nd_pos = *it;
106 NDDetectorSystPointers[func_it] = XsecCov->retPointer(em_escale_sqrt_nd_pos);
107 }
108 else if (name == "EMEScaleInvSqrtND") {
109 em_escale_invsqrt_nd_pos = *it;
110 NDDetectorSystPointers[func_it] = XsecCov->retPointer(em_escale_invsqrt_nd_pos);
111 }
112 else if (name == "HadResND") {
113 had_res_nd_pos = *it;
114 NDDetectorSystPointers[func_it] = XsecCov->retPointer(had_res_nd_pos);
115 }
116 else if (name == "MuResND") {
117 mu_res_nd_pos = *it;
118 NDDetectorSystPointers[func_it] = XsecCov->retPointer(mu_res_nd_pos);
119 }
120 else if (name == "NResND") {
121 n_res_nd_pos = *it;
122 NDDetectorSystPointers[func_it] = XsecCov->retPointer(n_res_nd_pos);
123 }
124 else if (name == "EMResND") {
125 em_res_nd_pos = *it;
126 NDDetectorSystPointers[func_it] = XsecCov->retPointer(em_res_nd_pos);
127 }
128 else {
129 MACH3LOG_ERROR("Found a functional parameter which wasn't specified in the xml | samplePDFDUNEBeamND: {}",name);
130 throw MaCh3Exception(__FILE__, __LINE__);
131 }
132 }
133 */
134
135 std::cout << "-------------------------------------------------------------------" <<std::endl;
136}

◆ ReturnKinematicParameter() [1/2]

double samplePDFDUNEBeamND::ReturnKinematicParameter ( double KinematicVariable,
int iSample,
int iEvent )
protected

Returns pointer to kinemtatic parameter for event in Structs DUNE.

Parameters
KinematicVariableKinematic parameter ID as double (gets cast -> int)
iSampleSample ID
iEventEvent ID
Returns
Value of kinetmatic parameter corresponding for a given event

Definition at line 351 of file samplePDFDUNEBeamND.cpp.

351 {
352 return *GetPointerToKinematicParameter(KinematicVariable, iSample, iEvent);
353}

◆ ReturnKinematicParameter() [2/2]

double samplePDFDUNEBeamND::ReturnKinematicParameter ( std::string KinematicParameter,
int iSample,
int iEvent )
protected

Returns pointer to kinemtatic parameter for event in Structs DUNE.

Parameters
KinematicParameterKinematic parameter name as string (gets cast -> int)
iSampleSample ID
iEventEvent ID
Returns
Value of kinematic parameter corresponding for a given event

Definition at line 355 of file samplePDFDUNEBeamND.cpp.

355 {
356 return *GetPointerToKinematicParameter(KinematicParameter, iSample, iEvent);
357}

◆ ReturnKinematicParameterBinning()

std::vector< double > samplePDFDUNEBeamND::ReturnKinematicParameterBinning ( std::string KinematicParameter)
protected

Gets binning for a given parameter.

Parameters
KinematicParameterStrParameter name
Returns
Vector containing parameter bins

Definition at line 437 of file samplePDFDUNEBeamND.cpp.

438{
439 std::vector<double> binningVector;
440 return binningVector;
441}

◆ ReturnKinematicParameterFromString()

int samplePDFDUNEBeamND::ReturnKinematicParameterFromString ( std::string KinematicParameterStr)
protected

Get kinematic parameter ID from string name.

Parameters
KinematicStr
Returns
Parameter ID

Definition at line 443 of file samplePDFDUNEBeamND.cpp.

443 {
444 return -1;
445}

◆ ReturnStringFromKinematicParameter()

std::string samplePDFDUNEBeamND::ReturnStringFromKinematicParameter ( int KinematicParameter)
protected

Gets name of kinematic parmaeter.

Parameters
KinParParameter ID
Returns
Name of parameter

Definition at line 447 of file samplePDFDUNEBeamND.cpp.

447 {
448 return "";
449}

◆ setupExperimentMC()

int samplePDFDUNEBeamND::setupExperimentMC ( int iSample)
protected

Function to setup MC from file.

Parameters
iSamplesample ID
Returns
Total number of events

Definition at line 168 of file samplePDFDUNEBeamND.cpp.

168 {
169
170 dunemc_base *duneobj = &(dunendmcSamples[iSample]);
171 int nupdgUnosc = sample_nupdgunosc[iSample];
172 int nupdg = sample_nupdg[iSample];
173
174 MACH3LOG_INFO("-------------------------------------------------------------------");
175 MACH3LOG_INFO("input file: {}", mc_files.at(iSample));
176
177 _sampleFile = TFile::Open(mc_files.at(iSample).c_str(), "READ");
178 _data = (TTree*)_sampleFile->Get("caf");
179
180 if(_data){
181 MACH3LOG_INFO("Found \"caf\" tree in {}", mc_files[iSample]);
182 MACH3LOG_INFO("With number of entries: {}", _data->GetEntries());
183 }
184 else{
185 MACH3LOG_ERROR("Could not find \"caf\" tree in {}", mc_files[iSample]);
186 throw MaCh3Exception(__FILE__, __LINE__);
187 }
188
189 _data->SetBranchStatus("*", 0);
190 _data->SetBranchStatus("Ev", 1);
191 _data->SetBranchAddress("Ev", &_ev);
192 _data->SetBranchStatus("Ev_reco", 1);
193 _data->SetBranchAddress("Ev_reco", &_erec);
194 _data->SetBranchStatus("Elep_reco", 1);
195 _data->SetBranchAddress("Elep_reco", &_erec_lep);
196 _data->SetBranchStatus("mode",1);
197 _data->SetBranchAddress("mode",&_mode);
198 _data->SetBranchStatus("isCC", 1);
199 _data->SetBranchAddress("isCC", &_isCC);
200 _data->SetBranchStatus("reco_q", 1);
201 _data->SetBranchAddress("reco_q", &_reco_q);
202 _data->SetBranchStatus("nuPDGunosc", 1);
203 _data->SetBranchAddress("nuPDGunosc", &_nuPDGunosc);
204 _data->SetBranchStatus("nuPDG", 1);
205 _data->SetBranchAddress("nuPDG", &_nuPDG);
206 _data->SetBranchStatus("BeRPA_A_cvwgt", 1);
207 _data->SetBranchAddress("BeRPA_A_cvwgt", &_BeRPA_cvwgt);
208
209 _data->SetBranchStatus("eRecoP", 1);
210 _data->SetBranchAddress("eRecoP", &_eRecoP);
211 _data->SetBranchStatus("eRecoPip", 1);
212 _data->SetBranchAddress("eRecoPip", &_eRecoPip);
213 _data->SetBranchStatus("eRecoPim", 1);
214 _data->SetBranchAddress("eRecoPim", &_eRecoPim);
215 _data->SetBranchStatus("eRecoPi0", 1);
216 _data->SetBranchAddress("eRecoPi0", &_eRecoPi0);
217 _data->SetBranchStatus("eRecoN", 1);
218 _data->SetBranchAddress("eRecoN", &_eRecoN);
219
220 _data->SetBranchStatus("LepE", 1);
221 _data->SetBranchAddress("LepE", &_LepE);
222 _data->SetBranchStatus("eP", 1);
223 _data->SetBranchAddress("eP", &_eP);
224 _data->SetBranchStatus("ePip", 1);
225 _data->SetBranchAddress("ePip", &_ePip);
226 _data->SetBranchStatus("ePim", 1);
227 _data->SetBranchAddress("ePim", &_ePim);
228 _data->SetBranchStatus("ePi0", 1);
229 _data->SetBranchAddress("ePi0", &_ePi0);
230 _data->SetBranchStatus("eN", 1);
231 _data->SetBranchAddress("eN", &_eN);
232
233 if (!IsRHC) {
234 duneobj->norm_s = (1e21/1.5e21);
235 } else {
236 duneobj->norm_s = (1e21/1.905e21);
237 }
238 duneobj->pot_s = (pot)/1e21;
239
240 duneobj->nEvents = _data->GetEntries();
241
242 duneobj->rw_yrec = new double[duneobj->nEvents];
243 duneobj->rw_erec_lep = new double[duneobj->nEvents];
244 duneobj->rw_erec_had = new double[duneobj->nEvents];
245 duneobj->rw_etru = new double[duneobj->nEvents];
246 duneobj->rw_erec = new double[duneobj->nEvents];
247 duneobj->rw_erec_shifted = new double[duneobj->nEvents];
248 duneobj->rw_theta = new double[duneobj->nEvents];
249 duneobj->flux_w = new double[duneobj->nEvents];
250 duneobj->rw_isCC = new int[duneobj->nEvents];
251 duneobj->rw_reco_q = new double[duneobj->nEvents];
252 duneobj->rw_nuPDGunosc = new int[duneobj->nEvents];
253 duneobj->rw_nuPDG = new int[duneobj->nEvents];
254 duneobj->rw_berpaacvwgt = new double[duneobj->nEvents];
255
256 duneobj->rw_eRecoP = new double[duneobj->nEvents];
257 duneobj->rw_eRecoPip = new double[duneobj->nEvents];
258 duneobj->rw_eRecoPim = new double[duneobj->nEvents];
259 duneobj->rw_eRecoPi0 = new double[duneobj->nEvents];
260 duneobj->rw_eRecoN = new double[duneobj->nEvents];
261
262 duneobj->rw_LepE = new double[duneobj->nEvents];
263 duneobj->rw_eP = new double[duneobj->nEvents];
264 duneobj->rw_ePip = new double[duneobj->nEvents];
265 duneobj->rw_ePim = new double[duneobj->nEvents];
266 duneobj->rw_ePi0 = new double[duneobj->nEvents];
267 duneobj->rw_eN = new double[duneobj->nEvents];
268
269 duneobj->nupdg = new int[duneobj->nEvents];
270 duneobj->nupdgUnosc = new int[duneobj->nEvents];
271 duneobj->mode = new double[duneobj->nEvents];
272 duneobj->Target = new int[duneobj->nEvents];
273
274 _data->GetEntry(0);
275
276 //FILL DUNE STRUCT
277 for (int i = 0; i < duneobj->nEvents; ++i) { // Loop through tree
278 _data->GetEntry(i);
279
280 duneobj->nupdg[i] = sample_nupdg[iSample];
281 duneobj->nupdgUnosc[i] = sample_nupdgunosc[iSample];
282
283 duneobj->rw_erec[i] = (double)_erec;
284 duneobj->rw_erec_shifted[i] = (double)_erec;
285 duneobj->rw_erec_lep[i] = (double)_erec_lep;
286 duneobj->rw_erec_had[i] = (double)(_erec - _erec_lep);
287 duneobj->rw_yrec[i] = (double)((_erec - _erec_lep)/_erec);
288 duneobj->rw_etru[i] = (double)_ev; // in GeV
289 duneobj->rw_theta[i] = (double)_LepNuAngle;
290 duneobj->rw_isCC[i] = _isCC;
291 duneobj->rw_reco_q[i] = _reco_q;
292 duneobj->rw_nuPDGunosc[i] = _nuPDGunosc;
293 duneobj->rw_nuPDG[i] = _nuPDG;
294 duneobj->rw_berpaacvwgt[i] = (double)_BeRPA_cvwgt;
295
296 duneobj->rw_eRecoP[i] = (double)_eRecoP;
297 duneobj->rw_eRecoPip[i] = (double)_eRecoPip;
298 duneobj->rw_eRecoPim[i] = (double)_eRecoPim;
299 duneobj->rw_eRecoPi0[i] = (double)_eRecoPi0;
300 duneobj->rw_eRecoN[i] = (double)_eRecoN;
301
302 duneobj->rw_LepE[i] = (double)_LepE;
303 duneobj->rw_eP[i] = (double)_eP;
304 duneobj->rw_ePip[i] = (double)_ePip;
305 duneobj->rw_ePim[i] = (double)_ePim;
306 duneobj->rw_ePi0[i] = (double)_ePi0;
307 duneobj->rw_eN[i] = (double)_eN;
308
309 //Assume everything is on Argon for now....
310 duneobj->Target[i] = 40;
311
312 int mode= TMath::Abs(_mode);
313 duneobj->mode[i]=(double)GENIEMode_ToMaCh3Mode(mode, _isCC);
314
315 duneobj->flux_w[i] = 1.0;
316 }
317
318 _sampleFile->Close();
319 return duneobj->nEvents;
320}
int GENIEMode_ToMaCh3Mode(int GENIE_mode, int isCC)
TTree * _data
TTree containing sample Data.
TFile * _sampleFile
File containing sample objects.
double pot
Value of POT used for sample.
double * rw_eRecoP
Definition StructsDUNE.h:22
double * rw_erec_had
Definition StructsDUNE.h:19
double * rw_eRecoPim
Definition StructsDUNE.h:24
int * nupdgUnosc
Definition StructsDUNE.h:15
double * rw_erec_lep
Definition StructsDUNE.h:20
double * mode
Definition StructsDUNE.h:65
double * rw_berpaacvwgt
Definition StructsDUNE.h:46
double * rw_yrec
Definition StructsDUNE.h:21
double * rw_erec_shifted
Definition StructsDUNE.h:18
double * rw_ePip
Definition StructsDUNE.h:30
double pot_s
Definition StructsDUNE.h:59
double * rw_erec
Definition StructsDUNE.h:17
int * rw_isCC
Definition StructsDUNE.h:47
double * rw_eP
Definition StructsDUNE.h:29
double * rw_ePim
Definition StructsDUNE.h:31
double * flux_w
Definition StructsDUNE.h:63
double * rw_eRecoN
Definition StructsDUNE.h:26
double * rw_eN
Definition StructsDUNE.h:33
double * rw_etru
Definition StructsDUNE.h:35
int * rw_nuPDGunosc
Definition StructsDUNE.h:48
double * rw_theta
Definition StructsDUNE.h:37
double * rw_eRecoPi0
Definition StructsDUNE.h:25
double norm_s
Definition StructsDUNE.h:60
double * rw_eRecoPip
Definition StructsDUNE.h:23
double * rw_LepE
Definition StructsDUNE.h:28
double * rw_reco_q
Definition StructsDUNE.h:56
int * rw_nuPDG
Definition StructsDUNE.h:49
double * rw_ePi0
Definition StructsDUNE.h:32

◆ setupFDMC()

void samplePDFDUNEBeamND::setupFDMC ( int iSample)
protected

Tells FD base which variables to point to/be set to.

Parameters
iSampleSample ID

Definition at line 359 of file samplePDFDUNEBeamND.cpp.

359 {
360 dunemc_base *duneobj = &(dunendmcSamples[iSample]);
361 FarDetectorCoreInfo *fdobj = &(MCSamples[iSample]);
362
363 for(int iEvent = 0 ;iEvent < fdobj->nEvents ; ++iEvent){
364 fdobj->rw_etru[iEvent] = &(duneobj->rw_etru[iEvent]);
365 fdobj->mode[iEvent] = &(duneobj->mode[iEvent]);
366 fdobj->Target[iEvent] = &(duneobj->Target[iEvent]);
367 fdobj->isNC[iEvent] = !(duneobj->rw_isCC[iEvent]);
368 fdobj->nupdgUnosc[iEvent] = &(duneobj->nupdgUnosc[iEvent]);
369 fdobj->nupdg[iEvent] = &(duneobj->nupdg[iEvent]);
370 }
371}

◆ SetupSplines()

void samplePDFDUNEBeamND::SetupSplines ( )
protected

Sets up splines.

Todo
move all of the spline setup into core

Definition at line 138 of file samplePDFDUNEBeamND.cpp.

138 {
140 if(XsecCov->GetNumParamsFromDetID(SampleDetID, kSpline) > 0){
141 MACH3LOG_INFO("Found {} splines for this sample so I will create a spline object", XsecCov->GetNumParamsFromDetID(SampleDetID, kSpline));
142 SplineHandler = std::unique_ptr<splineFDBase>(new splinesDUNE(XsecCov));
143 InitialiseSplineObject();
144 }
145 else{
146 MACH3LOG_INFO("Found {} splines for this sample so I will not load or evaluate splines", XsecCov->GetNumParamsFromDetID(SampleDetID, kSpline));
147 SplineHandler = nullptr;
148 }
149
150 return;
151}
Specialisation of FD (binned) spline class.
Definition splinesDUNE.h:8

◆ SetupWeightPointers()

void samplePDFDUNEBeamND::SetupWeightPointers ( )
protected

Sets up pointers weights for each event (oscillation/xsec/etc.)

Definition at line 153 of file samplePDFDUNEBeamND.cpp.

153 {
154 for (int i = 0; i < (int)dunendmcSamples.size(); ++i) {
155 for (int j = 0; j < dunendmcSamples[i].nEvents; ++j) {
156 MCSamples[i].ntotal_weight_pointers[j] = 6;
157 MCSamples[i].total_weight_pointers[j].resize(MCSamples[i].ntotal_weight_pointers[j]);
158 MCSamples[i].total_weight_pointers[j][0] = &(dunendmcSamples[i].pot_s);
159 MCSamples[i].total_weight_pointers[j][1] = &(dunendmcSamples[i].norm_s);
160 MCSamples[i].total_weight_pointers[j][2] = MCSamples[i].osc_w_pointer[j];
161 MCSamples[i].total_weight_pointers[j][3] = &(dunendmcSamples[i].rw_berpaacvwgt[j]);
162 MCSamples[i].total_weight_pointers[j][4] = &(dunendmcSamples[i].flux_w[j]);
163 MCSamples[i].total_weight_pointers[j][5] = &(MCSamples[i].xsec_w[j]);
164 }
165 }
166}

Member Data Documentation

◆ _BeRPA_cvwgt

double samplePDFDUNEBeamND::_BeRPA_cvwgt
protected

Definition at line 156 of file samplePDFDUNEBeamND.h.

◆ _data

TTree* samplePDFDUNEBeamND::_data
protected

TTree containing sample Data.

Definition at line 126 of file samplePDFDUNEBeamND.h.

◆ _eN

double samplePDFDUNEBeamND::_eN
protected

Definition at line 155 of file samplePDFDUNEBeamND.h.

◆ _eP

double samplePDFDUNEBeamND::_eP
protected

Definition at line 151 of file samplePDFDUNEBeamND.h.

◆ _ePi0

double samplePDFDUNEBeamND::_ePi0
protected

Definition at line 154 of file samplePDFDUNEBeamND.h.

◆ _ePim

double samplePDFDUNEBeamND::_ePim
protected

Definition at line 153 of file samplePDFDUNEBeamND.h.

◆ _ePip

double samplePDFDUNEBeamND::_ePip
protected

Definition at line 152 of file samplePDFDUNEBeamND.h.

◆ _erec

double samplePDFDUNEBeamND::_erec
protected

Definition at line 137 of file samplePDFDUNEBeamND.h.

◆ _erec_had

double samplePDFDUNEBeamND::_erec_had
protected

Definition at line 139 of file samplePDFDUNEBeamND.h.

◆ _erec_lep

double samplePDFDUNEBeamND::_erec_lep
protected

Definition at line 138 of file samplePDFDUNEBeamND.h.

◆ _eRecoN

double samplePDFDUNEBeamND::_eRecoN
protected

Definition at line 147 of file samplePDFDUNEBeamND.h.

◆ _eRecoP

double samplePDFDUNEBeamND::_eRecoP
protected

Definition at line 143 of file samplePDFDUNEBeamND.h.

◆ _eRecoPi0

double samplePDFDUNEBeamND::_eRecoPi0
protected

Definition at line 146 of file samplePDFDUNEBeamND.h.

◆ _eRecoPim

double samplePDFDUNEBeamND::_eRecoPim
protected

Definition at line 145 of file samplePDFDUNEBeamND.h.

◆ _eRecoPip

double samplePDFDUNEBeamND::_eRecoPip
protected

Definition at line 144 of file samplePDFDUNEBeamND.h.

◆ _ev

double samplePDFDUNEBeamND::_ev
protected

Definition at line 136 of file samplePDFDUNEBeamND.h.

◆ _isCC

int samplePDFDUNEBeamND::_isCC
protected

Definition at line 157 of file samplePDFDUNEBeamND.h.

◆ _isFHC

int samplePDFDUNEBeamND::_isFHC
protected

Definition at line 162 of file samplePDFDUNEBeamND.h.

◆ _isND

int samplePDFDUNEBeamND::_isND
protected

Definition at line 161 of file samplePDFDUNEBeamND.h.

◆ _LepE

double samplePDFDUNEBeamND::_LepE
protected

Definition at line 150 of file samplePDFDUNEBeamND.h.

◆ _LepNuAngle

double samplePDFDUNEBeamND::_LepNuAngle
protected

Definition at line 149 of file samplePDFDUNEBeamND.h.

◆ _LepTheta

double samplePDFDUNEBeamND::_LepTheta
protected

Definition at line 166 of file samplePDFDUNEBeamND.h.

◆ _mode

int samplePDFDUNEBeamND::_mode
protected

Definition at line 133 of file samplePDFDUNEBeamND.h.

◆ _nuPDG

int samplePDFDUNEBeamND::_nuPDG
protected

Definition at line 159 of file samplePDFDUNEBeamND.h.

◆ _nuPDGunosc

int samplePDFDUNEBeamND::_nuPDGunosc
protected

Definition at line 158 of file samplePDFDUNEBeamND.h.

◆ _nutype

TString samplePDFDUNEBeamND::_nutype
protected

Definition at line 131 of file samplePDFDUNEBeamND.h.

◆ _Q2

double samplePDFDUNEBeamND::_Q2
protected

Definition at line 167 of file samplePDFDUNEBeamND.h.

◆ _reco_nue

int samplePDFDUNEBeamND::_reco_nue
protected

Definition at line 141 of file samplePDFDUNEBeamND.h.

◆ _reco_numu

int samplePDFDUNEBeamND::_reco_numu
protected

Definition at line 140 of file samplePDFDUNEBeamND.h.

◆ _reco_q

int samplePDFDUNEBeamND::_reco_q
protected

Definition at line 168 of file samplePDFDUNEBeamND.h.

◆ _run

int samplePDFDUNEBeamND::_run
protected

Definition at line 160 of file samplePDFDUNEBeamND.h.

◆ _sampleFile

TFile* samplePDFDUNEBeamND::_sampleFile
protected

File containing sample objects.

Definition at line 123 of file samplePDFDUNEBeamND.h.

◆ _vtx_x

double samplePDFDUNEBeamND::_vtx_x
protected

Definition at line 163 of file samplePDFDUNEBeamND.h.

◆ _vtx_y

double samplePDFDUNEBeamND::_vtx_y
protected

Definition at line 164 of file samplePDFDUNEBeamND.h.

◆ _vtx_z

double samplePDFDUNEBeamND::_vtx_z
protected

Definition at line 165 of file samplePDFDUNEBeamND.h.

◆ dunendmcSamples

std::vector<struct dunemc_base> samplePDFDUNEBeamND::dunendmcSamples
protected

Array filled with MC samples for each oscillation channel.

Definition at line 120 of file samplePDFDUNEBeamND.h.

◆ em_escale_invsqrt_nd_pos

double samplePDFDUNEBeamND::em_escale_invsqrt_nd_pos
protected

Definition at line 190 of file samplePDFDUNEBeamND.h.

◆ em_escale_nd_pos

double samplePDFDUNEBeamND::em_escale_nd_pos
protected

Definition at line 188 of file samplePDFDUNEBeamND.h.

◆ em_escale_sqrt_nd_pos

double samplePDFDUNEBeamND::em_escale_sqrt_nd_pos
protected

Definition at line 189 of file samplePDFDUNEBeamND.h.

◆ em_res_nd_pos

double samplePDFDUNEBeamND::em_res_nd_pos
protected

Definition at line 194 of file samplePDFDUNEBeamND.h.

◆ had_escale_invsqrt_nd_pos

double samplePDFDUNEBeamND::had_escale_invsqrt_nd_pos
protected

Definition at line 181 of file samplePDFDUNEBeamND.h.

◆ had_escale_nd_pos

double samplePDFDUNEBeamND::had_escale_nd_pos
protected

Definition at line 179 of file samplePDFDUNEBeamND.h.

◆ had_escale_sqrt_nd_pos

double samplePDFDUNEBeamND::had_escale_sqrt_nd_pos
protected

Definition at line 180 of file samplePDFDUNEBeamND.h.

◆ had_res_nd_pos

double samplePDFDUNEBeamND::had_res_nd_pos
protected

Definition at line 191 of file samplePDFDUNEBeamND.h.

◆ iselike

bool samplePDFDUNEBeamND::iselike
protected

Definition at line 171 of file samplePDFDUNEBeamND.h.

◆ isND

bool samplePDFDUNEBeamND::isND
protected

Definition at line 172 of file samplePDFDUNEBeamND.h.

◆ IsRHC

bool samplePDFDUNEBeamND::IsRHC
protected

Definition at line 173 of file samplePDFDUNEBeamND.h.

◆ mu_escale_invsqrt_nd_pos

double samplePDFDUNEBeamND::mu_escale_invsqrt_nd_pos
protected

Definition at line 184 of file samplePDFDUNEBeamND.h.

◆ mu_escale_nd_pos

double samplePDFDUNEBeamND::mu_escale_nd_pos
protected

Definition at line 182 of file samplePDFDUNEBeamND.h.

◆ mu_escale_sqrt_nd_pos

double samplePDFDUNEBeamND::mu_escale_sqrt_nd_pos
protected

Definition at line 183 of file samplePDFDUNEBeamND.h.

◆ mu_res_nd_pos

double samplePDFDUNEBeamND::mu_res_nd_pos
protected

Definition at line 192 of file samplePDFDUNEBeamND.h.

◆ n_escale_invsqrt_nd_pos

double samplePDFDUNEBeamND::n_escale_invsqrt_nd_pos
protected

Definition at line 187 of file samplePDFDUNEBeamND.h.

◆ n_escale_nd_pos

double samplePDFDUNEBeamND::n_escale_nd_pos
protected

Definition at line 185 of file samplePDFDUNEBeamND.h.

◆ n_escale_sqrt_nd_pos

double samplePDFDUNEBeamND::n_escale_sqrt_nd_pos
protected

Definition at line 186 of file samplePDFDUNEBeamND.h.

◆ n_res_nd_pos

double samplePDFDUNEBeamND::n_res_nd_pos
protected

Definition at line 193 of file samplePDFDUNEBeamND.h.

◆ NDDetectorSystPointers

std::vector<const double*> samplePDFDUNEBeamND::NDDetectorSystPointers
protected

ND Detector Systematics.

Definition at line 197 of file samplePDFDUNEBeamND.h.

◆ nNDDetectorSystPointers

int samplePDFDUNEBeamND::nNDDetectorSystPointers
protected

Number of FD Detector Systematics.

Definition at line 200 of file samplePDFDUNEBeamND.h.

◆ pot

double samplePDFDUNEBeamND::pot
protected

Value of POT used for sample.

Definition at line 129 of file samplePDFDUNEBeamND.h.

◆ tot_escale_invsqrt_nd_pos

double samplePDFDUNEBeamND::tot_escale_invsqrt_nd_pos
protected

Definition at line 178 of file samplePDFDUNEBeamND.h.

◆ tot_escale_nd_pos

double samplePDFDUNEBeamND::tot_escale_nd_pos
protected

Definition at line 176 of file samplePDFDUNEBeamND.h.

◆ tot_escale_sqrt_nd_pos

double samplePDFDUNEBeamND::tot_escale_sqrt_nd_pos
protected

Definition at line 177 of file samplePDFDUNEBeamND.h.


The documentation for this class was generated from the following files: