FairShip
Loading...
Searching...
No Matches
EvtCalcGenerator Class Reference

#include <EvtCalcGenerator.h>

Inheritance diagram for EvtCalcGenerator:
Collaboration diagram for EvtCalcGenerator:

Public Member Functions

 EvtCalcGenerator ()
 
 ~EvtCalcGenerator () override
 
Bool_t ReadEvent (FairPrimaryGenerator *) override
 
Bool_t Init (const char *, int) override
 
Bool_t Init (const char *) override
 
Int_t GetNevents ()
 
void SetPositions (Double_t zTa, Double_t zDV)
 
Double_t GetNdaughters (const std::unique_ptr< TTree > &) const
 
Double_t GetMotherPx (const std::unique_ptr< TTree > &) const
 
Double_t GetMotherPy (const std::unique_ptr< TTree > &) const
 
Double_t GetMotherPz (const std::unique_ptr< TTree > &) const
 
Double_t GetMotherE (const std::unique_ptr< TTree > &) const
 
Double_t GetVx (const std::unique_ptr< TTree > &) const
 
Double_t GetVy (const std::unique_ptr< TTree > &) const
 
Double_t GetVz (const std::unique_ptr< TTree > &) const
 
Double_t GetDecayProb (const std::unique_ptr< TTree > &) const
 
Double_t GetDauPx (const std::unique_ptr< TTree > &, int) const
 
Double_t GetDauPy (const std::unique_ptr< TTree > &, int) const
 
Double_t GetDauPz (const std::unique_ptr< TTree > &, int) const
 
Double_t GetDauE (const std::unique_ptr< TTree > &, int) const
 
Double_t GetDauPDG (const std::unique_ptr< TTree > &, int) const
 
- Public Member Functions inherited from SHiP::Generator
 Generator ()
 
virtual ~Generator ()
 
virtual Bool_t Init (const char *, int)=0
 
virtual Bool_t Init (const char *)=0
 
virtual Bool_t Init (const std::vector< std::string > &inFiles, int startNumber)
 
virtual Bool_t Init (const std::vector< std::string > &inFiles)
 
virtual void UseExternalFile (std::string x, Int_t i)
 
virtual void UseExternalFile (std::vector< std::string > &inFiles, Int_t i)
 

Protected Types

enum class  BranchIndices {
  MotherPx = 0 , MotherPy = 1 , MotherPz = 2 , MotherE = 3 ,
  DecayProb = 6 , Vx = 7 , Vy = 8 , Vz = 9
}
 

Protected Member Functions

Double_t GetBranchValue (const std::unique_ptr< TTree > &, unsigned) const
 
Double_t GetDaughterValue (const std::unique_ptr< TTree > &tree, int, int) const
 

Protected Attributes

Double_t ztarget
 
Double_t zDecayVolume
 
std::unique_ptr< TFile > fInputFile
 
std::unique_ptr< TTree > fTree
 
std::vector< double > branchVars
 
int fNevents
 
int fn
 
int nBranches
 
int Ndau
 
- Protected Attributes inherited from SHiP::Generator
std::optional< std::string > fextFile
 
Int_t firstEvent = 0
 

Detailed Description

Definition at line 16 of file EvtCalcGenerator.h.

Member Enumeration Documentation

◆ BranchIndices

enum class EvtCalcGenerator::BranchIndices
strongprotected
Enumerator
MotherPx 
MotherPy 
MotherPz 
MotherE 
DecayProb 
Vx 
Vy 
Vz 

Definition at line 65 of file EvtCalcGenerator.h.

Constructor & Destructor Documentation

◆ EvtCalcGenerator()

EvtCalcGenerator::EvtCalcGenerator ( )

default constructor

Definition at line 17 of file EvtCalcGenerator.cxx.

17{}

◆ ~EvtCalcGenerator()

EvtCalcGenerator::~EvtCalcGenerator ( )
override

destructor

Definition at line 47 of file EvtCalcGenerator.cxx.

47{}

Member Function Documentation

◆ GetBranchValue()

