SND@LHC Software
Loading...
Searching...
No Matches
NtupleGenerator_FLUKA.cxx
Go to the documentation of this file.
1#include <math.h>
2#include "TROOT.h"
3#include "TSystem.h"
4#include "TFile.h"
5#include "FairPrimaryGenerator.h"
7#include "FairLogger.h"
8
9// read events from ntuples produced by FLUKA
10
11// ----- Default constructor -------------------------------------------
13// -------------------------------------------------------------------------
14// ----- Default constructor -------------------------------------------
15Bool_t NtupleGenerator_FLUKA::Init(const char* fileName) {
16 return Init(fileName, 0);
17}
18// ----- Default constructor -------------------------------------------
19Bool_t NtupleGenerator_FLUKA::Init(const char* fileName, const int firstEvent) {
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}
56// -------------------------------------------------------------------------
57
58
59// ----- Destructor ----------------------------------------------------
61{
62 // cout << "destroy Ntuple" << endl;
63 fInputFile->Close();
64 fInputFile->Delete();
65 delete fInputFile;
66}
67// -------------------------------------------------------------------------
68
69// ----- Passing the event ---------------------------------------------
70Bool_t NtupleGenerator_FLUKA::ReadEvent(FairPrimaryGenerator* cpg)
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}
104
105// -------------------------------------------------------------------------
107{
108 return fNevents;
109}
110
111
virtual Bool_t Init(const char *, int)
Bool_t ReadEvent(FairPrimaryGenerator *)
ClassImp(ecalContFact) ecalContFact