Fork me on GitHub

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

Last change on this file since 1136 was 1134, checked in by Pavel Demin, 11 years ago

fix particles status selection in DelphesCMSFWLite.cpp

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