Fork me on GitHub

source: git/readers/DelphesPythia8.cpp@ b94aacf

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

add event information

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