MaCh3 DUNE 1.0.0
Reference Guide
All Classes Files Functions Variables Enumerations Enumerator Pages
samplePDFDUNEBeamNDGar Class Reference

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

#include <samplePDFDUNE/samplePDFDUNEBeamNDGar.h>

Inheritance diagram for samplePDFDUNEBeamNDGar:
Collaboration diagram for samplePDFDUNEBeamNDGar:

Public Types

enum  KinematicTypes {
  kTrueNeutrinoEnergy , kRecoNeutrinoEnergy , kTrueXPos , kTrueYPos ,
  kTrueZPos , kTrueRad , kNMuonsRecoOverTruth , kRecoLepEnergy ,
  kTrueLepEnergy , kRecoXPos , kRecoYPos , kRecoZPos ,
  kRecoRad , kLepPT , kLepPZ , kPionMultiplicity ,
  kNRecoParticles , kInFDV , kTrueMinusRecoEnergyRatio , kTrueMinusRecoEnergy ,
  kNTrueMuons , kNRecoMuons
}
 Enum to identify kinematics. More...
 

Public Member Functions

 samplePDFDUNEBeamNDGar (std::string mc_version, covarianceXsec *xsec_cov)
 Constructor.
 
 ~samplePDFDUNEBeamNDGar ()
 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.
 
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.
 
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::vector< double > ReturnKinematicParameterBinning (KinematicTypes KinematicParameter)
 Gets binning for a given parameter.
 
int ReturnKinematicParameterFromString (std::string KinematicParameterStr)
 Get kinematic parameter ID from string name.
 
std::string ReturnStringFromKinematicParameter (int KinematicParameter)
 Gets name of kinematic parmaeter.
 

Protected Attributes

std::vector< struct dunemc_basedunendgarmcSamples
 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_nue
 
double _elep_reco
 
double _LepNuAngle
 
int _reco_numu
 
int _reco_nue
 
double _BeRPA_cvwgt = 1
 
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
 
bool iscalo_reco
 
float muonscore_threshold
 
caf::StandardRecord * sr = new caf::StandardRecord()
 

Detailed Description

Base class for handling beam ND GAR samples.

Definition at line 24 of file samplePDFDUNEBeamNDGar.h.

Member Enumeration Documentation

◆ KinematicTypes

Enum to identify kinematics.

Enumerator
kTrueNeutrinoEnergy 
kRecoNeutrinoEnergy 
kTrueXPos 
kTrueYPos 
kTrueZPos 
kTrueRad 
kNMuonsRecoOverTruth 
kRecoLepEnergy 
kTrueLepEnergy 
kRecoXPos 
kRecoYPos 
kRecoZPos 
kRecoRad 
kLepPT 
kLepPZ 
kPionMultiplicity 
kNRecoParticles 
kInFDV 
kTrueMinusRecoEnergyRatio 
kTrueMinusRecoEnergy 
kNTrueMuons 
kNRecoMuons 

Definition at line 36 of file samplePDFDUNEBeamNDGar.h.

Constructor & Destructor Documentation

◆ samplePDFDUNEBeamNDGar()

samplePDFDUNEBeamNDGar::samplePDFDUNEBeamNDGar ( std::string mc_version,
covarianceXsec * xsec_cov )

Constructor.

Parameters
mc_versionConfiguration file
xsec_covcross-section covariance matrix

Definition at line 11 of file samplePDFDUNEBeamNDGar.cpp.

11 : samplePDFFDBase(mc_version_, XsecCov_) {
12 Initialise();
13}

◆ ~samplePDFDUNEBeamNDGar()

samplePDFDUNEBeamNDGar::~samplePDFDUNEBeamNDGar ( )

destructor

Definition at line 15 of file samplePDFDUNEBeamNDGar.cpp.

15 {
16}

Member Function Documentation

◆ applyShifts()

void samplePDFDUNEBeamNDGar::applyShifts ( int iSample,
int iEvent )
inlineprotected

NOT IMPLEMENTED: Apply kinematic shifts.

Parameters
iSampleSample Number
iEventEvent number

Definition at line 65 of file samplePDFDUNEBeamNDGar.h.

65{}

◆ CalcXsecWeightFunc()

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

NOT IMPLEMENTED: Dunder method to calculate xsec weights.

Parameters
iSamplesample ID
iEventEvent number

Definition at line 60 of file samplePDFDUNEBeamNDGar.h.

60{return 1.;}

