Fork me on GitHub

source: git/readers/DelphesProMC.cpp@ d352618

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

fix appName in DelphesProMC.cpp

  • Property mode set to 100644
File size: 7.2 KB
Line 
1#include <stdexcept>
2#include <iostream>
3#include <sstream>
4#include <memory>
5
6#include <map>
7
8#include <stdlib.h>
9#include <signal.h>
10#include <stdio.h>
11
12#include "TROOT.h"
13#include "TApplication.h"
14
15#include "TFile.h"
16#include "TObjArray.h"
17#include "TStopwatch.h"
18#include "TDatabasePDG.h"
19#include "TParticlePDG.h"
20#include "TLorentzVector.h"
21
22#include "modules/Delphes.h"
23#include "classes/DelphesStream.h"
24#include "classes/DelphesClasses.h"
25#include "classes/DelphesFactory.h"
26
27#include "ExRootAnalysis/ExRootTreeWriter.h"
28#include "ExRootAnalysis/ExRootTreeBranch.h"
29#include "ExRootAnalysis/ExRootProgressBar.h"
30
31#include "ProMC/ProMC.pb.h"
32#include "ProMC/ProMCBook.h"
33#include "ProMC/ProMCHeader.pb.h"
34
35using namespace std;
36
37//---------------------------------------------------------------------------
38
39void ConvertInput(ProMCEvent &event, ExRootTreeBranch *branch, DelphesFactory *factory,
40 TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray,
41 TStopwatch *readStopWatch, TStopwatch *procStopWatch)
42{
43 Int_t i;
44
45 ProMCEvent_Event *mutableEvent;
46 ProMCEvent_Particles *mutableParticles;
47
48 HepMCEvent *element;
49 Candidate *candidate;
50 TDatabasePDG *pdg;
51 TParticlePDG *pdgParticle;
52 Int_t pdgCode;
53
54 Int_t pid, status;
55 Double_t px, py, pz, mass;
56 Double_t x, y, z, t;
57
58 pdg = TDatabasePDG::Instance();
59
60 // event information
61 mutableEvent = event.mutable_event();
62
63 element = static_cast<HepMCEvent *>(branch->NewEntry());
64
65 element->Number = mutableEvent->number();
66
67 element->ProcessID = mutableEvent->process_id();
68 element->MPI = mutableEvent->mpi();
69 element->Weight = mutableEvent->weight();
70 element->Scale = mutableEvent->scale();
71 element->AlphaQED = mutableEvent->alpha_qed();
72 element->AlphaQCD = mutableEvent->alpha_qcd();
73
74 element->ID1 = mutableEvent->id1();
75 element->ID2 = mutableEvent->id2();
76 element->X1 = mutableEvent->x1();
77 element->X2 = mutableEvent->x2();
78 element->ScalePDF = mutableEvent->scale_pdf();
79 element->PDF1 = mutableEvent->pdf1();
80 element->PDF2 = mutableEvent->pdf2();
81
82 element->ReadTime = readStopWatch->RealTime();
83 element->ProcTime = procStopWatch->RealTime();
84
85 mutableParticles = event.mutable_particles();
86
87 for(i = 0; i < mutableParticles->pdg_id_size(); ++i)
88 {
89 pid = mutableParticles->pdg_id(i);
90 status = mutableParticles->status(i);
91 px = mutableParticles->px(i); py = mutableParticles->py(i); pz = mutableParticles->pz(i); mass = mutableParticles->mass(i);
92 x = mutableParticles->x(i); y = mutableParticles->y(i); z = mutableParticles->z(i); t = mutableParticles->t(i);
93
94 candidate = factory->NewCandidate();
95
96 candidate->PID = pid;
97 pdgCode = TMath::Abs(candidate->PID);
98
99 candidate->Status = status;
100
101 candidate->M1 = mutableParticles->mother1(i);
102 candidate->M2 = mutableParticles->mother2(i);
103
104 candidate->D1 = mutableParticles->daughter1(i);
105 candidate->D2 = mutableParticles->daughter2(i);
106
107 pdgParticle = pdg->GetParticle(pid);
108 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
109 candidate->Mass = mass;
110
111 candidate->Momentum.SetXYZM(px, py, pz, mass);
112
113 candidate->Position.SetXYZT(x, y, z, t);
114
115 allParticleOutputArray->Add(candidate);
116
117 if(!pdgParticle) continue;
118
119 if(status == 1)
120 {
121 stableParticleOutputArray->Add(candidate);
122 }
123 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
124 {
125 partonOutputArray->Add(candidate);
126 }
127 }
128}
129
130//---------------------------------------------------------------------------
131
132static bool interrupted = false;
133
134void SignalHandler(int sig)
135{
136 interrupted = true;
137}
138
139//---------------------------------------------------------------------------
140
141int main(int argc, char *argv[])
142{
143 char appName[] = "DelphesProMC";
144 stringstream message;
145 ProMCBook *inputFile = 0;
146 TFile *outputFile = 0;
147 TStopwatch readStopWatch, procStopWatch;
148 ExRootTreeWriter *treeWriter = 0;
149 ExRootTreeBranch *branchEvent = 0;
150 ExRootConfReader *confReader = 0;
151 Delphes *modularDelphes = 0;
152 DelphesFactory *factory = 0;
153 TObjArray *allParticleOutputArray = 0, *stableParticleOutputArray = 0, *partonOutputArray = 0;
154 Int_t i;
155 Long64_t eventCounter, numberOfEvents;
156
157 if(argc < 4)
158 {
159 cout << " Usage: " << appName << " config_file" << " output_file" << " input_file(s)" << endl;
160 cout << " config_file - configuration file in Tcl format," << endl;
161 cout << " output_file - output file in ROOT format," << endl;
162 cout << " input_file(s) - input file(s) in ProMC format." << endl;
163 return 1;
164 }
165
166 signal(SIGINT, SignalHandler);
167
168 gROOT->SetBatch();
169
170 int appargc = 1;
171 char *appargv[] = {appName};
172 TApplication app(appName, &appargc, appargv);
173
174 try
175 {
176 outputFile = TFile::Open(argv[2], "CREATE");
177
178 if(outputFile == NULL)
179 {
180 message << "can't open " << argv[2] << endl;
181 throw runtime_error(message.str());
182 }
183
184 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
185
186 branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
187
188 confReader = new ExRootConfReader;
189 confReader->ReadFile(argv[1]);
190
191 modularDelphes = new Delphes("Delphes");
192 modularDelphes->SetConfReader(confReader);
193 modularDelphes->SetTreeWriter(treeWriter);
194
195 factory = modularDelphes->GetFactory();
196 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
197 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
198 partonOutputArray = modularDelphes->ExportArray("partons");
199
200 modularDelphes->InitTask();
201
202 for(i = 3; i < argc && !interrupted; ++i)
203 {
204 cout << "** Reading " << argv[i] << endl;
205
206 inputFile = new ProMCBook(argv[i], "r");
207
208 if(inputFile == NULL)
209 {
210 message << "can't open " << argv[i] << endl;
211 throw runtime_error(message.str());
212 }
213
214 numberOfEvents = inputFile->getEvents();
215
216 if(numberOfEvents <= 0) continue;
217
218 ExRootProgressBar progressBar(numberOfEvents - 1);
219
220 // Loop over all objects
221 modularDelphes->Clear();
222 treeWriter->Clear();
223 readStopWatch.Start();
224 for(eventCounter = 0; eventCounter < numberOfEvents && !interrupted; ++eventCounter)
225 {
226 if(inputFile->next() != 0) continue;
227 ProMCEvent event = inputFile->get();
228
229 readStopWatch.Stop();
230
231 procStopWatch.Start();
232 ConvertInput(event, branchEvent, factory,
233 allParticleOutputArray, stableParticleOutputArray, partonOutputArray,
234 &readStopWatch, &procStopWatch);
235 modularDelphes->ProcessTask();
236 procStopWatch.Stop();
237
238 treeWriter->Fill();
239
240 modularDelphes->Clear();
241 treeWriter->Clear();
242
243 readStopWatch.Start();
244 progressBar.Update(eventCounter);
245 }
246
247 progressBar.Update(eventCounter, eventCounter, kTRUE);
248 progressBar.Finish();
249
250 inputFile->close();
251 delete inputFile;
252 }
253
254 modularDelphes->FinishTask();
255 treeWriter->Write();
256
257 cout << "** Exiting..." << endl;
258
259 delete modularDelphes;
260 delete confReader;
261 delete treeWriter;
262 delete outputFile;
263
264 return 0;
265 }
266 catch(runtime_error &e)
267 {
268 if(treeWriter) delete treeWriter;
269 if(outputFile) delete outputFile;
270 cerr << "** ERROR: " << e.what() << endl;
271 return 1;
272 }
273}
Note: See TracBrowser for help on using the repository browser.