All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BooNEInterface.h
Go to the documentation of this file.
1 /**
2  * @file BooNEInterface.h
3  * @author Marco Del Tutto
4  * @author Zarko Pavlovic
5  * @date 3 Dec 2021
6  * @brief Interface to read flux files in BooNE tuple format.
7  */
8 
9 #ifndef _BOONEINTERFACE_H_
10 #define _BOONEINTERFACE_H_
11 
12 #include "FluxInterface.h"
13 
14 //
15 // boone ntuple structure
16 //
17 struct BooNENtuple {
18  float beamwgt; /// Magic Weight: CrossSect * BooBeamNT
19  int ntp; /// 1,2,3,4: nue, nuebar, numu, numubar
20  int npart; /// number of particles in the chain
21  /// npart-1 == proton
22  /// 0 == neutrino
23  int id[20]; /// id of each particle in before chain, array length 'npart'
24  float ini_pos[20][3]; /// 3-pos of particle in before chain, array length 'npart'
25  float ini_mom[20][3]; /// 3-mom of particle in before chain, array length 'npart'
26  float ini_eng[20]; /// E of particle in before chain, array length 'npart'
27  float ini_t[20]; /// "decay" time of particle (starting from proton)
28  /// in before chain, array length 'npart'
29  float fin_mom[20][3]; /// final 3-mom of particle in before chain, array length 'npart'
30  float fin_pol[20][3]; /// final 3-polarization of particle in before chain, array length 'npart'
31 
32 };
33 
34 //
35 // boone beam ntuple structure
36 //
37 struct BeamNtuple {
38  float tank_pos_beam[3]; /// 3-position of the flux window;
39  float targ_pos_beam[3]; /// Frame conversion from beam to flux frame
40  float pot; /// pot in the file
41 
42  float windowbase[3]; /// Flux window base vector
43  float windowdir1[3]; /// Flux window direction vector 1
44  float windowdir2[3]; /// Flux window direction vector 2
45 };
46 
47 class TTree;
48 class TFile;
49 
50 namespace fluxr {
52  {
53  public:
56 
57  const Long64_t GetEntries() {return fNEntries;};
58  const int GetRun() {return fRun;};
59  const void SetRun(int run) {fRun = run;};
60  const float GetPOT() {return fPOT;};
61  const TLorentzVector GetNuPosition() {return fNuPos;};
62  const TLorentzVector GetNuMomentum() {return fNuMom;};
63 
64  void SetRootFile(TFile* rootFileName);
65  bool FillMCFlux(Long64_t ientry, simb::MCFlux& mcflux);
66 
67  private:
68  TTree* fBNBtree;
69  TTree* fWindowtree;
72  double fMaxWeight;
73  double fMinWeight;
74  double fSumWeights;
75 
76  Long64_t fNEntries;
77  int fRun;
78  float fPOT;
79  TLorentzVector fNuPos;
80  TLorentzVector fNuMom;
81  };
82 }
83 
84 #endif // _BOONEINTERFACE_H_
float ini_t[20]
E of particle in before chain, array length 'npart'.
int npart
1,2,3,4: nue, nuebar, numu, numubar
float ini_pos[20][3]
id of each particle in before chain, array length 'npart'
float windowbase[3]
pot in the file
float ini_eng[20]
3-mom of particle in before chain, array length 'npart'
const TLorentzVector GetNuPosition()
float pot
Frame conversion from beam to flux frame.
const Long64_t GetEntries()
bool FillMCFlux(Long64_t ientry, simb::MCFlux &mcflux)
const void SetRun(int run)
void SetRootFile(TFile *rootFileName)
const TLorentzVector GetNuMomentum()
float fin_mom[20][3]
float windowdir2[3]
Flux window direction vector 1.
float windowdir1[3]
Flux window base vector.
float ini_mom[20][3]
3-pos of particle in before chain, array length 'npart'
float tank_pos_beam[3]
float targ_pos_beam[3]
3-position of the flux window;
const float GetPOT()
float fin_pol[20][3]
final 3-mom of particle in before chain, array length 'npart'
TLorentzVector fNuMom
int ntp
Magic Weight: CrossSect * BooBeamNT.
TLorentzVector fNuPos