◆ GetPointerToKinematicParameter() [1/3]

const double * samplePDFDUNEBeamNDGar::GetPointerToKinematicParameter ( 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 309 of file samplePDFDUNEBeamNDGar.cpp.

309 {
310 KinematicTypes KinPar = (KinematicTypes) std::round(KinematicVariable);
311 return GetPointerToKinematicParameter(KinPar,iSample,iEvent);
312}
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 * samplePDFDUNEBeamNDGar::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 252 of file samplePDFDUNEBeamNDGar.cpp.

252 {
253 double* KinematicValue;
254
255 switch(KinematicParameter) {
257 KinematicValue = &dunendgarmcSamples[iSample].rw_etru[iEvent];
258 break;
260 KinematicValue = &dunendgarmcSamples[iSample].rw_erec[iEvent];
261 break;
262 case kTrueXPos:
263 KinematicValue = &dunendgarmcSamples[iSample].rw_vtx_x[iEvent];
264 break;
265 case kTrueYPos:
266 KinematicValue = &dunendgarmcSamples[iSample].rw_vtx_y[iEvent];
267 break;
268 case kTrueZPos:
269 KinematicValue = &dunendgarmcSamples[iSample].rw_vtx_z[iEvent];
270 break;
271 case kTrueRad:
272 KinematicValue = &dunendgarmcSamples[iSample].rw_rad[iEvent];
273 break;
275 KinematicValue = &dunendgarmcSamples[iSample].nmuonsratio[iEvent];
276 break;
277 case kRecoLepEnergy:
278 KinematicValue = &dunendgarmcSamples[iSample].rw_elep_reco[iEvent];
279 break;
280 case kTrueLepEnergy:
281 KinematicValue = &dunendgarmcSamples[iSample].rw_elep_true[iEvent];
282 break;
283 case kRecoXPos:
284 KinematicValue = &dunendgarmcSamples[iSample].rw_reco_vtx_x[iEvent];
285 break;
286 case kRecoYPos:
287 KinematicValue = &dunendgarmcSamples[iSample].rw_reco_vtx_y[iEvent];
288 break;
289 case kRecoZPos:
290 KinematicValue = &dunendgarmcSamples[iSample].rw_reco_vtx_z[iEvent];
291 break;
292 case kRecoRad:
293 KinematicValue = &dunendgarmcSamples[iSample].rw_reco_rad[iEvent];
294 break;
295 case kLepPT:
296 KinematicValue = &dunendgarmcSamples[iSample].rw_lep_pT[iEvent];
297 break;
298 case kLepPZ:
299 KinematicValue = &dunendgarmcSamples[iSample].rw_lep_pZ[iEvent];
300 break;
301 default:
302 MACH3LOG_ERROR("Did not recognise Kinematic Parameter type...");
303 throw MaCh3Exception(__FILE__, __LINE__);
304 }
305
306 return KinematicValue;
307}
std::vector< struct dunemc_base > dunendgarmcSamples
Array filled with MC samples for each oscillation channel.

◆ GetPointerToKinematicParameter() [3/3]

const double * samplePDFDUNEBeamNDGar::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
Value of kinematic parameter corresponding for a given event

Definition at line 314 of file samplePDFDUNEBeamNDGar.cpp.

314 {
315 KinematicTypes KinPar = static_cast<KinematicTypes>(ReturnKinematicParameterFromString(KinematicParameter));
316 return GetPointerToKinematicParameter(KinPar,iSample,iEvent);
317}
int ReturnKinematicParameterFromString(std::string KinematicParameterStr)
Get kinematic parameter ID from string name.

◆ Init()

void samplePDFDUNEBeamNDGar::Init ( )
protected

Initialises object.

Definition at line 18 of file samplePDFDUNEBeamNDGar.cpp.

18 {
19 dunendgarmcSamples.resize(nSamples,dunemc_base());
20
21 iscalo_reco = SampleManager->raw()["SampleBools"]["iscalo_reco"].as<bool>(); //NK determine what reco used
22}

◆ ReturnKinematicParameter() [1/2]

double samplePDFDUNEBeamNDGar::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 319 of file samplePDFDUNEBeamNDGar.cpp.

319 {
320 return *GetPointerToKinematicParameter(KinematicVariable, iSample, iEvent);
321}

◆ ReturnKinematicParameter() [2/2]

double samplePDFDUNEBeamNDGar::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 323 of file samplePDFDUNEBeamNDGar.cpp.

323 {
324 return *GetPointerToKinematicParameter(KinematicParameter, iSample, iEvent);
325}

◆ ReturnKinematicParameterBinning() [1/2]

std::vector< double > samplePDFDUNEBeamNDGar::ReturnKinematicParameterBinning ( KinematicTypes KinematicParameter)
protected

Gets binning for a given parameter.

Parameters
KinParParameter ID
Returns
Vector containing parameter bins

Definition at line 350 of file samplePDFDUNEBeamNDGar.cpp.

350 {
351 std::vector<double> binningVector;
352 switch(KinPar){
354 for(double ibins =0; ibins<10*10; ibins++){
355 double binval = ibins/10;
356 binningVector.push_back(binval);
357 }
358 break;
360 for(double ibins =0; ibins<10*10; ibins++){
361 double binval = ibins/10;
362 binningVector.push_back(binval);
363 }
364 break;
365 case kRecoXPos:
366 case kTrueXPos:
367 for(double ibins =0; ibins<259*2; ibins++){
368 binningVector.push_back(ibins-259);
369 }
370 break;
371 case kRecoYPos:
372 case kTrueYPos:
373 for(double ibins =0; ibins<277*2; ibins++){
374 binningVector.push_back(ibins-277-150);
375 }
376 break;
377 case kRecoZPos:
378 case kTrueZPos:
379 for(double ibins =0; ibins<277*2; ibins++){
380 binningVector.push_back(ibins-277+1486);
381 }
382 break;
384 for(double ibins =0; ibins<10; ibins++){
385 binningVector.push_back(ibins);
386 }
387 break;
388 case kNRecoParticles:
389 for(double ibins =0; ibins<50; ibins++){
390 binningVector.push_back(ibins);
391 }
392 break;
393 case kInFDV:
394 for(double ibins =0; ibins<3; ibins++){
395 binningVector.push_back(ibins);
396 }
397 break;
400 for(double ibins =0; ibins<20*10; ibins++){
401 binningVector.push_back(-10+(double)(ibins)/10);
402 }
403 break;
405 for(double ibins =0; ibins<20*10; ibins++){
406 binningVector.push_back(-10+(double)(ibins)/10);
407 }
408 break;
409 case kNTrueMuons:
410 case kNRecoMuons:
411 for(double ibins =0; ibins<10; ibins++){
412 binningVector.push_back(ibins);
413 }
414 break;
415 case kRecoLepEnergy:
416 for(double ibins =0; ibins<10*10; ibins++){
417 binningVector.push_back((double)(ibins)/10);
418 }
419 break;
420 case kTrueLepEnergy:
421 for(double ibins =0; ibins<10*10; ibins++){
422 binningVector.push_back((double)(ibins)/10);
423 }
424 break;
425 case kTrueRad:
426 case kRecoRad:
427 for(double ibins =0; ibins<300; ibins++){
428 binningVector.push_back(ibins);
429 }
430 break;
431 case kLepPT:
432 case kLepPZ:
433 for(double ibins =0; ibins<10*10; ibins++){
434 binningVector.push_back((double)(ibins)/10);
435 }
436 break;
437 default:
438 for(double ibins =0; ibins<10*100; ibins++){
439 binningVector.push_back(ibins/100);
440 }
441 break;
442 }
443
444 return binningVector;
445}

◆ ReturnKinematicParameterBinning() [2/2]

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

Gets binning for a given parameter.

Parameters
KinematicParameterStrParameter name
Returns
Vector containing parameter bins

Definition at line 345 of file samplePDFDUNEBeamNDGar.cpp.

345 {
346 KinematicTypes KinPar = static_cast<KinematicTypes>(ReturnKinematicParameterFromString(KinematicParameterStr));
347 return ReturnKinematicParameterBinning(KinPar);
348}
std::vector< double > ReturnKinematicParameterBinning(std::string KinematicParameter)
Gets binning for a given parameter.

◆ ReturnKinematicParameterFromString()

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

Get kinematic parameter ID from string name.

Parameters
KinematicStr
Returns
Parameter ID

Definition at line 447 of file samplePDFDUNEBeamNDGar.cpp.

447 {
448 return -1;
449}

◆ ReturnStringFromKinematicParameter()

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

Gets name of kinematic parmaeter.

Parameters
KinParParameter ID
Returns
Name of parameter

Definition at line 451 of file samplePDFDUNEBeamNDGar.cpp.

451 {
452 return "";
453}

◆ setupExperimentMC()

int samplePDFDUNEBeamNDGar::setupExperimentMC ( int iSample)
protected

Function to setup MC from file.

Parameters
iSamplesample ID
Returns
Total number of events

Definition at line 56 of file samplePDFDUNEBeamNDGar.cpp.

56 {
57
58 dunemc_base *duneobj = &(dunendgarmcSamples[iSample]);
59 int nutype = sample_nutype[iSample];
60 int oscnutype = sample_oscnutype[iSample];
61 bool signal = sample_signal[iSample];
62
63 MACH3LOG_INFO("-------------------------------------------------------------------");
64 MACH3LOG_INFO("Input File: {}", mc_files.at(iSample));
65
66 _sampleFile = TFile::Open(mc_files.at(iSample).c_str(), "READ");
67 _data = (TTree*)_sampleFile->Get("cafTree");
68
69 if(_data){
70 MACH3LOG_INFO("Found \"caf\" tree in {}", mc_files[iSample]);
71 MACH3LOG_INFO("With number of entries: {}", _data->GetEntries());
72 }
73 else{
74 MACH3LOG_ERROR("Could not find \"caf\" tree in {}", mc_files[iSample]);
75 throw MaCh3Exception(__FILE__, __LINE__);
76 }
77
78 _data->SetBranchStatus("*", 1);
79 _data->SetBranchAddress("rec", &sr);
80
81 duneobj->norm_s = 1.0;
82 duneobj->pot_s = (pot)/1e21;
83
84 duneobj->nEvents = _data->GetEntries();
85 duneobj->nutype = nutype;
86 duneobj->oscnutype = oscnutype;
87 duneobj->signal = signal;
88
89 // allocate memory for dunendgarmc variables
90 duneobj->rw_yrec = new double[duneobj->nEvents];
91 duneobj->rw_elep_reco = new double[duneobj->nEvents];
92 duneobj->rw_etru = new double[duneobj->nEvents];
93 duneobj->rw_erec = new double[duneobj->nEvents];
94 duneobj->flux_w = new double[duneobj->nEvents];
95 duneobj->rw_isCC = new int[duneobj->nEvents];
96 duneobj->rw_nuPDGunosc = new int[duneobj->nEvents];
97 duneobj->rw_nuPDG = new int[duneobj->nEvents];
98 duneobj->rw_berpaacvwgt = new double[duneobj->nEvents];
99
100 duneobj->mode = new double[duneobj->nEvents];
101
102 duneobj->nproton = new int[duneobj->nEvents];
103 duneobj->nneutron = new int[duneobj->nEvents];
104 duneobj->npip = new int[duneobj->nEvents];
105 duneobj->npim = new int[duneobj->nEvents];
106 duneobj->npi0 = new int[duneobj->nEvents];
107
108 duneobj->nrecomuon = new int[duneobj->nEvents];
109 duneobj->ntruemuon = new int[duneobj->nEvents];
110 duneobj->nmuonsratio = new double[duneobj->nEvents];
111 duneobj->ntruemuonprim = new int[duneobj->nEvents];
112
113 duneobj->nrecoparticles = new int[duneobj->nEvents];
114 duneobj->in_fdv = new bool[duneobj->nEvents];
115 duneobj->rw_elep_true = new double[duneobj->nEvents];
116
117 duneobj->rw_vtx_x = new double[duneobj->nEvents];
118 duneobj->rw_vtx_y = new double[duneobj->nEvents];
119 duneobj->rw_vtx_z = new double[duneobj->nEvents];
120 duneobj->rw_rad = new double[duneobj->nEvents];
121
122 duneobj->rw_lep_pT = new double[duneobj->nEvents];
123 duneobj->rw_lep_pZ = new double[duneobj->nEvents];
124
125 duneobj->rw_reco_vtx_x = new double[duneobj->nEvents];
126 duneobj->rw_reco_vtx_y = new double[duneobj->nEvents];
127 duneobj->rw_reco_vtx_z = new double[duneobj->nEvents];
128 duneobj->rw_reco_rad = new double[duneobj->nEvents];
129
130 duneobj->Target = new int[duneobj->nEvents];
131
132 int num_no_ixns =0;
133 int num_no_recparticles = 0;
134 int num_in_fdv = 0;
135 int num_in_fdv_noreco = 0;
136 int num_notin_fdv =0;
137 int num_nanenergy =0;
138 int num_nanparticles =0;
139
140 //FILL DUNE STRUCT
141 for (int i = 0; i < (duneobj->nEvents); ++i) { // Loop through tree
142 _data->GetEntry(i);
143 double radius = pow((pow((sr->mc.nu[0].vtx.y+150),2) + pow((sr->mc.nu[0].vtx.z-1486),2)),0.5);
144 if(std::abs(sr->mc.nu[0].vtx.x)<=209.0 && radius<=227.02){
145 num_in_fdv++;
146 duneobj->in_fdv[i] = 1;
147 } else{
148 num_notin_fdv++;
149 duneobj->in_fdv[i] = 0;
150 }
151
152 if(sr->common.ixn.ngsft == 0){
153 //duneobj->rw_erec[i] = (double)(0);
154 float erec_total =0;
155 duneobj->rw_elep_reco[i] = double(0);
156 duneobj->rw_yrec[i] = (double)(0);
157 num_no_ixns++;
158 duneobj->nrecoparticles[i] = (int)(0);
159 } else{
160 duneobj->nrecoparticles[i] = (int)(0);
161 float erec_total =0;
162 float elep_reco =0;
163 float muonscore = muonscore_threshold;
164 int nixns = (int)(sr->common.ixn.ngsft);
165 for(int i_ixn =0; i_ixn<nixns; i_ixn++) {
166 int nrecoparticles = (int)(sr->common.ixn.gsft[i_ixn].part.ngsft);
167 duneobj->nrecoparticles[i] += (int)(sr->common.ixn.gsft[i_ixn].part.ngsft);
168 int nanparticles = 0;
169 if(nrecoparticles ==0){
170 double radius = pow((pow((sr->mc.nu[0].vtx.y+150.),2) + pow((sr->mc.nu[0].vtx.z-1486.),2)),0.5);
171 if(std::abs(sr->mc.nu[0].vtx.x)<=209.0 || radius<=227.02) {
172 num_in_fdv_noreco++;
173 }
174 num_no_recparticles++;}
175 for(int i_part =0; i_part<nrecoparticles; i_part++) {
176 float erec_part = (float)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].E);
177 if(std::isnan(erec_part)){nanparticles++;}
178 erec_total+=erec_part;
179 if((float)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].score.gsft_pid.muon_score>muonscore)){
180 if(erec_part>elep_reco){
181 elep_reco = erec_part;
182 duneobj->rw_reco_vtx_x[i] = (double)((float)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].start.x));
183 duneobj->rw_reco_vtx_y[i] = (double)((float)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].start.y));
184 duneobj->rw_reco_vtx_z[i] = (double)((float)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].start.z));
185 duneobj->rw_lep_pT[i] = (double)(pow(pow((float)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].p.x), 2) + pow((float)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].p.y), 2), 0.5));
186 duneobj->rw_lep_pZ[i] = (double)(sr->common.ixn.gsft[i_ixn].part.gsft[i_part].p.z);
187 }
188 duneobj->nrecomuon[i]++;
189 }
190 }
191 num_nanparticles = num_nanparticles + (nanparticles/nrecoparticles);
192 } //ADD PRIMARY LEPTON ENERGY ELEP_RECO
193 if(std::isnan(erec_total)){num_nanenergy++; erec_total = (float)(sr->common.ixn.gsft[0].Enu.lep_calo);}
194 if(iscalo_reco){duneobj->rw_erec[i]=(double)(sr->common.ixn.gsft[0].Enu.lep_calo);}
195 else{duneobj->rw_erec[i]=(double)(erec_total);}
196 duneobj->rw_elep_reco[i] = (double)(elep_reco);
197 }
198
199 if(duneobj->rw_erec[i] != 0){duneobj->rw_yrec[i] = (double)(((duneobj->rw_erec[i])-(duneobj->rw_elep_reco[i]))/(duneobj->rw_erec[i]));}
200 else{duneobj->rw_yrec[i] = (double)(0);}
201 duneobj->rw_etru[i] = (double)(sr->mc.nu[0].E); // in GeV
202 duneobj->rw_isCC[i] = (int)(sr->mc.nu[0].iscc);
203 duneobj->rw_nuPDGunosc[i] = sr->mc.nu[0].pdgorig;
204 duneobj->rw_nuPDG[i] = sr->mc.nu[0].pdg;
205 duneobj->rw_berpaacvwgt[i] = _BeRPA_cvwgt;
206
207 int ntrueparticles = (int)(sr->mc.nu[0].nprim);
208 for(int i_truepart =0; i_truepart<ntrueparticles; i_truepart++){
209 if(std::abs(sr->mc.nu[0].prim[i_truepart].pdg) == 13){
210 duneobj->ntruemuon[i]++;
211 duneobj->ntruemuonprim[i]++;
212 }
213 }
214 int ntruesecparticles = (int)(sr->mc.nu[0].nsec);
215 for(int i_truepart =0; i_truepart<ntruesecparticles; i_truepart++){
216 if(std::abs(sr->mc.nu[0].sec[i_truepart].pdg) == 13){
217 duneobj->ntruemuon[i]++;
218 }
219 }
220
221 duneobj->nproton[i] = sr->mc.nu[0].nproton;
222 duneobj->nneutron[i] = sr->mc.nu[0].nneutron;
223 duneobj->npip[i] = sr->mc.nu[0].npip;
224 duneobj->npim[i] = sr->mc.nu[0].npim;
225 duneobj->npi0[i] = sr->mc.nu[0].npi0;
226
227 duneobj->nmuonsratio[i] = (double)(duneobj->nrecomuon[i])/(double)(duneobj->ntruemuonprim[i]);
228 duneobj->rw_vtx_x[i] = (double)(sr->mc.nu[0].vtx.x);
229 duneobj->rw_vtx_y[i] = (double)(sr->mc.nu[0].vtx.y);
230 duneobj->rw_vtx_z[i] = (double)(sr->mc.nu[0].vtx.z);
231
232 duneobj->rw_rad[i] = (double)(pow((pow((duneobj->rw_vtx_y[i]+150),2) + pow((duneobj->rw_vtx_z[i]-1486),2)),0.5));
233 duneobj->rw_reco_rad[i] = (double)(pow(pow((duneobj->rw_reco_vtx_y[i]+150),2) + pow((duneobj->rw_reco_vtx_z[i]-1486), 2), 0.5));
234 duneobj->rw_elep_true[i] = (double)(sr->mc.nu[0].prim[0].p.E);
235
236 //Assume everything is on Argon for now....
237 duneobj->Target[i] = 40;
238
239 _mode = sr->mc.nu[0].mode;
240 _isCC = (int)(sr->mc.nu[0].iscc);
241
242 int mode= TMath::Abs(_mode);
243 duneobj->mode[i]=(double)GENIEMode_ToMaCh3Mode(mode, _isCC);
244
245 duneobj->flux_w[i] = 1.0;
246 }
247
248 _sampleFile->Close();
249 return duneobj->nEvents;
250}
int GENIEMode_ToMaCh3Mode(int GENIE_mode, int isCC)
TFile * _sampleFile
File containing sample objects.
TTree * _data
TTree containing sample Data.
double pot
Value of POT used for sample.
double * rw_vtx_x
Definition StructsDUNE.h:52
double * rw_vtx_y
Definition StructsDUNE.h:53
double * rw_reco_rad
Definition StructsDUNE.h:86
double * mode
Definition StructsDUNE.h:65
double * rw_berpaacvwgt
Definition StructsDUNE.h:46
int * npip
number of (post-FSI) primary pi+
Definition StructsDUNE.h:72
double * rw_yrec
Definition StructsDUNE.h:21
int * ntruemuon
Definition StructsDUNE.h:76
double pot_s
Definition StructsDUNE.h:59
double * nmuonsratio
Definition StructsDUNE.h:79
double * rw_reco_vtx_z
Definition StructsDUNE.h:85
int * nrecoparticles
Definition StructsDUNE.h:92
double * rw_erec
Definition StructsDUNE.h:17
int * rw_isCC
Definition StructsDUNE.h:47
int * ntruemuonprim
Definition StructsDUNE.h:77
double * flux_w
Definition StructsDUNE.h:63
int * npim
number of (post-FSI) primary pi-
Definition StructsDUNE.h:73
double * rw_etru
Definition StructsDUNE.h:35
int * rw_nuPDGunosc
Definition StructsDUNE.h:48
int * npi0
number of (post-FSI) primary pi0
Definition StructsDUNE.h:74
double * rw_elep_reco
Definition StructsDUNE.h:89
double * rw_reco_vtx_y
Definition StructsDUNE.h:84
int * nrecomuon
Definition StructsDUNE.h:78
double norm_s
Definition StructsDUNE.h:60
double * rw_reco_vtx_x
Definition StructsDUNE.h:83
int * nproton
number of (post-FSI) primary protons
Definition StructsDUNE.h:70
double * rw_lep_pT
Definition StructsDUNE.h:81
double * rw_lep_pZ
Definition StructsDUNE.h:82
double * rw_vtx_z
Definition StructsDUNE.h:54
int * nneutron
number of (post-FSI) primary neutrons
Definition StructsDUNE.h:71
double * rw_elep_true
Definition StructsDUNE.h:90
int * rw_nuPDG
Definition StructsDUNE.h:49
double * rw_rad
Definition StructsDUNE.h:87
bool * in_fdv
Definition StructsDUNE.h:93

