Fork me on GitHub

source: git/readers/DelphesPythia8.cpp@ 971db5b

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 971db5b was 971db5b, checked in by Alexandre Mertens <alexandre.mertens@…>, 8 years ago

first version working with particles (not jet) and fixed value of Energy

  • Property mode set to 100644
File size: 11.1 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// Single-particle gun. The particle must be a colour singlet.
149// Input: flavour, energy, direction (theta, phi).
150// If theta < 0 then random choice over solid angle.
151// Optional final argument to put particle at rest => E = m.
152// from pythia8 example 21
153
154void fillParticle(int id, double ee, double thetaIn, double phiIn,
155 Pythia8::Event& event, Pythia8::ParticleData& pdt, Pythia8::Rndm& rndm, bool atRest = false) {
156
157 // Reset event record to allow for new event.
158 event.reset();
159
160 // Select particle mass; where relevant according to Breit-Wigner.
161 double mm = pdt.mSel(id);
162 double pp = Pythia8::sqrtpos(ee*ee - mm*mm);
163
164 // Special case when particle is supposed to be at rest.
165 if (atRest) {
166 ee = mm;
167 pp = 0.;
168 }
169
170 // Angles as input or uniform in solid angle.
171 double cThe, sThe, phi;
172 if (thetaIn >= 0.) {
173 cThe = cos(thetaIn);
174 sThe = sin(thetaIn);
175 phi = phiIn;
176 } else {
177 cThe = 2. * rndm.flat() - 1.;
178 sThe = Pythia8::sqrtpos(1. - cThe * cThe);
179 phi = 2. * M_PI * rndm.flat();
180 }
181
182 // Store the particle in the event record.
183 event.append( id, 1, 0, 0, pp * sThe * cos(phi), pp * sThe * sin(phi), pp * cThe, ee, mm);
184
185}
186
187
188
189//---------------------------------------------------------------------------
190
191int main(int argc, char *argv[])
192{
193 char appName[] = "DelphesPythia8";
194 stringstream message;
195 FILE *inputFile = 0;
196 TFile *outputFile = 0;
197 TStopwatch readStopWatch, procStopWatch;
198 ExRootTreeWriter *treeWriter = 0;
199 ExRootTreeBranch *branchEvent = 0;
200 ExRootTreeBranch *branchEventLHEF = 0, *branchWeightLHEF = 0;
201 ExRootConfReader *confReader = 0;
202 Delphes *modularDelphes = 0;
203 DelphesFactory *factory = 0;
204 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
205 TObjArray *stableParticleOutputArrayLHEF = 0, *allParticleOutputArrayLHEF = 0, *partonOutputArrayLHEF = 0;
206 DelphesLHEFReader *reader = 0;
207 Long64_t eventCounter, errorCounter;
208 Long64_t numberOfEvents, timesAllowErrors;
209
210 Pythia8::Pythia *pythia = 0;
211
212 if(argc != 4)
213 {
214 cout << " Usage: " << appName << " config_file" << " pythia_card" << " output_file" << endl;
215 cout << " config_file - configuration file in Tcl format," << endl;
216 cout << " pythia_card - Pythia8 configuration file," << endl;
217 cout << " output_file - output file in ROOT format." << endl;
218 return 1;
219 }
220
221 signal(SIGINT, SignalHandler);
222
223 gROOT->SetBatch();
224
225 int appargc = 1;
226 char *appargv[] = {appName};
227 TApplication app(appName, &appargc, appargv);
228
229 try
230 {
231 outputFile = TFile::Open(argv[3], "CREATE");
232
233 if(outputFile == NULL)
234 {
235 message << "can't create output file " << argv[3];
236 throw runtime_error(message.str());
237 }
238
239 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
240
241 branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
242
243 confReader = new ExRootConfReader;
244 confReader->ReadFile(argv[1]);
245
246 modularDelphes = new Delphes("Delphes");
247 modularDelphes->SetConfReader(confReader);
248 modularDelphes->SetTreeWriter(treeWriter);
249
250 factory = modularDelphes->GetFactory();
251 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
252 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
253 partonOutputArray = modularDelphes->ExportArray("partons");
254
255 // Initialize Pythia
256 pythia = new Pythia8::Pythia;
257 //Pythia8::Event& event = pythia->event;
258 //Pythia8::ParticleData& pdt = pythia->particleData;
259
260 if(pythia == NULL)
261 {
262 throw runtime_error("can't create Pythia instance");
263 }
264
265 // Read in commands from configuration file
266 if(!pythia->readFile(argv[2]))
267 {
268 message << "can't read Pythia8 configuration file " << argv[2] << endl;
269 throw runtime_error(message.str());
270 }
271
272 // Extract settings to be used in the main program
273 numberOfEvents = pythia->mode("Main:numberOfEvents");
274 timesAllowErrors = pythia->mode("Main:timesAllowErrors");
275
276 // Check if particle gun
277 if (!pythia->flag("Main:spareFlag1"))
278 {
279 inputFile = fopen(pythia->word("Beams:LHEF").c_str(), "r");
280 if(inputFile)
281 {
282 reader = new DelphesLHEFReader;
283 reader->SetInputFile(inputFile);
284
285 branchEventLHEF = treeWriter->NewBranch("EventLHEF", LHEFEvent::Class());
286 branchWeightLHEF = treeWriter->NewBranch("WeightLHEF", LHEFWeight::Class());
287
288 allParticleOutputArrayLHEF = modularDelphes->ExportArray("allParticlesLHEF");
289 stableParticleOutputArrayLHEF = modularDelphes->ExportArray("stableParticlesLHEF");
290 partonOutputArrayLHEF = modularDelphes->ExportArray("partonsLHEF");
291 }
292 }
293
294 modularDelphes->InitTask();
295
296 pythia->init();
297
298 // ExRootProgressBar progressBar(numberOfEvents - 1);
299 ExRootProgressBar progressBar(-1);
300
301 // Loop over all events
302 errorCounter = 0;
303 treeWriter->Clear();
304 modularDelphes->Clear();
305 readStopWatch.Start();
306 for(eventCounter = 0; eventCounter < numberOfEvents && !interrupted; ++eventCounter)
307 {
308 while(reader && reader->ReadBlock(factory, allParticleOutputArrayLHEF,
309 stableParticleOutputArrayLHEF, partonOutputArrayLHEF) && !reader->EventReady());
310
311 if (pythia->flag("Main:spareFlag1"))
312 {
313 fillParticle( pythia->mode("Main:spareMode1"), 30, -1., 0.,pythia->event, pythia->particleData, pythia->rndm, 0);
314 }
315
316 if(!pythia->next())
317 {
318 // If failure because reached end of file then exit event loop
319 if(pythia->info.atEndOfFile())
320 {
321 cerr << "Aborted since reached end of Les Houches Event File" << endl;
322 break;
323 }
324
325 // First few failures write off as "acceptable" errors, then quit
326 if(++errorCounter > timesAllowErrors)
327 {
328 cerr << "Event generation aborted prematurely, owing to error!" << endl;
329 break;
330 }
331
332 modularDelphes->Clear();
333 reader->Clear();
334 continue;
335 }
336
337 readStopWatch.Stop();
338
339 procStopWatch.Start();
340 ConvertInput(eventCounter, pythia, branchEvent, factory,
341 allParticleOutputArray, stableParticleOutputArray, partonOutputArray,
342 &readStopWatch, &procStopWatch);
343 modularDelphes->ProcessTask();
344 procStopWatch.Stop();
345
346 if(reader)
347 {
348 reader->AnalyzeEvent(branchEventLHEF, eventCounter, &readStopWatch, &procStopWatch);
349 reader->AnalyzeWeight(branchWeightLHEF);
350 }
351
352 treeWriter->Fill();
353
354 treeWriter->Clear();
355 modularDelphes->Clear();
356 if(reader) reader->Clear();
357
358 readStopWatch.Start();
359 progressBar.Update(eventCounter, eventCounter);
360 }
361
362 progressBar.Update(eventCounter, eventCounter, kTRUE);
363 progressBar.Finish();
364
365 pythia->stat();
366
367 modularDelphes->FinishTask();
368 treeWriter->Write();
369
370 cout << "** Exiting..." << endl;
371
372 delete reader;
373 delete pythia;
374 delete modularDelphes;
375 delete confReader;
376 delete treeWriter;
377 delete outputFile;
378
379 return 0;
380 }
381 catch(runtime_error &e)
382 {
383 if(treeWriter) delete treeWriter;
384 if(outputFile) delete outputFile;
385 cerr << "** ERROR: " << e.what() << endl;
386 return 1;
387 }
388}
Note: See TracBrowser for help on using the repository browser.