Double_t EvtCalcGenerator::GetBranchValue ( const std::unique_ptr< TTree > &  tree,
unsigned  index 
) const
protected

Definition at line 50 of file EvtCalcGenerator.cxx.

51 {
52 if (index < branchVars.size()) {
53 return branchVars[index];
54 } else {
55 throw std::out_of_range("Branch index out of range");
56 }
57}
std::vector< double > branchVars

◆ GetDauE()

Double_t EvtCalcGenerator::GetDauE ( const std::unique_ptr< TTree > &  tree,
int  dauID 
) const

Definition at line 120 of file EvtCalcGenerator.cxx.

121 {
122 return GetDaughterValue(tree, dauID, 3);
123}
Double_t GetDaughterValue(const std::unique_ptr< TTree > &tree, int, int) const

◆ GetDaughterValue()

Double_t EvtCalcGenerator::GetDaughterValue ( const std::unique_ptr< TTree > &  tree,
int  dauID,
int  offset 
) const
protected

Definition at line 59 of file EvtCalcGenerator.cxx.

60 {
61 int baseIndex = 10 + (dauID * 6);
62 return GetBranchValue(tree, baseIndex + offset);
63}
Double_t GetBranchValue(const std::unique_ptr< TTree > &, unsigned) const

◆ GetDauPDG()

Double_t EvtCalcGenerator::GetDauPDG ( const std::unique_ptr< TTree > &  tree,
int  dauID 
) const

Definition at line 124 of file EvtCalcGenerator.cxx.

125 {
126 return GetDaughterValue(tree, dauID, 5);
127}

◆ GetDauPx()

Double_t EvtCalcGenerator::GetDauPx ( const std::unique_ptr< TTree > &  tree,
int  dauID 
) const

Definition at line 108 of file EvtCalcGenerator.cxx.

109 {
110 return GetDaughterValue(tree, dauID, 0);
111}

◆ GetDauPy()

Double_t EvtCalcGenerator::GetDauPy ( const std::unique_ptr< TTree > &  tree,
int  dauID 
) const

Definition at line 112 of file EvtCalcGenerator.cxx.

113 {
114 return GetDaughterValue(tree, dauID, 1);
115}

◆ GetDauPz()

Double_t EvtCalcGenerator::GetDauPz ( const std::unique_ptr< TTree > &  tree,
int  dauID 
) const

Definition at line 116 of file EvtCalcGenerator.cxx.

117 {
118 return GetDaughterValue(tree, dauID, 2);
119}

◆ GetDecayProb()

Double_t EvtCalcGenerator::GetDecayProb ( const std::unique_ptr< TTree > &  tree) const

Definition at line 102 of file EvtCalcGenerator.cxx.

103 {
104 return GetBranchValue(tree, static_cast<int>(BranchIndices::DecayProb));
105}

◆ GetMotherE()

Double_t EvtCalcGenerator::GetMotherE ( const std::unique_ptr< TTree > &  tree) const

Definition at line 85 of file EvtCalcGenerator.cxx.

86 {
87 return GetBranchValue(tree, static_cast<int>(BranchIndices::MotherE));
88}

◆ GetMotherPx()

Double_t EvtCalcGenerator::GetMotherPx ( const std::unique_ptr< TTree > &  tree) const

Definition at line 73 of file EvtCalcGenerator.cxx.

74 {
75 return GetBranchValue(tree, static_cast<int>(BranchIndices::MotherPx));
76}

◆ GetMotherPy()

Double_t EvtCalcGenerator::GetMotherPy ( const std::unique_ptr< TTree > &  tree) const

Definition at line 77 of file EvtCalcGenerator.cxx.

78 {
79 return GetBranchValue(tree, static_cast<int>(BranchIndices::MotherPy));
80}

◆ GetMotherPz()

Double_t EvtCalcGenerator::GetMotherPz ( const std::unique_ptr< TTree > &  tree) const

Definition at line 81 of file EvtCalcGenerator.cxx.

82 {
83 return GetBranchValue(tree, static_cast<int>(BranchIndices::MotherPz));
84}