◆ setupFDMC()

void samplePDFDUNEBeamNDGar::setupFDMC ( int iSample)
protected

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

Parameters
iSampleSample ID

Definition at line 327 of file samplePDFDUNEBeamNDGar.cpp.

327 {
328 dunemc_base *duneobj = &(dunendgarmcSamples[iSample]);
329 fdmc_base *fdobj = &(MCSamples[iSample]);
330
331 fdobj->nutype = duneobj->nutype;
332 fdobj->oscnutype = duneobj->oscnutype;
333 fdobj->signal = duneobj->signal;
334 fdobj->SampleDetID = SampleDetID;
335
336 for(int iEvent = 0 ;iEvent < fdobj->nEvents ; ++iEvent){
337 fdobj->rw_etru[iEvent] = &(duneobj->rw_etru[iEvent]);
338 fdobj->mode[iEvent] = &(duneobj->mode[iEvent]);
339 fdobj->Target[iEvent] = &(duneobj->Target[iEvent]);
340
341 }
342
343}

◆ SetupSplines()

void samplePDFDUNEBeamNDGar::SetupSplines ( )
protected

Sets up splines.

Todo
move all of the spline setup into core

Definition at line 24 of file samplePDFDUNEBeamNDGar.cpp.

24 {
26 if(XsecCov->GetNumParamsFromDetID(SampleDetID, kSpline) > 0){
27 MACH3LOG_INFO("Found {} splines for this sample so I will create a spline object", XsecCov->GetNumParamsFromDetID(SampleDetID, kSpline));
28 splinesDUNE* DUNESplines = new splinesDUNE(XsecCov);
29 splineFile = (splineFDBase*)DUNESplines;
30 InitialiseSplineObject();
31 }
32 else{
33 MACH3LOG_INFO("Found {} splines for this sample so I will not load or evaluate splines", XsecCov->GetNumParamsFromDetID(SampleDetID, kSpline));
34 splineFile = nullptr;
35 }
36
37 return;
38}
Specialisation of FD (binned) spline class.
Definition splinesDUNE.h:8

