Fork me on GitHub

source: git/readers/DelphesPythia8.cpp@ 239e1d0

ImprovedOutputFile Timing dual_readout llp 3.3.3pre11
Last change on this file since 239e1d0 was 2acc2e9, checked in by Alexandre Mertens <alexandre.mertens@…>, 9 years ago

working for di-jet and single e, mu, gamma

  • Property mode set to 100644
File size: 12.1 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[5150e39]19#include <stdexcept>
20#include <iostream>
21#include <sstream>
[cd616b9]22#include <string>
[5150e39]23
24#include <signal.h>
25
[280667b]26#include "Pythia.h"
[5150e39]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"
[cd616b9]41#include "classes/DelphesLHEFReader.h"
[5150e39]42
43#include "ExRootAnalysis/ExRootTreeWriter.h"
44#include "ExRootAnalysis/ExRootTreeBranch.h"
45#include "ExRootAnalysis/ExRootProgressBar.h"
46
47using namespace std;
48
49//---------------------------------------------------------------------------
50
[b94aacf]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)
[5150e39]55{
56 int i;
57
[ad78364]58 HepMCEvent *element;
[5150e39]59 Candidate *candidate;
60 TDatabasePDG *pdg;
61 TParticlePDG *pdgParticle;
62 Int_t pdgCode;
63
64 Int_t pid, status;
[80d4a34]65 Double_t px, py, pz, e, mass;
[474eb76]66 Double_t x, y, z, t;
[5150e39]67
[b94aacf]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
[ad78364]80 element->ID1 = pythia->info.id1();
81 element->ID2 = pythia->info.id2();
82 element->X1 = pythia->info.x1();
83 element->X2 = pythia->info.x2();
[b94aacf]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
[5150e39]91 pdg = TDatabasePDG::Instance();
92
[d8b1858]93 for(i = 1; i < pythia->event.size(); ++i)
[5150e39]94 {
95 Pythia8::Particle &particle = pythia->event[i];
96
97 pid = particle.id();
[3e276d4]98 status = particle.statusHepMC();
[80d4a34]99 px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.e(); mass = particle.m();
[5150e39]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
[474eb76]109 candidate->M1 = particle.mother1() - 1;
110 candidate->M2 = particle.mother2() - 1;
[5150e39]111
[474eb76]112 candidate->D1 = particle.daughter1() - 1;
113 candidate->D2 = particle.daughter2() - 1;
[5150e39]114
115 pdgParticle = pdg->GetParticle(pid);
116 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
[bba4646]117 candidate->Mass = mass;
[5150e39]118
119 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
120
[474eb76]121 candidate->Position.SetXYZT(x, y, z, t);
[5150e39]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
[971db5b]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
[2acc2e9]154void fillParticle(int id, double ee_max, double thetaIn, double phiIn,
[971db5b]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
[2acc2e9]160 // Angles uniform in solid angle.
161 double cThe, sThe, phi, ee;
162 cThe = 2. * rndm.flat() - 1.;
163 sThe = Pythia8::sqrtpos(1. - cThe * cThe);
164 phi = 2. * M_PI * rndm.flat();
165 ee = pow(10,1+(log10(ee_max)-1)*rndm.flat());
[971db5b]166 double mm = pdt.mSel(id);
167 double pp = Pythia8::sqrtpos(ee*ee - mm*mm);
168
169 // Store the particle in the event record.
170 event.append( id, 1, 0, 0, pp * sThe * cos(phi), pp * sThe * sin(phi), pp * cThe, ee, mm);
171
172}
173
[2acc2e9]174void fillPartons(int type, double ee_max, Pythia8::Event& event, Pythia8::ParticleData& pdt,
175 Pythia8::Rndm& rndm) {
176
177 // Reset event record to allow for new event.
178 event.reset();
179
180 // Angles uniform in solid angle.
181 double cThe, sThe, phi, ee;
182
183 // Information on a q qbar system, to be hadronized.
184
185 cThe = 2. * rndm.flat() - 1.;
186 sThe = Pythia8::sqrtpos(1. - cThe * cThe);
187 phi = 2. * M_PI * rndm.flat();
188 ee = pow(10,1+(log10(ee_max)-1)*rndm.flat());
189 double mm = pdt.m0(type);
190 double pp = Pythia8::sqrtpos(ee*ee - mm*mm);
191 if (type == 21)
192 {
193 event.append( 21, 23, 101, 102, pp * sThe * cos(phi), pp * sThe * sin(phi), pp * cThe, ee);
194 event.append( 21, 23, 102, 101, -pp * sThe * cos(phi), -pp * sThe * sin(phi), -pp * cThe, ee);
195 }
196 else
197 {
198 event.append( type, 23, 101, 0, pp * sThe * cos(phi), pp * sThe * sin(phi), pp * cThe, ee, mm);
199 event.append( -type, 23, 0, 101, -pp * sThe * cos(phi), -pp * sThe * sin(phi), -pp * cThe, ee, mm);
200 }
201}
[971db5b]202
203
[5150e39]204//---------------------------------------------------------------------------
205
206int main(int argc, char *argv[])
207{
208 char appName[] = "DelphesPythia8";
209 stringstream message;
[cd616b9]210 FILE *inputFile = 0;
[5150e39]211 TFile *outputFile = 0;
212 TStopwatch readStopWatch, procStopWatch;
213 ExRootTreeWriter *treeWriter = 0;
214 ExRootTreeBranch *branchEvent = 0;
[1e8afcc]215 ExRootTreeBranch *branchEventLHEF = 0, *branchWeightLHEF = 0;
[5150e39]216 ExRootConfReader *confReader = 0;
217 Delphes *modularDelphes = 0;
218 DelphesFactory *factory = 0;
219 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
[cd616b9]220 TObjArray *stableParticleOutputArrayLHEF = 0, *allParticleOutputArrayLHEF = 0, *partonOutputArrayLHEF = 0;
221 DelphesLHEFReader *reader = 0;
[5150e39]222 Long64_t eventCounter, errorCounter;
223 Long64_t numberOfEvents, timesAllowErrors;
224
225 Pythia8::Pythia *pythia = 0;
226
227 if(argc != 4)
228 {
229 cout << " Usage: " << appName << " config_file" << " pythia_card" << " output_file" << endl;
230 cout << " config_file - configuration file in Tcl format," << endl;
231 cout << " pythia_card - Pythia8 configuration file," << endl;
232 cout << " output_file - output file in ROOT format." << endl;
233 return 1;
234 }
235
236 signal(SIGINT, SignalHandler);
237
238 gROOT->SetBatch();
239
240 int appargc = 1;
241 char *appargv[] = {appName};
242 TApplication app(appName, &appargc, appargv);
243
244 try
245 {
246 outputFile = TFile::Open(argv[3], "CREATE");
247
248 if(outputFile == NULL)
249 {
250 message << "can't create output file " << argv[3];
251 throw runtime_error(message.str());
252 }
253
254 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
255
[b94aacf]256 branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
[5150e39]257
258 confReader = new ExRootConfReader;
259 confReader->ReadFile(argv[1]);
260
261 modularDelphes = new Delphes("Delphes");
262 modularDelphes->SetConfReader(confReader);
263 modularDelphes->SetTreeWriter(treeWriter);
264
265 factory = modularDelphes->GetFactory();
266 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
267 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
268 partonOutputArray = modularDelphes->ExportArray("partons");
269
[cd616b9]270 // Initialize Pythia
[5150e39]271 pythia = new Pythia8::Pythia;
[971db5b]272 //Pythia8::Event& event = pythia->event;
273 //Pythia8::ParticleData& pdt = pythia->particleData;
[5150e39]274
275 if(pythia == NULL)
276 {
277 throw runtime_error("can't create Pythia instance");
278 }
279
280 // Read in commands from configuration file
[c0f9b5b]281 if(!pythia->readFile(argv[2]))
282 {
283 message << "can't read Pythia8 configuration file " << argv[2] << endl;
284 throw runtime_error(message.str());
285 }
[5150e39]286
287 // Extract settings to be used in the main program
288 numberOfEvents = pythia->mode("Main:numberOfEvents");
289 timesAllowErrors = pythia->mode("Main:timesAllowErrors");
290
[971db5b]291 // Check if particle gun
292 if (!pythia->flag("Main:spareFlag1"))
[cd616b9]293 {
[971db5b]294 inputFile = fopen(pythia->word("Beams:LHEF").c_str(), "r");
295 if(inputFile)
296 {
297 reader = new DelphesLHEFReader;
298 reader->SetInputFile(inputFile);
[cd616b9]299
[971db5b]300 branchEventLHEF = treeWriter->NewBranch("EventLHEF", LHEFEvent::Class());
301 branchWeightLHEF = treeWriter->NewBranch("WeightLHEF", LHEFWeight::Class());
[cd616b9]302
[971db5b]303 allParticleOutputArrayLHEF = modularDelphes->ExportArray("allParticlesLHEF");
304 stableParticleOutputArrayLHEF = modularDelphes->ExportArray("stableParticlesLHEF");
305 partonOutputArrayLHEF = modularDelphes->ExportArray("partonsLHEF");
306 }
[cd616b9]307 }
308
[a446115]309 modularDelphes->InitTask();
310
[5150e39]311 pythia->init();
312
[a0538b9]313 // ExRootProgressBar progressBar(numberOfEvents - 1);
314 ExRootProgressBar progressBar(-1);
[5150e39]315
316 // Loop over all events
317 errorCounter = 0;
318 treeWriter->Clear();
319 modularDelphes->Clear();
320 readStopWatch.Start();
321 for(eventCounter = 0; eventCounter < numberOfEvents && !interrupted; ++eventCounter)
322 {
[cd616b9]323 while(reader && reader->ReadBlock(factory, allParticleOutputArrayLHEF,
324 stableParticleOutputArrayLHEF, partonOutputArrayLHEF) && !reader->EventReady());
325
[971db5b]326 if (pythia->flag("Main:spareFlag1"))
327 {
[2acc2e9]328 if (pythia->mode("Main:spareMode1") == 11 || pythia->mode("Main:spareMode1") == 13 || pythia->mode("Main:spareMode1") == 22)
329 {
330 fillParticle( pythia->mode("Main:spareMode1"), pythia->parm("Main:spareParm1"), -1., 0.,pythia->event, pythia->particleData, pythia->rndm, 0);
331 }
332 else fillPartons( pythia->mode("Main:spareMode1"), pythia->parm("Main:spareParm1"), pythia->event, pythia->particleData, pythia->rndm);
[971db5b]333 }
334
[5150e39]335 if(!pythia->next())
336 {
337 // If failure because reached end of file then exit event loop
[cd616b9]338 if(pythia->info.atEndOfFile())
[5150e39]339 {
340 cerr << "Aborted since reached end of Les Houches Event File" << endl;
341 break;
342 }
343
344 // First few failures write off as "acceptable" errors, then quit
[cd616b9]345 if(++errorCounter > timesAllowErrors)
346 {
347 cerr << "Event generation aborted prematurely, owing to error!" << endl;
348 break;
349 }
350
351 modularDelphes->Clear();
352 reader->Clear();
353 continue;
[5150e39]354 }
355
356 readStopWatch.Stop();
357
358 procStopWatch.Start();
[b94aacf]359 ConvertInput(eventCounter, pythia, branchEvent, factory,
360 allParticleOutputArray, stableParticleOutputArray, partonOutputArray,
361 &readStopWatch, &procStopWatch);
[5150e39]362 modularDelphes->ProcessTask();
363 procStopWatch.Stop();
364
[cd616b9]365 if(reader)
366 {
367 reader->AnalyzeEvent(branchEventLHEF, eventCounter, &readStopWatch, &procStopWatch);
[1e8afcc]368 reader->AnalyzeWeight(branchWeightLHEF);
[cd616b9]369 }
370
[5150e39]371 treeWriter->Fill();
372
373 treeWriter->Clear();
374 modularDelphes->Clear();
[cd616b9]375 if(reader) reader->Clear();
[5150e39]376
377 readStopWatch.Start();
[a0538b9]378 progressBar.Update(eventCounter, eventCounter);
[5150e39]379 }
380
[a0538b9]381 progressBar.Update(eventCounter, eventCounter, kTRUE);
[5150e39]382 progressBar.Finish();
383
[3e276d4]384 pythia->stat();
[5150e39]385
386 modularDelphes->FinishTask();
387 treeWriter->Write();
388
389 cout << "** Exiting..." << endl;
390
[cd616b9]391 delete reader;
[5150e39]392 delete pythia;
393 delete modularDelphes;
394 delete confReader;
395 delete treeWriter;
396 delete outputFile;
397
398 return 0;
399 }
400 catch(runtime_error &e)
401 {
402 if(treeWriter) delete treeWriter;
403 if(outputFile) delete outputFile;
404 cerr << "** ERROR: " << e.what() << endl;
405 return 1;
406 }
407}
Note: See TracBrowser for help on using the repository browser.