◆ GetNdaughters()

Double_t EvtCalcGenerator::GetNdaughters ( const std::unique_ptr< TTree > &  tree) const

Definition at line 67 of file EvtCalcGenerator.cxx.

68 {
69 return GetBranchValue(tree, nBranches - 1);
70}

◆ GetNevents()

Int_t EvtCalcGenerator::GetNevents ( )
inline

Definition at line 30 of file EvtCalcGenerator.h.

30{ return fNevents; }

◆ GetVx()

Double_t EvtCalcGenerator::GetVx ( const std::unique_ptr< TTree > &  tree) const

Definition at line 91 of file EvtCalcGenerator.cxx.

91 {
92 return GetBranchValue(tree, static_cast<int>(BranchIndices::Vx));
93}

◆ GetVy()

Double_t EvtCalcGenerator::GetVy ( const std::unique_ptr< TTree > &  tree) const

Definition at line 94 of file EvtCalcGenerator.cxx.

94 {
95 return GetBranchValue(tree, static_cast<int>(BranchIndices::Vy));
96}

◆ GetVz()

Double_t EvtCalcGenerator::GetVz ( const std::unique_ptr< TTree > &  tree) const

Definition at line 97 of file EvtCalcGenerator.cxx.

97 {
98 return GetBranchValue(tree, static_cast<int>(BranchIndices::Vz));
99}

◆ Init() [1/2]

Bool_t EvtCalcGenerator::Init ( const char *  fileName)
overridevirtual

Implements SHiP::Generator.

Definition at line 20 of file EvtCalcGenerator.cxx.

20 {
21 return Init(fileName, 0);
22}
Bool_t Init(const char *, int) override

◆ Init() [2/2]

Bool_t EvtCalcGenerator::Init ( const char *  fileName,
int  startEvent 
)
overridevirtual

Implements SHiP::Generator.

Definition at line 24 of file EvtCalcGenerator.cxx.

24 {
25 fInputFile = std::unique_ptr<TFile>(TFile::Open(fileName, "read"));
26 LOGF(info, "Info EvtCalcGenerator: Opening input file %s", fileName);
27
28 fTree =
29 std::unique_ptr<TTree>(dynamic_cast<TTree*>(fInputFile->Get("LLP_tree")));
30 fNevents = fTree->GetEntries();
31 fn = startEvent;
32
33 auto* branches = fTree->GetListOfBranches();
34 nBranches = branches->GetEntries();
35 branchVars.resize(nBranches);
36
37 for (int i = 0; i < nBranches; ++i) {
38 auto* branch = dynamic_cast<TBranch*>(branches->At(i));
39 if (fTree->FindBranch(branch->GetName())) {
40 fTree->SetBranchAddress(branch->GetName(), &branchVars[i]);
41 }
42 }
43
44 return kTRUE;
45}
std::unique_ptr< TFile > fInputFile
std::unique_ptr< TTree > fTree
int i
Definition: ShipAna.py:97

◆ ReadEvent()

Bool_t EvtCalcGenerator::ReadEvent ( FairPrimaryGenerator *  cpg)
override

Definition at line 130 of file EvtCalcGenerator.cxx.