◆ SetupWeightPointers()

void samplePDFDUNEBeamNDGar::SetupWeightPointers ( )
protected

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

Definition at line 41 of file samplePDFDUNEBeamNDGar.cpp.

41 {
42 for (int i = 0; i < (int)dunendgarmcSamples.size(); ++i) {
43 for (int j = 0; j < dunendgarmcSamples[i].nEvents; ++j) {
44 MCSamples[i].ntotal_weight_pointers[j] = 6;
45 MCSamples[i].total_weight_pointers[j] = new const double*[MCSamples[i].ntotal_weight_pointers[j]];
46 MCSamples[i].total_weight_pointers[j][0] = &(dunendgarmcSamples[i].pot_s);
47 MCSamples[i].total_weight_pointers[j][1] = &(dunendgarmcSamples[i].norm_s);
48 MCSamples[i].total_weight_pointers[j][2] = MCSamples[i].osc_w_pointer[j];
49 MCSamples[i].total_weight_pointers[j][3] = &(dunendgarmcSamples[i].rw_berpaacvwgt[j]);
50 MCSamples[i].total_weight_pointers[j][4] = &(dunendgarmcSamples[i].flux_w[j]);
51 MCSamples[i].total_weight_pointers[j][5] = &(MCSamples[i].xsec_w[j]);
52 }
53 }
54}

