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 <algorithm>
|
---|
20 | #include <stdexcept>
|
---|
21 | #include <iostream>
|
---|
22 | #include <sstream>
|
---|
23 | #include <memory>
|
---|
24 |
|
---|
25 | #include <map>
|
---|
26 | #include <vector>
|
---|
27 |
|
---|
28 | #include <stdlib.h>
|
---|
29 | #include <signal.h>
|
---|
30 | #include <stdio.h>
|
---|
31 |
|
---|
32 | #include "TROOT.h"
|
---|
33 | #include "TApplication.h"
|
---|
34 |
|
---|
35 | #include "TFile.h"
|
---|
36 | #include "TObjArray.h"
|
---|
37 | #include "TStopwatch.h"
|
---|
38 | #include "TDatabasePDG.h"
|
---|
39 | #include "TParticlePDG.h"
|
---|
40 | #include "TLorentzVector.h"
|
---|
41 |
|
---|
42 | #include "modules/Delphes.h"
|
---|
43 | #include "classes/DelphesStream.h"
|
---|
44 | #include "classes/DelphesClasses.h"
|
---|
45 | #include "classes/DelphesFactory.h"
|
---|
46 |
|
---|
47 | #include "ExRootAnalysis/ExRootTreeWriter.h"
|
---|
48 | #include "ExRootAnalysis/ExRootTreeBranch.h"
|
---|
49 | #include "ExRootAnalysis/ExRootProgressBar.h"
|
---|
50 |
|
---|
51 | #include "FWCore/FWLite/interface/FWLiteEnabler.h"
|
---|
52 | #include "DataFormats/FWLite/interface/Event.h"
|
---|
53 | #include "DataFormats/FWLite/interface/Handle.h"
|
---|
54 | #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
|
---|
55 | #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
|
---|
56 | #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
|
---|
57 | #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
|
---|
58 | #include "SimDataFormats/GeneratorProducts/interface/WeightsInfo.h"
|
---|
59 |
|
---|
60 | using namespace std;
|
---|
61 |
|
---|
62 | //---------------------------------------------------------------------------
|
---|
63 |
|
---|
64 | void ConvertInput(fwlite::Event &event, Long64_t eventCounter,
|
---|
65 | ExRootTreeBranch *branchEvent, ExRootTreeBranch *branchRwgt,
|
---|
66 | DelphesFactory *factory, TObjArray *allParticleOutputArray,
|
---|
67 | TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray)
|
---|
68 | {
|
---|
69 | fwlite::Handle< GenEventInfoProduct > handleGenEventInfo;
|
---|
70 |
|
---|
71 | fwlite::Handle< LHEEventProduct > handleLHEEvent;
|
---|
72 |
|
---|
73 | fwlite::Handle< vector< reco::GenParticle > > handleParticle;
|
---|
74 | vector< reco::GenParticle >::const_iterator itParticle;
|
---|
75 |
|
---|
76 | vector< const reco::Candidate * > vectorCandidate;
|
---|
77 | vector< const reco::Candidate * >::iterator itCandidate;
|
---|
78 |
|
---|
79 | handleGenEventInfo.getByLabel(event, "generator");
|
---|
80 | handleLHEEvent.getByLabel(event, "externalLHEProducer");
|
---|
81 | handleParticle.getByLabel(event, "genParticles");
|
---|
82 |
|
---|
83 | HepMCEvent *element;
|
---|
84 | Weight *weight;
|
---|
85 | Candidate *candidate;
|
---|
86 | TDatabasePDG *pdg;
|
---|
87 | TParticlePDG *pdgParticle;
|
---|
88 | Int_t pdgCode;
|
---|
89 |
|
---|
90 | Int_t pid, status;
|
---|
91 | Double_t px, py, pz, e, mass;
|
---|
92 | Double_t x, y, z;
|
---|
93 |
|
---|
94 | const vector< gen::WeightsInfo > &vectorWeightsInfo = handleLHEEvent->weights();
|
---|
95 | vector< gen::WeightsInfo >::const_iterator itWeightsInfo;
|
---|
96 |
|
---|
97 | element = static_cast<HepMCEvent *>(branchEvent->NewEntry());
|
---|
98 |
|
---|
99 | element->Number = eventCounter;
|
---|
100 |
|
---|
101 | element->ProcessID = handleGenEventInfo->signalProcessID();
|
---|
102 | element->MPI = 1;
|
---|
103 | element->Weight = handleGenEventInfo->weight();
|
---|
104 | element->Scale = handleGenEventInfo->qScale();
|
---|
105 | element->AlphaQED = handleGenEventInfo->alphaQED();
|
---|
106 | element->AlphaQCD = handleGenEventInfo->alphaQCD();
|
---|
107 |
|
---|
108 | element->ID1 = 0;
|
---|
109 | element->ID2 = 0;
|
---|
110 | element->X1 = 0.0;
|
---|
111 | element->X2 = 0.0;
|
---|
112 | element->ScalePDF = 0.0;
|
---|
113 | element->PDF1 = 0.0;
|
---|
114 | element->PDF2 = 0.0;
|
---|
115 |
|
---|
116 | element->ReadTime = 0.0;
|
---|
117 | element->ProcTime = 0.0;
|
---|
118 |
|
---|
119 | for(itWeightsInfo = vectorWeightsInfo.begin(); itWeightsInfo != vectorWeightsInfo.end(); ++itWeightsInfo)
|
---|
120 | {
|
---|
121 | weight = static_cast<Weight *>(branchRwgt->NewEntry());
|
---|
122 | weight->Weight = itWeightsInfo->wgt;
|
---|
123 | }
|
---|
124 |
|
---|
125 | pdg = TDatabasePDG::Instance();
|
---|
126 |
|
---|
127 | for(itParticle = handleParticle->begin(); itParticle != handleParticle->end(); ++itParticle)
|
---|
128 | {
|
---|
129 | vectorCandidate.push_back(&*itParticle);
|
---|
130 | }
|
---|
131 |
|
---|
132 | for(itParticle = handleParticle->begin(); itParticle != handleParticle->end(); ++itParticle)
|
---|
133 | {
|
---|
134 | const reco::GenParticle &particle = *itParticle;
|
---|
135 |
|
---|
136 | pid = particle.pdgId();
|
---|
137 | status = particle.status();
|
---|
138 | px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.energy(); mass = particle.mass();
|
---|
139 | x = particle.vx(); y = particle.vy(); z = particle.vz();
|
---|
140 |
|
---|
141 | candidate = factory->NewCandidate();
|
---|
142 |
|
---|
143 | candidate->PID = pid;
|
---|
144 | pdgCode = TMath::Abs(candidate->PID);
|
---|
145 |
|
---|
146 | candidate->Status = status;
|
---|
147 |
|
---|
148 | if(particle.mother())
|
---|
149 | {
|
---|
150 | itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.mother());
|
---|
151 | if(itCandidate != vectorCandidate.end()) candidate->M1 = distance(vectorCandidate.begin(), itCandidate);
|
---|
152 | }
|
---|
153 |
|
---|
154 | itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(0));
|
---|
155 | if(itCandidate != vectorCandidate.end()) candidate->D1 = distance(vectorCandidate.begin(), itCandidate);
|
---|
156 |
|
---|
157 | itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(particle.numberOfDaughters() - 1));
|
---|
158 | if(itCandidate != vectorCandidate.end()) candidate->D2 = distance(vectorCandidate.begin(), itCandidate);
|
---|
159 |
|
---|
160 | pdgParticle = pdg->GetParticle(pid);
|
---|
161 | candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
|
---|
162 | candidate->Mass = mass;
|
---|
163 |
|
---|
164 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
165 |
|
---|
166 | candidate->Position.SetXYZT(x*10.0, y*10.0, z*10.0, 0.0);
|
---|
167 |
|
---|
168 | allParticleOutputArray->Add(candidate);
|
---|
169 |
|
---|
170 | if(!pdgParticle) continue;
|
---|
171 |
|
---|
172 | if(status == 1)
|
---|
173 | {
|
---|
174 | stableParticleOutputArray->Add(candidate);
|
---|
175 | }
|
---|
176 | else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
|
---|
177 | {
|
---|
178 | partonOutputArray->Add(candidate);
|
---|
179 | }
|
---|
180 | }
|
---|
181 | }
|
---|
182 |
|
---|
183 | //---------------------------------------------------------------------------
|
---|
184 |
|
---|
185 | static bool interrupted = false;
|
---|
186 |
|
---|
187 | void SignalHandler(int sig)
|
---|
188 | {
|
---|
189 | interrupted = true;
|
---|
190 | }
|
---|
191 |
|
---|
192 | //---------------------------------------------------------------------------
|
---|
193 |
|
---|
194 | int main(int argc, char *argv[])
|
---|
195 | {
|
---|
196 | char appName[] = "DelphesCMSFWLite";
|
---|
197 | stringstream message;
|
---|
198 | TFile *inputFile = 0;
|
---|
199 | TFile *outputFile = 0;
|
---|
200 | TStopwatch eventStopWatch;
|
---|
201 | ExRootTreeWriter *treeWriter = 0;
|
---|
202 | ExRootTreeBranch *branchEvent = 0, *branchRwgt = 0;
|
---|
203 | ExRootConfReader *confReader = 0;
|
---|
204 | Delphes *modularDelphes = 0;
|
---|
205 | DelphesFactory *factory = 0;
|
---|
206 | TObjArray *allParticleOutputArray = 0, *stableParticleOutputArray = 0, *partonOutputArray = 0;
|
---|
207 | Int_t i;
|
---|
208 | Long64_t eventCounter, numberOfEvents;
|
---|
209 |
|
---|
210 | if(argc < 4)
|
---|
211 | {
|
---|
212 | cout << " Usage: " << appName << " config_file" << " output_file" << " input_file(s)" << endl;
|
---|
213 | cout << " config_file - configuration file in Tcl format," << endl;
|
---|
214 | cout << " output_file - output file in ROOT format," << endl;
|
---|
215 | cout << " input_file(s) - input file(s) in ROOT format." << endl;
|
---|
216 | return 1;
|
---|
217 | }
|
---|
218 |
|
---|
219 | signal(SIGINT, SignalHandler);
|
---|
220 |
|
---|
221 | gROOT->SetBatch();
|
---|
222 |
|
---|
223 | int appargc = 1;
|
---|
224 | char *appargv[] = {appName};
|
---|
225 | TApplication app(appName, &appargc, appargv);
|
---|
226 |
|
---|
227 | FWLiteEnabler::enable();
|
---|
228 |
|
---|
229 | try
|
---|
230 | {
|
---|
231 | outputFile = TFile::Open(argv[2], "CREATE");
|
---|
232 |
|
---|
233 | if(outputFile == NULL)
|
---|
234 | {
|
---|
235 | message << "can't open " << argv[2] << endl;
|
---|
236 | throw runtime_error(message.str());
|
---|
237 | }
|
---|
238 |
|
---|
239 | treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
|
---|
240 |
|
---|
241 | branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
|
---|
242 | branchRwgt = treeWriter->NewBranch("Rwgt", Weight::Class());
|
---|
243 |
|
---|
244 | confReader = new ExRootConfReader;
|
---|
245 | confReader->ReadFile(argv[1]);
|
---|
246 |
|
---|
247 | modularDelphes = new Delphes("Delphes");
|
---|
248 | modularDelphes->SetConfReader(confReader);
|
---|
249 | modularDelphes->SetTreeWriter(treeWriter);
|
---|
250 |
|
---|
251 | factory = modularDelphes->GetFactory();
|
---|
252 | allParticleOutputArray = modularDelphes->ExportArray("allParticles");
|
---|
253 | stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
|
---|
254 | partonOutputArray = modularDelphes->ExportArray("partons");
|
---|
255 |
|
---|
256 | modularDelphes->InitTask();
|
---|
257 |
|
---|
258 | for(i = 3; i < argc && !interrupted; ++i)
|
---|
259 | {
|
---|
260 | cout << "** Reading " << argv[i] << endl;
|
---|
261 |
|
---|
262 | inputFile = TFile::Open(argv[i]);
|
---|
263 |
|
---|
264 | if(inputFile == NULL)
|
---|
265 | {
|
---|
266 | message << "can't open " << argv[i] << endl;
|
---|
267 | throw runtime_error(message.str());
|
---|
268 | }
|
---|
269 |
|
---|
270 | fwlite::Event event(inputFile);
|
---|
271 |
|
---|
272 | numberOfEvents = event.size();
|
---|
273 |
|
---|
274 | if(numberOfEvents <= 0) continue;
|
---|
275 |
|
---|
276 | // ExRootProgressBar progressBar(numberOfEvents - 1);
|
---|
277 | ExRootProgressBar progressBar(-1);
|
---|
278 |
|
---|
279 | // Loop over all objects
|
---|
280 | eventCounter = 0;
|
---|
281 | modularDelphes->Clear();
|
---|
282 | treeWriter->Clear();
|
---|
283 | for(event.toBegin(); !event.atEnd() && !interrupted; ++event)
|
---|
284 | {
|
---|
285 | ConvertInput(event, eventCounter, branchEvent, branchRwgt, factory,
|
---|
286 | allParticleOutputArray, stableParticleOutputArray, partonOutputArray);
|
---|
287 | modularDelphes->ProcessTask();
|
---|
288 |
|
---|
289 | treeWriter->Fill();
|
---|
290 |
|
---|
291 | modularDelphes->Clear();
|
---|
292 | treeWriter->Clear();
|
---|
293 |
|
---|
294 | progressBar.Update(eventCounter, eventCounter);
|
---|
295 | ++eventCounter;
|
---|
296 | }
|
---|
297 |
|
---|
298 | progressBar.Update(eventCounter, eventCounter, kTRUE);
|
---|
299 | progressBar.Finish();
|
---|
300 |
|
---|
301 | inputFile->Close();
|
---|
302 | }
|
---|
303 |
|
---|
304 | modularDelphes->FinishTask();
|
---|
305 | treeWriter->Write();
|
---|
306 |
|
---|
307 | cout << "** Exiting..." << endl;
|
---|
308 |
|
---|
309 | delete modularDelphes;
|
---|
310 | delete confReader;
|
---|
311 | delete treeWriter;
|
---|
312 | delete outputFile;
|
---|
313 |
|
---|
314 | return 0;
|
---|
315 | }
|
---|
316 | catch(runtime_error &e)
|
---|
317 | {
|
---|
318 | if(treeWriter) delete treeWriter;
|
---|
319 | if(outputFile) delete outputFile;
|
---|
320 | cerr << "** ERROR: " << e.what() << endl;
|
---|
321 | return 1;
|
---|
322 | }
|
---|
323 | }
|
---|