Fork me on GitHub

source: svn/trunk/examples/DelphesCMSFWLite.cpp@ 1106

Last change on this file since 1106 was 1051, checked in by Pavel Demin, 12 years ago
replace (fStatus == 2) with (pdgCode <= 5
pdgCode == 21 pdgCode == 15) for fPartonOutputArray
  • Property svn:eol-style set to native
  • Property svn:keywords set to Id Revision Date
File size: 5.8 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 "FWCore/FWLite/interface/AutoLibraryLoader.h"
32
33#include "DataFormats/FWLite/interface/Event.h"
34#include "DataFormats/FWLite/interface/Handle.h"
35#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
36
37using namespace std;
38
39//---------------------------------------------------------------------------
40
41void ConvertInput(fwlite::Event &event, DelphesFactory *factory, TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray)
42{
43 size_t i;
44 fwlite::Handle< vector< reco::GenParticle > > handleParticle;
45 handleParticle.getByLabel(event, "genParticles");
46
47 Candidate *candidate;
48 TDatabasePDG *pdg;
49 TParticlePDG *pdgParticle;
50 Int_t pdgCode;
51
52 Int_t pid, status;
53 Double_t px, py, pz, e;
54 Double_t x, y, z;
55
56 pdg = TDatabasePDG::Instance();
57
58 for(i = 0; i < handleParticle->size(); ++i)
59 {
60 const reco::GenParticle &particle = (*handleParticle)[i];
61
62 pid = particle.pdgId();
63 status = particle.status();
64 px = particle.px(), py = particle.py(), pz = particle.phi(), e = particle.energy();
65 x = particle.vx(), y = particle.vy(), z = particle.vz();
66
67 if(status == 1 || status == 2)
68 {
69 candidate = factory->NewCandidate();
70
71 candidate->PID = pid;
72 pdgCode = TMath::Abs(candidate->PID);
73
74 candidate->Status = status;
75
76 pdgParticle = pdg->GetParticle(pid);
77 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
78 candidate->Mass = pdgParticle ? pdgParticle->Mass() : -999.9;
79
80 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
81
82 candidate->Position.SetXYZT(x, y, z, 0.0);
83
84 allParticleOutputArray->Add(candidate);
85
86 if(!pdgParticle) return;
87
88 if(status == 1)
89 {
90 stableParticleOutputArray->Add(candidate);
91 }
92 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
93 {
94 partonOutputArray->Add(candidate);
95 }
96 }
97 }
98}
99
100//---------------------------------------------------------------------------
101
102static bool interrupted = false;
103
104void SignalHandler(int sig)
105{
106 interrupted = true;
107}
108
109//---------------------------------------------------------------------------
110
111int main(int argc, char *argv[])
112{
113 char appName[] = "DelphesCMSFWLite";
114 stringstream message;
115 TFile *inputFile = 0;
116 TFile *outputFile = 0;
117 TStopwatch eventStopWatch;
118 ExRootTreeWriter *treeWriter = 0;
119 ExRootConfReader *confReader = 0;
120 Delphes *modularDelphes = 0;
121 DelphesFactory *factory = 0;
122 TObjArray *allParticleOutputArray = 0, *stableParticleOutputArray = 0, *partonOutputArray = 0;
123 Int_t i;
124 Long64_t entry, allEntries;
125
126 if(argc < 4)
127 {
128 cout << " Usage: " << appName << " config_file" << " output_file" << " input_file(s)" << endl;
129 cout << " config_file - configuration file in Tcl format," << endl;
130 cout << " output_file - output file in ROOT format," << endl;
131 cout << " input_file(s) - input file(s) in ROOT format." << endl;
132 return 1;
133 }
134
135 signal(SIGINT, SignalHandler);
136
137 gROOT->SetBatch();
138
139 int appargc = 1;
140 char *appargv[] = {appName};
141 TApplication app(appName, &appargc, appargv);
142
143 AutoLibraryLoader::enable();
144
145 try
146 {
147 outputFile = TFile::Open(argv[2], "CREATE");
148
149 if(outputFile == NULL)
150 {
151 message << "can't open " << argv[2] << endl;
152 throw runtime_error(message.str());
153 }
154
155 treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
156
157 confReader = new ExRootConfReader;
158 confReader->ReadFile(argv[1]);
159
160 modularDelphes = new Delphes("Delphes");
161 modularDelphes->SetConfReader(confReader);
162 modularDelphes->SetTreeWriter(treeWriter);
163
164 factory = modularDelphes->GetFactory();
165 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
166 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
167 partonOutputArray = modularDelphes->ExportArray("partons");
168
169 modularDelphes->InitTask();
170
171 for(i = 3; i < argc && !interrupted; ++i)
172 {
173 cout << "** Reading " << argv[i] << endl;
174
175 inputFile = TFile::Open(argv[i]);
176
177 if(inputFile == NULL)
178 {
179 message << "can't open " << argv[i] << endl;
180 throw runtime_error(message.str());
181 }
182
183 fwlite::Event event(inputFile);
184
185 allEntries = event.size();
186
187 if(allEntries <= 0) continue;
188
189 ExRootProgressBar progressBar(allEntries - 1);
190
191 // Loop over all objects
192 entry = 0;
193 modularDelphes->Clear();
194 treeWriter->Clear();
195 for(event.toBegin(); !event.atEnd() && !interrupted; ++event)
196 {
197 ConvertInput(event, factory, allParticleOutputArray, stableParticleOutputArray, partonOutputArray);
198 modularDelphes->ProcessTask();
199
200 treeWriter->Fill();
201
202 modularDelphes->Clear();
203 treeWriter->Clear();
204
205 progressBar.Update(entry);
206 ++entry;
207 }
208 progressBar.Finish();
209
210 inputFile->Close();
211 }
212
213 modularDelphes->FinishTask();
214 treeWriter->Write();
215
216 cout << "** Exiting..." << endl;
217
218 delete modularDelphes;
219 delete confReader;
220 delete treeWriter;
221 delete outputFile;
222
223 return 0;
224 }
225 catch(runtime_error &e)
226 {
227 if(treeWriter) delete treeWriter;
228 if(outputFile) delete outputFile;
229 cerr << "** ERROR: " << e.what() << endl;
230 return 1;
231 }
232}
Note: See TracBrowser for help on using the repository browser.