Member Data Documentation

◆ _BeRPA_cvwgt

double samplePDFDUNEBeamNDGar::_BeRPA_cvwgt = 1
protected

Definition at line 146 of file samplePDFDUNEBeamNDGar.h.

◆ _data

TTree* samplePDFDUNEBeamNDGar::_data
protected

TTree containing sample Data.

Definition at line 129 of file samplePDFDUNEBeamNDGar.h.

◆ _elep_reco

double samplePDFDUNEBeamNDGar::_elep_reco
protected

Definition at line 142 of file samplePDFDUNEBeamNDGar.h.

◆ _erec

double samplePDFDUNEBeamNDGar::_erec
protected

Definition at line 140 of file samplePDFDUNEBeamNDGar.h.

◆ _erec_nue

double samplePDFDUNEBeamNDGar::_erec_nue
protected

Definition at line 141 of file samplePDFDUNEBeamNDGar.h.

◆ _ev

double samplePDFDUNEBeamNDGar::_ev
protected

Definition at line 139 of file samplePDFDUNEBeamNDGar.h.

◆ _isCC

int samplePDFDUNEBeamNDGar::_isCC
protected

Definition at line 147 of file samplePDFDUNEBeamNDGar.h.

◆ _isFHC

int samplePDFDUNEBeamNDGar::_isFHC
protected

