Fork me on GitHub

source: git/examples/ExternalFastJetHepMC.cpp@ a8782e8

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a8782e8 was a8782e8, checked in by Michele Selvaggi <michele.selvaggi@…>, 10 years ago

replace examples by cards

  • Property mode set to 100644
File size: 8.0 KB
Line 
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/*
20
21First, you need to set DELPHES_DIR, FASTJET_DIR and ROOT_DIR:
22
23DELPHES_DIR=<path to Delphes installation>
24FASTJET_DIR=<path to FastJet installation>
25ROOT_DIR=<path to ROOT installation>
26
27Then run the following commands to build ExternalFastJetHepMC:
28
29DELPHES_LIB="-L$DELPHES_DIR -lDelphesNoFastJet"
30
31FASTJET_INC=`$FASTJET_DIR/bin/fastjet-config --cxxflags`
32FASTJET_LIB=`$FASTJET_DIR/bin/fastjet-config --libs`
33
34ROOT_INC=`$ROOT_DIR/bin/root-config --incdir`
35ROOT_LIB=`$ROOT_DIR/bin/root-config --libs`
36
37CXXFLAGS="$FASTJET_INC -I$ROOT_INC -I$DELPHES_DIR -I$DELPHES_DIR/external"
38LDFLAGS="$FASTJET_LIB $ROOT_LIB $DELPHES_LIB"
39
40g++ $CXXFLAGS $LDFLAGS ExternalFastJetHepMC.cpp -o ExternalFastJetHepMC
41
42*/
43
44#include <stdexcept>
45#include <iostream>
46#include <sstream>
47
48#include <vector>
49
50#include <signal.h>
51
52#include "TROOT.h"
53#include "TApplication.h"
54
55#include "TFile.h"
56#include "TObjArray.h"
57#include "TDatabasePDG.h"
58#include "TParticlePDG.h"
59#include "TLorentzVector.h"
60
61#include "modules/Delphes.h"
62#include "classes/DelphesClasses.h"
63#include "classes/DelphesFactory.h"
64#include "classes/DelphesHepMCReader.h"
65
66#include "ExRootAnalysis/ExRootTreeWriter.h"
67#include "ExRootAnalysis/ExRootTreeBranch.h"
68#include "ExRootAnalysis/ExRootProgressBar.h"
69
70#include "fastjet/PseudoJet.hh"
71#include "fastjet/JetDefinition.hh"
72#include "fastjet/ClusterSequence.hh"
73
74// #include "fastjet/contrib/Nsubjettiness.hh"
75// #include "fastjet/contrib/Njettiness.hh"
76// #include "fastjet/contrib/NjettinessPlugin.hh"
77// #include "fastjet/contrib/WinnerTakeAllRecombiner.hh"
78
79using namespace std;
80using namespace fastjet;
81// using namespace fastjet::contrib;
82
83//---------------------------------------------------------------------------
84
85static bool interrupted = false;
86
87void SignalHandler(int sig)
88{
89 interrupted = true;
90}
91
92//---------------------------------------------------------------------------
93
94int main(int argc, char *argv[])
95{
96 char appName[] = "ExternalFastJetHepMC";
97 stringstream message;
98 FILE *inputFile = 0;
99 ExRootConfReader *confReader = 0;
100 Delphes *modularDelphes = 0;
101 DelphesFactory *factory = 0;
102 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
103 DelphesHepMCReader *reader = 0;
104 Int_t i, maxEvents, skipEvents;
105 Long64_t length, eventCounter;
106
107 TObjArray *inputArray = 0;
108 TIterator *inputIterator = 0;
109 Candidate *candidate = 0;
110 TLorentzVector momentum;
111
112 JetDefinition *definition = 0;
113// JetDefinition::Recombiner *recomb = 0;
114 vector<PseudoJet> inputList, outputList;
115 PseudoJet jet;
116
117 if(argc < 2)
118 {
119 cout << " Usage: " << appName << " config_file" << " [input_file(s)]" << endl;
120 cout << " config_file - configuration file in Tcl format," << endl;
121 cout << " input_file(s) - input file(s) in HepMC format," << endl;
122 cout << " with no input_file, or when input_file is -, read standard input." << endl;
123 return 1;
124 }
125
126 signal(SIGINT, SignalHandler);
127
128 gROOT->SetBatch();
129
130 int appargc = 1;
131 char *appargv[] = {appName};
132 TApplication app(appName, &appargc, appargv);
133
134 try
135 {
136 confReader = new ExRootConfReader;
137 confReader->ReadFile(argv[1]);
138
139 maxEvents = confReader->GetInt("::MaxEvents", 0);
140 skipEvents = confReader->GetInt("::SkipEvents", 0);
141
142 if(maxEvents < 0)
143 {
144 throw runtime_error("MaxEvents must be zero or positive");
145 }
146
147 if(skipEvents < 0)
148 {
149 throw runtime_error("SkipEvents must be zero or positive");
150 }
151
152 modularDelphes = new Delphes("Delphes");
153 modularDelphes->SetConfReader(confReader);
154
155 factory = modularDelphes->GetFactory();
156 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
157 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
158 partonOutputArray = modularDelphes->ExportArray("partons");
159
160 reader = new DelphesHepMCReader;
161
162 modularDelphes->InitTask();
163
164 ClusterSequence::print_banner();
165// recomb = new WinnerTakeAllRecombiner();
166// definition = new JetDefinition(antikt_algorithm, 0.5, recomb, Best);
167 definition = new JetDefinition(antikt_algorithm, 0.5);
168
169 inputArray = modularDelphes->ImportArray("Calorimeter/towers");
170 inputIterator = inputArray->MakeIterator();
171
172 i = 2;
173 do
174 {
175 if(interrupted) break;
176
177 if(i == argc || strncmp(argv[i], "-", 2) == 0)
178 {
179 cout << "** Reading standard input" << endl;
180 inputFile = stdin;
181 length = -1;
182 }
183 else
184 {
185 cout << "** Reading " << argv[i] << endl;
186 inputFile = fopen(argv[i], "r");
187
188 if(inputFile == NULL)
189 {
190 message << "can't open " << argv[i];
191 throw runtime_error(message.str());
192 }
193
194 fseek(inputFile, 0L, SEEK_END);
195 length = ftello(inputFile);
196 fseek(inputFile, 0L, SEEK_SET);
197
198 if(length <= 0)
199 {
200 fclose(inputFile);
201 ++i;
202 continue;
203 }
204 }
205
206 reader->SetInputFile(inputFile);
207
208 // Loop over all objects
209 eventCounter = 0;
210 modularDelphes->Clear();
211 reader->Clear();
212 while((maxEvents <= 0 || eventCounter - skipEvents < maxEvents) &&
213 reader->ReadBlock(factory, allParticleOutputArray,
214 stableParticleOutputArray, partonOutputArray) && !interrupted)
215 {
216 if(reader->EventReady())
217 {
218 ++eventCounter;
219
220 if(eventCounter > skipEvents)
221 {
222 modularDelphes->ProcessTask();
223
224 inputList.clear();
225 inputIterator->Reset();
226 while((candidate = static_cast<Candidate*>(inputIterator->Next())))
227 {
228 momentum = candidate->Momentum;
229 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
230 inputList.push_back(jet);
231 }
232 ClusterSequence sequence(inputList, *definition);
233 outputList.clear();
234 outputList = sorted_by_pt(sequence.inclusive_jets(0.0));
235
236 // tell the user what was done
237 // - the description of the algorithm used
238 // - show the inclusive jets as
239 // {index, rapidity, phi, pt}
240 //----------------------------------------------------------
241 if(eventCounter == skipEvents + 1)
242 {
243 cout << "Ran " << definition->description() << endl;
244
245 // label the columns
246 printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
247
248 // print out the details for each jet
249 for (unsigned int i = 0; i < outputList.size(); i++) {
250 printf("%5u %15.8f %15.8f %15.8f\n",
251 i, outputList[i].rap(), outputList[i].phi(),
252 outputList[i].perp());
253 }
254 }
255 }
256
257 modularDelphes->Clear();
258 reader->Clear();
259 }
260 }
261
262 if(inputFile != stdin) fclose(inputFile);
263
264 ++i;
265 }
266 while(i < argc);
267
268 modularDelphes->FinishTask();
269
270 cout << "** Exiting..." << endl;
271
272 delete definition;
273// delete recomb;
274
275 delete reader;
276 delete modularDelphes;
277 delete confReader;
278
279 return 0;
280 }
281 catch(runtime_error &e)
282 {
283 cerr << "** ERROR: " << e.what() << endl;
284 return 1;
285 }
286}
Note: See TracBrowser for help on using the repository browser.