SND@LHC Software
Loading...
Searching...
No Matches
NtupleGenerator_FLUKA Class Reference

#include <NtupleGenerator_FLUKA.h>

Inheritance diagram for NtupleGenerator_FLUKA:
Collaboration diagram for NtupleGenerator_FLUKA:

Public Member Functions

 NtupleGenerator_FLUKA ()
 
virtual ~NtupleGenerator_FLUKA ()
 
Bool_t ReadEvent (FairPrimaryGenerator *)
 
virtual Bool_t Init (const char *, int)
 
virtual Bool_t Init (const char *)
 
Int_t GetNevents ()
 
void SetZ (Double_t X)
 

Protected Member Functions

 ClassDef (NtupleGenerator_FLUKA, 3)
 

Protected Attributes

Double_t id [1]
 
Double_t generation [1]
 
Double_t E [1]
 
Double_t t [1]
 
Double_t px [1]
 
Double_t py [1]
 
Double_t pz [1]
 
Double_t x [1]
 
Double_t y [1]
 
Double_t z [1]
 
Double_t w [1]
 
Double_t SND_Z
 
TFile * fInputFile
 
TTree * fTree
 
int fNevents
 
int fn
 
TClonesArray * primaries
 

Detailed Description

Definition at line 16 of file NtupleGenerator_FLUKA.h.

Constructor & Destructor Documentation

◆ NtupleGenerator_FLUKA()

NtupleGenerator_FLUKA::NtupleGenerator_FLUKA ( )

default constructor

Definition at line 12 of file NtupleGenerator_FLUKA.cxx.

12{}

◆ ~NtupleGenerator_FLUKA()

NtupleGenerator_FLUKA::~NtupleGenerator_FLUKA ( )
virtual

destructor

Definition at line 60 of file NtupleGenerator_FLUKA.cxx.

61{
62 // cout << "destroy Ntuple" << endl;
63 fInputFile->Close();
64 fInputFile->Delete();
65 delete fInputFile;
66}

Member Function Documentation

◆ ClassDef()

NtupleGenerator_FLUKA::ClassDef ( NtupleGenerator_FLUKA  ,
 
)
protected

◆ GetNevents()

Int_t NtupleGenerator_FLUKA::GetNevents ( )

Definition at line 106 of file NtupleGenerator_FLUKA.cxx.

107{
108 return fNevents;
109}

◆ Init() [1/2]

Bool_t NtupleGenerator_FLUKA::Init ( const char *  fileName)
virtual

Definition at line 15 of file NtupleGenerator_FLUKA.cxx.

15 {
16 return Init(fileName, 0);
17}
virtual Bool_t Init(const char *, int)

◆ Init() [2/2]

Bool_t NtupleGenerator_FLUKA::Init ( const char *  fileName,
int  firstEvent 
)
virtual

Definition at line 19 of file NtupleGenerator_FLUKA.cxx.

19 {
20 TString fN = "";
21 if (0 == strncmp("/eos",fileName,4) ) {
22 fN = gSystem->Getenv("EOSSHIP");
23 }
24 fN+=fileName;
25 fInputFile = TFile::Open(fN);
26 if (fInputFile->IsZombie()) {
27 LOG(FATAL) << "NtupleGenerator_FLUKA: Error opening the Signal file" << fN;
28 return kFALSE;
29 }
30 LOG(INFO) << "NtupleGenerator_FLUKA: Opening input file " << fN;
31
32 fTree = (TTree *)fInputFile->Get("nt");
33
34 fNevents = fTree->GetEntries();
35 fn = firstEvent;
36 primaries = 0;
37 // Check the input TTree structure
38 if (fTree->FindLeaf("primaries_") != nullptr) {
39 fTree->SetBranchAddress("primaries",&primaries);
40 }
41 else {
42 fTree->SetBranchAddress("id",&id); // particle id
43 fTree->SetBranchAddress("generation",&generation); // origin generation number
44 fTree->SetBranchAddress("t",&t); // time of flight
45 fTree->SetBranchAddress("E",&E); // incoming muon energy
46 fTree->SetBranchAddress("w",&w); // weight of event
47 fTree->SetBranchAddress("x",&x); // position
48 fTree->SetBranchAddress("y",&y);
49 fTree->SetBranchAddress("z",&z);
50 fTree->SetBranchAddress("px",&px); // momentum
51 fTree->SetBranchAddress("py",&py);
52 fTree->SetBranchAddress("pz",&pz);
53 }
54 return kTRUE;
55}
int firstEvent
Definition MufluxReco.py:13