Definition at line 152 of file samplePDFDUNEBeamNDGar.h.

◆ _isND

int samplePDFDUNEBeamNDGar::_isND
protected

Definition at line 151 of file samplePDFDUNEBeamNDGar.h.

◆ _LepNuAngle

double samplePDFDUNEBeamNDGar::_LepNuAngle
protected

Definition at line 143 of file samplePDFDUNEBeamNDGar.h.

◆ _LepTheta

double samplePDFDUNEBeamNDGar::_LepTheta
protected

Definition at line 156 of file samplePDFDUNEBeamNDGar.h.

◆ _mode

int samplePDFDUNEBeamNDGar::_mode
protected

Definition at line 135 of file samplePDFDUNEBeamNDGar.h.

◆ _nuPDG

int samplePDFDUNEBeamNDGar::_nuPDG
protected

Definition at line 149 of file samplePDFDUNEBeamNDGar.h.

◆ _nuPDGunosc

int samplePDFDUNEBeamNDGar::_nuPDGunosc
protected

Definition at line 148 of file samplePDFDUNEBeamNDGar.h.

◆ _nutype

TString samplePDFDUNEBeamNDGar::_nutype
protected

Definition at line 134 of file samplePDFDUNEBeamNDGar.h.

◆ _Q2

double samplePDFDUNEBeamNDGar::_Q2
protected