130 {
131 if (fn == fNevents) {
132 LOG(warning) << "End of input file. Rewind.";
133 fn = 0;
134 }
135
136 fTree->GetEntry(fn);
137 fn++;
138 if (fn % 100 == 0) LOGF(info, "Info EvtCalcGenerator: event nr %d", fn);
139
141 // Vertex coordinates in the SHiP reference frame, expressed in [cm]
142 Double_t space_unit_conv = 100.; // m to cm
143 Double_t coord_shift = (zDecayVolume - ztarget) / space_unit_conv; // units m
144 Double_t vx_transf = GetVx(fTree) * space_unit_conv; // units cm
145 Double_t vy_transf = GetVy(fTree) * space_unit_conv; // units cm
146 Double_t vz_transf =
147 (GetVz(fTree) - coord_shift) * space_unit_conv; // units cm
148
149 Double_t c = 2.99792458e+10; // speed of light [cm/s]
150 Double_t tof = TMath::Sqrt(vx_transf * vx_transf + vy_transf * vy_transf +
151 vz_transf * vz_transf) /
152 c;
153 Double_t decay_prob = GetDecayProb(fTree);
154 Double_t pdg_dau = 0;
155
156 // Mother LLP
157 Bool_t wanttracking = false;
158 Double_t pdg_llp = 999; // Geantino, placeholder
159
160 cpg->AddTrack(pdg_llp, GetMotherPx(fTree), GetMotherPy(fTree),
161 GetMotherPz(fTree), vx_transf, vy_transf, vz_transf, -1.,
162 wanttracking, GetMotherE(fTree), tof, decay_prob);
163
164 wanttracking = true;
165
166 // Secondaries
167 for (int iPart = 0; iPart < Ndau; ++iPart) {
168 pdg_dau = GetDauPDG(fTree, iPart);
169 if (pdg_dau != -999) {
170 cpg->AddTrack(pdg_dau, GetDauPx(fTree, iPart), GetDauPy(fTree, iPart),
171 GetDauPz(fTree, iPart), vx_transf, vy_transf, vz_transf, 0.,
172 wanttracking, GetDauE(fTree, iPart), tof, decay_prob);
173 }
174 }
175
176 return kTRUE;
177}
Double_t GetDauPDG(const std::unique_ptr< TTree > &, int) const
Double_t GetNdaughters(const std::unique_ptr< TTree > &) const
Double_t GetDauE(const std::unique_ptr< TTree > &, int) const
Double_t GetMotherE(const std::unique_ptr< TTree > &) const
Double_t GetDauPx(const std::unique_ptr< TTree > &, int) const
Double_t GetDecayProb(const std::unique_ptr< TTree > &) const
Double_t GetDauPz(const std::unique_ptr< TTree > &, int) const
Double_t GetMotherPy(const std::unique_ptr< TTree > &) const
Double_t GetMotherPz(const std::unique_ptr< TTree > &) const
Double_t GetMotherPx(const std::unique_ptr< TTree > &) const
Double_t GetVz(const std::unique_ptr< TTree > &) const
Double_t GetDauPy(const std::unique_ptr< TTree > &, int) const
Double_t GetVx(const std::unique_ptr< TTree > &) const
Double_t GetVy(const std::unique_ptr< TTree > &) const
constants c
Definition: hnl.py:115

◆ SetPositions()

void EvtCalcGenerator::SetPositions ( Double_t  zTa,
Double_t  zDV 
)
inline

Definition at line 31 of file EvtCalcGenerator.h.

31 {
32 ztarget = zTa; // units cm (midpoint)
33 zDecayVolume = zDV; // units cm (midpoint)
34 }

Member Data Documentation

◆ branchVars

std::vector<double> EvtCalcGenerator::branchVars
protected

Definition at line 63 of file EvtCalcGenerator.h.

◆ fInputFile

std::unique_ptr<TFile> EvtCalcGenerator::fInputFile
protected

Definition at line 61 of file EvtCalcGenerator.h.

◆ fn

int EvtCalcGenerator::fn
protected

Definition at line 77 of file EvtCalcGenerator.h.

◆ fNevents

int EvtCalcGenerator::fNevents
protected

Definition at line 76 of file EvtCalcGenerator.h.

◆ fTree

std::unique_ptr<TTree> EvtCalcGenerator::fTree
protected

Definition at line 62 of file EvtCalcGenerator.h.

◆ nBranches

int EvtCalcGenerator::nBranches
protected

Definition at line 78 of file EvtCalcGenerator.h.

◆ Ndau

int EvtCalcGenerator::Ndau
protected

Definition at line 79 of file EvtCalcGenerator.h.

◆ zDecayVolume

Double_t EvtCalcGenerator::zDecayVolume
protected

Definition at line 60 of file EvtCalcGenerator.h.

◆ ztarget

Double_t EvtCalcGenerator::ztarget
protected

Definition at line 60 of file EvtCalcGenerator.h.


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