Fork me on GitHub

source: git/readers/DelphesPythia8.cpp@ 31def62

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 31def62 was c0f9b5b, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

check if can read Pythia8 configuration file

  • Property mode set to 100644
File size: 9.6 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19#include <stdexcept>
20#include <iostream>
21#include <sstream>
22#include <string>
23
24#include <signal.h>
25
26#include "Pythia.h"
27
28#include "TROOT.h"
29#include "TApplication.h"
30
31#include "TFile.h"
32#include "TObjArray.h"
33#include "TStopwatch.h"
34#include "TDatabasePDG.h"
35#include "TParticlePDG.h"
36#include "TLorentzVector.h"
37
38#include "modules/Delphes.h"
39#include "classes/DelphesClasses.h"
40#include "classes/DelphesFactory.h"
41#include "classes/DelphesLHEFReader.h"
42
43#include "ExRootAnalysis/ExRootTreeWriter.h"
44#include "ExRootAnalysis/ExRootTreeBranch.h"
45#include "ExRootAnalysis/ExRootProgressBar.h"
46
47using namespace std;
48
49//---------------------------------------------------------------------------
50
51void ConvertInput(Long64_t eventCounter, Pythia8::Pythia *pythia,
52 ExRootTreeBranch *branch, DelphesFactory *factory,
53 TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray,
54 TStopwatch *readStopWatch, TStopwatch *procStopWatch)
55{
56 int i;
57
58 HepMCEvent *element;
59 Candidate *candidate;
60 TDatabasePDG *pdg;
61 TParticlePDG *pdgParticle;
62 Int_t pdgCode;
63
64 Int_t pid, status;
65 Double_t px, py, pz, e, mass;
66 Double_t x, y, z, t;
67
68 // event information
69 element = static_cast<HepMCEvent *>(branch->NewEntry());
70
71 element->Number = eventCounter;
72
73 element->ProcessID = pythia->info.code();
74 element->MPI = 1;
75 element->Weight = pythia->info.weight();
76 element->Scale = pythia->info.QRen();
77 element->AlphaQED = pythia->info.alphaEM();
78 element->AlphaQCD = pythia->info.alphaS();
79
80 element->ID1 = pythia->info.id1();
81 element->ID2 = pythia->info.id2();
82 element->X1 = pythia->info.x1();
83 element->X2 = pythia->info.x2();
84 element->ScalePDF = pythia->info.QFac();
85 element->PDF1 = pythia->info.pdf1();
86 element->PDF2 = pythia->info.pdf2();
87
88 element->ReadTime = readStopWatch->RealTime();
89 element->ProcTime = procStopWatch->RealTime();
90
91 pdg = TDatabasePDG::Instance();
92
93 for(i = 1; i < pythia->event.size(); ++i)
94 {
95 Pythia8::Particle &particle = pythia->event[i];
96
97 pid = particle.id();
98 status = particle.statusHepMC();
99 px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.e(); mass = particle.m();
100 x = particle.xProd(); y = particle.yProd(); z = particle.zProd(); t = particle.tProd();
101
102 candidate = factory->NewCandidate();
103
104 candidate->PID = pid;
105 pdgCode = TMath::Abs(candidate->PID);
106
107 candidate->Status = status;
108
109 candidate->M1 = particle.mother1() - 1;
110 candidate->M2 = particle.mother2() - 1;
111
112 candidate->D1 = particle.daughter1() - 1;
113 candidate->D2 = particle.daughter2() - 1;
114
115 pdgParticle = pdg->GetParticle(pid);
116 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
117 candidate->Mass = mass;
118
119 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
120
121 candidate->Position.SetXYZT(x, y, z, t);
122
123 allParticleOutputArray->Add(candidate);
124
125 if(!pdgParticle) continue;
126
127 if(status == 1)
128 {
129 stableParticleOutputArray->Add(candidate);
130 }
131 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
132 {
133 partonOutputArray->Add(candidate);
134 }
135 }
136}
137
138//---------------------------------------------------------------------------
139
140static bool interrupted = false;
141
142void SignalHandler(int sig)
143{
144 interrupted = true;
145}
146
147//---------------------------------------------------------------------------
148
149int main(int argc, char *argv[])
150{
151 char appName[] = "DelphesPythia8";
152 stringstream message;
153 FILE *inputFile = 0;
154 TFile *outputFile = 0;
155 TStopwatch readStopWatch, procStopWatch;
156 ExRootTreeWriter *treeWriter = 0;
157 ExRootTreeBranch *branchEvent = 0;
158 ExRootTreeBranch *branchEventLHEF = 0, *branchWeightLHEF = 0;
159 ExRootConfReader *confReader = 0;
160 Delphes *modularDelphes = 0;
161 DelphesFactory *factory = 0;
162 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
163 TObjArray *stableParticleOutputArrayLHEF = 0, *allParticleOutputArrayLHEF = 0, *partonOutputArrayLHEF = 0;
164 DelphesLHEFReader *reader = 0;
165 Long64_t eventCounter, errorCounter;
166 Long64_t numberOfEvents, timesAllowErrors;
167
168 Pythia8::Pythia *pythia = 0;
169
170 if(argc != 4)
171 {
172 cout << " Usage: " << appName << " config_file" << " pythia_card" << " output_file" << endl;
173 cout << " config_file - configuration file in Tcl format," << endl;
174 cout << " pythia_card - Pythia8 configuration file," << endl;
175 cout << " output_file - output file in ROOT format." << endl;
176 return 1;
177 }
178
179 signal(SIGINT, SignalHandler);
180
181 gROOT->SetBatch();
182
183 int appargc = 1;
184 char *appargv[] = {appName};
185 TApplication app(appName, &appargc, appargv);
186
187 try
188 {
189 outputFile = TFile::Open(argv[3], "CREATE");
190
191 if(outputFile == NULL)
192 {
193 message << "can't create output file " << argv[3];
194 throw runtime_error(message.str());
195 }
196
197 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
198
199 branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
200
201 confReader = new ExRootConfReader;
202 confReader->ReadFile(argv[1]);
203
204 modularDelphes = new Delphes("Delphes");
205 modularDelphes->SetConfReader(confReader);
206 modularDelphes->SetTreeWriter(treeWriter);
207
208 factory = modularDelphes->GetFactory();
209 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
210 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
211 partonOutputArray = modularDelphes->ExportArray("partons");
212
213 // Initialize Pythia
214 pythia = new Pythia8::Pythia;
215
216 if(pythia == NULL)
217 {
218 throw runtime_error("can't create Pythia instance");
219 }
220
221 // Read in commands from configuration file
222 if(!pythia->readFile(argv[2]))
223 {
224 message << "can't read Pythia8 configuration file " << argv[2] << endl;
225 throw runtime_error(message.str());
226 }
227
228 // Extract settings to be used in the main program
229 numberOfEvents = pythia->mode("Main:numberOfEvents");
230 timesAllowErrors = pythia->mode("Main:timesAllowErrors");
231
232 inputFile = fopen(pythia->word("Beams:LHEF").c_str(), "r");
233 if(inputFile)
234 {
235 reader = new DelphesLHEFReader;
236 reader->SetInputFile(inputFile);
237
238 branchEventLHEF = treeWriter->NewBranch("EventLHEF", LHEFEvent::Class());
239 branchWeightLHEF = treeWriter->NewBranch("WeightLHEF", LHEFWeight::Class());
240
241 allParticleOutputArrayLHEF = modularDelphes->ExportArray("allParticlesLHEF");
242 stableParticleOutputArrayLHEF = modularDelphes->ExportArray("stableParticlesLHEF");
243 partonOutputArrayLHEF = modularDelphes->ExportArray("partonsLHEF");
244 }
245
246 modularDelphes->InitTask();
247
248 pythia->init();
249
250 // ExRootProgressBar progressBar(numberOfEvents - 1);
251 ExRootProgressBar progressBar(-1);
252
253 // Loop over all events
254 errorCounter = 0;
255 treeWriter->Clear();
256 modularDelphes->Clear();
257 readStopWatch.Start();
258 for(eventCounter = 0; eventCounter < numberOfEvents && !interrupted; ++eventCounter)
259 {
260 while(reader && reader->ReadBlock(factory, allParticleOutputArrayLHEF,
261 stableParticleOutputArrayLHEF, partonOutputArrayLHEF) && !reader->EventReady());
262
263 if(!pythia->next())
264 {
265 // If failure because reached end of file then exit event loop
266 if(pythia->info.atEndOfFile())
267 {
268 cerr << "Aborted since reached end of Les Houches Event File" << endl;
269 break;
270 }
271
272 // First few failures write off as "acceptable" errors, then quit
273 if(++errorCounter > timesAllowErrors)
274 {
275 cerr << "Event generation aborted prematurely, owing to error!" << endl;
276 break;
277 }
278
279 modularDelphes->Clear();
280 reader->Clear();
281 continue;
282 }
283
284 readStopWatch.Stop();
285
286 procStopWatch.Start();
287 ConvertInput(eventCounter, pythia, branchEvent, factory,
288 allParticleOutputArray, stableParticleOutputArray, partonOutputArray,
289 &readStopWatch, &procStopWatch);
290 modularDelphes->ProcessTask();
291 procStopWatch.Stop();
292
293 if(reader)
294 {
295 reader->AnalyzeEvent(branchEventLHEF, eventCounter, &readStopWatch, &procStopWatch);
296 reader->AnalyzeWeight(branchWeightLHEF);
297 }
298
299 treeWriter->Fill();
300
301 treeWriter->Clear();
302 modularDelphes->Clear();
303 if(reader) reader->Clear();
304
305 readStopWatch.Start();
306 progressBar.Update(eventCounter, eventCounter);
307 }
308
309 progressBar.Update(eventCounter, eventCounter, kTRUE);
310 progressBar.Finish();
311
312 pythia->stat();
313
314 modularDelphes->FinishTask();
315 treeWriter->Write();
316
317 cout << "** Exiting..." << endl;
318
319 delete reader;
320 delete pythia;
321 delete modularDelphes;
322 delete confReader;
323 delete treeWriter;
324 delete outputFile;
325
326 return 0;
327 }
328 catch(runtime_error &e)
329 {
330 if(treeWriter) delete treeWriter;
331 if(outputFile) delete outputFile;
332 cerr << "** ERROR: " << e.what() << endl;
333 return 1;
334 }
335}
Note: See TracBrowser for help on using the repository browser.