◆ ReadEvent()

Bool_t NtupleGenerator_FLUKA::ReadEvent ( FairPrimaryGenerator *  cpg)

public method ReadEvent

Definition at line 70 of file NtupleGenerator_FLUKA.cxx.

71{
72 fTree->GetEntry(fn);
73 fn++;
74 if (fn %10000==0) {LOG(INFO)<< "reading event "<<fn;}
75 if (fn > fNevents) {
76 LOG(WARNING) << "No more input events";
77 return kFALSE; }
78 /* Add all primary tracks in the event. Typically, muon bkg MC has single muon per event.
79 Several muon primaries per event exist in the multi-muon MC.
80 */
81 if (fTree->FindLeaf("primaries_") != nullptr) {
82 int nf = primaries->GetEntries();
83 for(int i=0; i<nf; i++){
84 //casting
85 struct PrimaryTrack* primary = dynamic_cast<struct PrimaryTrack*>(primaries->AddrAt(i));
86 // what to do with generation info? convert time from ns to sec
87 cpg->AddTrack(int(primary->id),primary->px,primary->py,primary->pz,
88 primary->x,primary->y,primary->z-SND_Z,-1,true,primary->E,
89 primary->t/1E9,primary->w,(TMCProcess)primary->generation);
90 LOG(DEBUG) << "NtupleGenerator_FLUKA: add muon " << i << "," << int(primary->id)
91 <<","<<primary->px<<","<<primary->py<<","<<primary->pz<<","<<primary->x
92 <<","<<primary->y<<","<<primary->z-SND_Z<<","<<primary->generation
93 <<","<<primary->E<<","<<primary->t<<"ns,"<<primary->w;
94 }
95 }
96 else {
97 // what to do with generation info? convert time from ns to sec
98 cpg->AddTrack(int(id[0]),px[0],py[0],pz[0],x[0],y[0],z[0]-SND_Z,-1,true,E[0],t[0]/1E9,w[0],(TMCProcess)generation[0]);
99 LOG(DEBUG) << "NtupleGenerator_FLUKA: add muon " << id[0]<<","<<px[0]<<","<<py[0]<<","<<pz[0]<<","<<x[0]<<","<<y[0]<<","<<z[0]-SND_Z<<","<<generation[0]<<","<<E[0]<<","<<t[0]<<"ns,"<<w[0];
100 }
101
102 return kTRUE;
103}
int i
Definition ShipAna.py:86

◆ SetZ()

void NtupleGenerator_FLUKA::SetZ ( Double_t  X)
inline

Definition at line 31 of file NtupleGenerator_FLUKA.h.

Member Data Documentation

◆ E

Double_t NtupleGenerator_FLUKA::E[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ fInputFile

TFile* NtupleGenerator_FLUKA::fInputFile
protected

Definition at line 36 of file NtupleGenerator_FLUKA.h.

◆ fn

int NtupleGenerator_FLUKA::fn
protected

Definition at line 39 of file NtupleGenerator_FLUKA.h.

◆ fNevents

int NtupleGenerator_FLUKA::fNevents
protected

Definition at line 38 of file NtupleGenerator_FLUKA.h.

◆ fTree

TTree* NtupleGenerator_FLUKA::fTree
protected

Definition at line 37 of file NtupleGenerator_FLUKA.h.

◆ generation

Double_t NtupleGenerator_FLUKA::generation[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ id

Double_t NtupleGenerator_FLUKA::id[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ primaries

TClonesArray* NtupleGenerator_FLUKA::primaries
protected

Definition at line 40 of file NtupleGenerator_FLUKA.h.

◆ px

Double_t NtupleGenerator_FLUKA::px[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ py

Double_t NtupleGenerator_FLUKA::py[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ pz

Double_t NtupleGenerator_FLUKA::pz[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ SND_Z

Double_t NtupleGenerator_FLUKA::SND_Z
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ t

Double_t NtupleGenerator_FLUKA::t[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ w

Double_t NtupleGenerator_FLUKA::w[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ x

Double_t NtupleGenerator_FLUKA::x[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ y

Double_t NtupleGenerator_FLUKA::y[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.

◆ z

Double_t NtupleGenerator_FLUKA::z[1]
protected

Definition at line 35 of file NtupleGenerator_FLUKA.h.


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