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, 4)
 

Protected Attributes

Double_t runN
 
Double_t eventN
 
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 13 of file NtupleGenerator_FLUKA.cxx.

13{}

◆ ~NtupleGenerator_FLUKA()

NtupleGenerator_FLUKA::~NtupleGenerator_FLUKA ( )
virtual

destructor

Definition at line 63 of file NtupleGenerator_FLUKA.cxx.

64{
65 // cout << "destroy Ntuple" << endl;
66 fInputFile->Close();
67 fInputFile->Delete();
68 delete fInputFile;
69}

Member Function Documentation

◆ ClassDef()

NtupleGenerator_FLUKA::ClassDef ( NtupleGenerator_FLUKA  ,
 
)
protected

◆ GetNevents()

Int_t NtupleGenerator_FLUKA::GetNevents ( )

Definition at line 114 of file NtupleGenerator_FLUKA.cxx.

115{
116 return fNevents;
117}

◆ Init() [1/2]

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

Definition at line 16 of file NtupleGenerator_FLUKA.cxx.

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

◆ Init() [2/2]

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

Definition at line 20 of file NtupleGenerator_FLUKA.cxx.

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

◆ ReadEvent()

Bool_t NtupleGenerator_FLUKA::ReadEvent ( FairPrimaryGenerator *  cpg)

public method ReadEvent

Definition at line 73 of file NtupleGenerator_FLUKA.cxx.

74{
75 fTree->GetEntry(fn);
76 fn++;
77 if (fn %10000==0) {LOG(INFO)<< "reading event "<<fn;}
78 if (fn > fNevents) {
79 LOG(WARNING) << "No more input events";
80 return kFALSE; }
81 /* Add all primary tracks in the event. Typically, muon bkg MC has single muon per event.
82 Several muon primaries per event exist in the multi-muon MC.
83 */
84 if (fTree->FindLeaf("primaries_") != nullptr) {
85 int nf = primaries->GetEntries();
86 for(int i=0; i<nf; i++){
87 //casting
88 struct PrimaryTrack* primary = dynamic_cast<struct PrimaryTrack*>(primaries->AddrAt(i));
89 //convert time from ns to sec
90 cpg->AddTrack(int(primary->id),primary->px,primary->py,primary->pz,
91 primary->x,primary->y,primary->z-SND_Z,-1,true,primary->E,
92 primary->t/1E9,primary->w,(TMCProcess)primary->generation);
93 LOG(DEBUG) << "NtupleGenerator_FLUKA: add muon " << i << "," << int(primary->id)
94 <<","<<primary->px<<","<<primary->py<<","<<primary->pz<<","<<primary->x
95 <<","<<primary->y<<","<<primary->z-SND_Z<<","<<primary->generation
96 <<","<<primary->E<<","<<primary->t<<"ns,"<<primary->w;
97 }
98 }
99 else {
100 // convert time from ns to sec
101 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]);
102 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];
103 }
104
105 // Set the generation run and event numbers through the MC event header
106 FairMCEventHeader* header = cpg->GetEvent();
107 header->SetEventID(int(eventN));
108 header->SetRunID(int(runN));
109
110 return kTRUE;
111}
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.

◆ eventN

Double_t NtupleGenerator_FLUKA::eventN
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.

◆ runN

Double_t NtupleGenerator_FLUKA::runN
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: