Fork me on GitHub

source: git/readers/DelphesPythia8.cpp@ 6e0bbe3

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 6e0bbe3 was 474eb76, checked in by pavel <pavel@…>, 11 years ago

fix DelphesPythia8.cpp

  • Property mode set to 100644
File size: 6.3 KB
RevLine 
[5150e39]1#include <stdexcept>
2#include <iostream>
3#include <sstream>
4
5#include <signal.h>
6
[474eb76]7#include "Pythia.h"
[5150e39]8
9#include "TROOT.h"
10#include "TApplication.h"
11
12#include "TFile.h"
13#include "TObjArray.h"
14#include "TStopwatch.h"
15#include "TDatabasePDG.h"
16#include "TParticlePDG.h"
17#include "TLorentzVector.h"
18
19#include "modules/Delphes.h"
20#include "classes/DelphesClasses.h"
21#include "classes/DelphesFactory.h"
22
23#include "ExRootAnalysis/ExRootTreeWriter.h"
24#include "ExRootAnalysis/ExRootTreeBranch.h"
25#include "ExRootAnalysis/ExRootProgressBar.h"
26
27using namespace std;
28
29//---------------------------------------------------------------------------
30
31void ConvertInput(Pythia8::Pythia *pythia, DelphesFactory *factory, TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray)
32{
33 int i;
34
35 Candidate *candidate;
36 TDatabasePDG *pdg;
37 TParticlePDG *pdgParticle;
38 Int_t pdgCode;
39
40 Int_t pid, status;
41 Double_t px, py, pz, e;
[474eb76]42 Double_t x, y, z, t;
[5150e39]43
44 pdg = TDatabasePDG::Instance();
45
46 for(i = 0; i < pythia->event.size(); ++i)
47 {
48 Pythia8::Particle &particle = pythia->event[i];
49
50 pid = particle.id();
51 status = particle.status();
52 px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.e();
53 x = particle.xProd(); y = particle.yProd(); z = particle.zProd(); t = particle.tProd();
54
55 candidate = factory->NewCandidate();
56
57 candidate->PID = pid;
58 pdgCode = TMath::Abs(candidate->PID);
59
60 candidate->Status = status;
61
[474eb76]62 candidate->M1 = particle.mother1() - 1;
63 candidate->M2 = particle.mother2() - 1;
[5150e39]64
[474eb76]65 candidate->D1 = particle.daughter1() - 1;
66 candidate->D2 = particle.daughter2() - 1;
[5150e39]67
68 pdgParticle = pdg->GetParticle(pid);
69 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
70 candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
71
72 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
73
[474eb76]74 candidate->Position.SetXYZT(x, y, z, t);
[5150e39]75
76 allParticleOutputArray->Add(candidate);
77
78 if(!pdgParticle) continue;
79
80 if(status == 1)
81 {
82 stableParticleOutputArray->Add(candidate);
83 }
84 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
85 {
86 partonOutputArray->Add(candidate);
87 }
88 }
89}
90
91//---------------------------------------------------------------------------
92
93static bool interrupted = false;
94
95void SignalHandler(int sig)
96{
97 interrupted = true;
98}
99
100//---------------------------------------------------------------------------
101
102int main(int argc, char *argv[])
103{
104 char appName[] = "DelphesPythia8";
105 stringstream message;
106 TFile *outputFile = 0;
107 TStopwatch readStopWatch, procStopWatch;
108 ExRootTreeWriter *treeWriter = 0;
109 ExRootTreeBranch *branchEvent = 0;
110 ExRootConfReader *confReader = 0;
111 Delphes *modularDelphes = 0;
112 DelphesFactory *factory = 0;
113 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
114 Long64_t eventCounter, errorCounter;
115 Long64_t numberOfEvents, timesAllowErrors;
116
117 Pythia8::Pythia *pythia = 0;
118
119 if(argc != 4)
120 {
121 cout << " Usage: " << appName << " config_file" << " pythia_card" << " output_file" << endl;
122 cout << " config_file - configuration file in Tcl format," << endl;
123 cout << " pythia_card - Pythia8 configuration file," << endl;
124 cout << " output_file - output file in ROOT format." << endl;
125 return 1;
126 }
127
128 signal(SIGINT, SignalHandler);
129
130 gROOT->SetBatch();
131
132 int appargc = 1;
133 char *appargv[] = {appName};
134 TApplication app(appName, &appargc, appargv);
135
136 try
137 {
138 outputFile = TFile::Open(argv[3], "CREATE");
139
140 if(outputFile == NULL)
141 {
142 message << "can't create output file " << argv[3];
143 throw runtime_error(message.str());
144 }
145
146 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
147
148 branchEvent = treeWriter->NewBranch("Event", LHEFEvent::Class());
149
150 confReader = new ExRootConfReader;
151 confReader->ReadFile(argv[1]);
152
153 modularDelphes = new Delphes("Delphes");
154 modularDelphes->SetConfReader(confReader);
155 modularDelphes->SetTreeWriter(treeWriter);
156
157 factory = modularDelphes->GetFactory();
158 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
159 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
160 partonOutputArray = modularDelphes->ExportArray("partons");
161
162 modularDelphes->InitTask();
163
164 // Initialize pythia
165 pythia = new Pythia8::Pythia;
166
167 if(pythia == NULL)
168 {
169 throw runtime_error("can't create Pythia instance");
170 }
171
172 // Read in commands from configuration file
173 pythia->readFile(argv[2]);
174
175 // Extract settings to be used in the main program
176 numberOfEvents = pythia->mode("Main:numberOfEvents");
177 timesAllowErrors = pythia->mode("Main:timesAllowErrors");
178
179 pythia->init();
180
181 ExRootProgressBar progressBar(numberOfEvents - 1);
182
183 // Loop over all events
184 errorCounter = 0;
185 treeWriter->Clear();
186 modularDelphes->Clear();
187 readStopWatch.Start();
188 for(eventCounter = 0; eventCounter < numberOfEvents && !interrupted; ++eventCounter)
189 {
190 if(!pythia->next())
191 {
192 // If failure because reached end of file then exit event loop
193 if (pythia->info.atEndOfFile())
194 {
195 cerr << "Aborted since reached end of Les Houches Event File" << endl;
196 break;
197 }
198
199 // First few failures write off as "acceptable" errors, then quit
200 if (++errorCounter < timesAllowErrors) continue;
201 cerr << "Event generation aborted prematurely, owing to error!" << endl;
202 break;
203 }
204
205 readStopWatch.Stop();
206
207 procStopWatch.Start();
208 ConvertInput(pythia, factory, allParticleOutputArray, stableParticleOutputArray, partonOutputArray);
209 modularDelphes->ProcessTask();
210 procStopWatch.Stop();
211
212 treeWriter->Fill();
213
214 treeWriter->Clear();
215 modularDelphes->Clear();
216
217 readStopWatch.Start();
218 progressBar.Update(eventCounter);
219 }
220
221 progressBar.Finish();
222
223 pythia->statistics();
224
225 modularDelphes->FinishTask();
226 treeWriter->Write();
227
228 cout << "** Exiting..." << endl;
229
230 delete pythia;
231 delete modularDelphes;
232 delete confReader;
233 delete treeWriter;
234 delete outputFile;
235
236 return 0;
237 }
238 catch(runtime_error &e)
239 {
240 if(treeWriter) delete treeWriter;
241 if(outputFile) delete outputFile;
242 cerr << "** ERROR: " << e.what() << endl;
243 return 1;
244 }
245}
Note: See TracBrowser for help on using the repository browser.