Definition at line 157 of file samplePDFDUNEBeamNDGar.h.

◆ _reco_nue

int samplePDFDUNEBeamNDGar::_reco_nue
protected

Definition at line 145 of file samplePDFDUNEBeamNDGar.h.

◆ _reco_numu

int samplePDFDUNEBeamNDGar::_reco_numu
protected

Definition at line 144 of file samplePDFDUNEBeamNDGar.h.

◆ _run

int samplePDFDUNEBeamNDGar::_run
protected

Definition at line 150 of file samplePDFDUNEBeamNDGar.h.

◆ _sampleFile

TFile* samplePDFDUNEBeamNDGar::_sampleFile
protected

File containing sample objects.

Definition at line 127 of file samplePDFDUNEBeamNDGar.h.

◆ _vtx_x

double samplePDFDUNEBeamNDGar::_vtx_x
protected

Definition at line 153 of file samplePDFDUNEBeamNDGar.h.

◆ _vtx_y

double samplePDFDUNEBeamNDGar::_vtx_y
protected

Definition at line 154 of file samplePDFDUNEBeamNDGar.h.

◆ _vtx_z

double samplePDFDUNEBeamNDGar::_vtx_z
protected

Definition at line 155 of file samplePDFDUNEBeamNDGar.h.

◆ dunendgarmcSamples

std::vector<struct dunemc_base> samplePDFDUNEBeamNDGar::dunendgarmcSamples
protected

Array filled with MC samples for each oscillation channel.

Definition at line 124 of file samplePDFDUNEBeamNDGar.h.

◆ iscalo_reco

bool samplePDFDUNEBeamNDGar::iscalo_reco
protected

Definition at line 159 of file samplePDFDUNEBeamNDGar.h.

◆ muonscore_threshold

float samplePDFDUNEBeamNDGar::muonscore_threshold
protected

Definition at line 160 of file samplePDFDUNEBeamNDGar.h.

◆ pot

double samplePDFDUNEBeamNDGar::pot
protected

Value of POT used for sample.

Definition at line 131 of file samplePDFDUNEBeamNDGar.h.

◆ sr

caf::StandardRecord* samplePDFDUNEBeamNDGar::sr = new caf::StandardRecord()
protected

Definition at line 162 of file samplePDFDUNEBeamNDGar.h.


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