Fork me on GitHub

source: git/examples/StandaloneHepMC.cpp@ 2eb25b1

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 2eb25b1 was 2eb25b1, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

add examples/StandaloneHepMC.cpp and libDelphesNoFastJet

  • Property mode set to 100644
File size: 6.7 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#include <stdexcept>
20#include <iostream>
21#include <sstream>
22
23#include <vector>
24
25#include <signal.h>
26
27#include "TROOT.h"
28#include "TApplication.h"
29
30#include "TFile.h"
31#include "TObjArray.h"
32#include "TDatabasePDG.h"
33#include "TParticlePDG.h"
34#include "TLorentzVector.h"
35
36#include "modules/Delphes.h"
37#include "classes/DelphesClasses.h"
38#include "classes/DelphesFactory.h"
39#include "classes/DelphesHepMCReader.h"
40
41#include "ExRootAnalysis/ExRootTreeWriter.h"
42#include "ExRootAnalysis/ExRootTreeBranch.h"
43#include "ExRootAnalysis/ExRootProgressBar.h"
44
45#include "fastjet/PseudoJet.hh"
46#include "fastjet/JetDefinition.hh"
47#include "fastjet/ClusterSequence.hh"
48#include "fastjet/Selector.hh"
49#include "fastjet/ClusterSequenceArea.hh"
50#include "fastjet/tools/JetMedianBackgroundEstimator.hh"
51
52#include "fastjet/plugins/SISCone/fastjet/SISConePlugin.hh"
53#include "fastjet/plugins/CDFCones/fastjet/CDFMidPointPlugin.hh"
54#include "fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh"
55
56#include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh"
57#include "fastjet/contribs/Nsubjettiness/Njettiness.hh"
58#include "fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh"
59#include "fastjet/contribs/Nsubjettiness/WinnerTakeAllRecombiner.hh"
60
61using namespace std;
62using namespace fastjet;
63using namespace fastjet::contrib;
64
65//---------------------------------------------------------------------------
66
67static bool interrupted = false;
68
69void SignalHandler(int sig)
70{
71 interrupted = true;
72}
73
74//---------------------------------------------------------------------------
75
76int main(int argc, char *argv[])
77{
78 char appName[] = "StandaloneHepMC";
79 stringstream message;
80 FILE *inputFile = 0;
81 ExRootConfReader *confReader = 0;
82 Delphes *modularDelphes = 0;
83 DelphesFactory *factory = 0;
84 TObjArray *stableParticleOutputArray = 0, *allParticleOutputArray = 0, *partonOutputArray = 0;
85 DelphesHepMCReader *reader = 0;
86 Int_t i, maxEvents, skipEvents;
87 Long64_t length, eventCounter;
88
89 TObjArray *inputArray = 0;
90 TIterator *inputIterator = 0;
91 Candidate *candidate = 0;
92 TLorentzVector momentum;
93
94 JetDefinition *definition = 0;
95 JetDefinition::Recombiner *recomb = 0;
96 vector<PseudoJet> inputList, outputList;
97 PseudoJet jet;
98
99 if(argc < 2)
100 {
101 cout << " Usage: " << appName << " config_file" << " [input_file(s)]" << endl;
102 cout << " config_file - configuration file in Tcl format," << endl;
103 cout << " input_file(s) - input file(s) in HepMC format," << endl;
104 cout << " with no input_file, or when input_file is -, read standard input." << endl;
105 return 1;
106 }
107
108 signal(SIGINT, SignalHandler);
109
110 gROOT->SetBatch();
111
112 int appargc = 1;
113 char *appargv[] = {appName};
114 TApplication app(appName, &appargc, appargv);
115
116 try
117 {
118 confReader = new ExRootConfReader;
119 confReader->ReadFile(argv[1]);
120
121 maxEvents = confReader->GetInt("::MaxEvents", 0);
122 skipEvents = confReader->GetInt("::SkipEvents", 0);
123
124 if(maxEvents < 0)
125 {
126 throw runtime_error("MaxEvents must be zero or positive");
127 }
128
129 if(skipEvents < 0)
130 {
131 throw runtime_error("SkipEvents must be zero or positive");
132 }
133
134 modularDelphes = new Delphes("Delphes");
135 modularDelphes->SetConfReader(confReader);
136
137 factory = modularDelphes->GetFactory();
138 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
139 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
140 partonOutputArray = modularDelphes->ExportArray("partons");
141
142 reader = new DelphesHepMCReader;
143
144 modularDelphes->InitTask();
145
146 ClusterSequence::print_banner();
147 recomb = new WinnerTakeAllRecombiner();
148 definition = new JetDefinition(antikt_algorithm, 0.5, recomb, Best);
149
150 inputArray = modularDelphes->ImportArray("Calorimeter/towers");
151 inputIterator = inputArray->MakeIterator();
152
153 i = 2;
154 do
155 {
156 if(interrupted) break;
157
158 if(i == argc || strncmp(argv[i], "-", 2) == 0)
159 {
160 cout << "** Reading standard input" << endl;
161 inputFile = stdin;
162 length = -1;
163 }
164 else
165 {
166 cout << "** Reading " << argv[i] << endl;
167 inputFile = fopen(argv[i], "r");
168
169 if(inputFile == NULL)
170 {
171 message << "can't open " << argv[i];
172 throw runtime_error(message.str());
173 }
174
175 fseek(inputFile, 0L, SEEK_END);
176 length = ftello(inputFile);
177 fseek(inputFile, 0L, SEEK_SET);
178
179 if(length <= 0)
180 {
181 fclose(inputFile);
182 ++i;
183 continue;
184 }
185 }
186
187 reader->SetInputFile(inputFile);
188
189 // Loop over all objects
190 eventCounter = 0;
191 modularDelphes->Clear();
192 reader->Clear();
193 while((maxEvents <= 0 || eventCounter - skipEvents < maxEvents) &&
194 reader->ReadBlock(factory, allParticleOutputArray,
195 stableParticleOutputArray, partonOutputArray) && !interrupted)
196 {
197 if(reader->EventReady())
198 {
199 ++eventCounter;
200
201 if(eventCounter > skipEvents)
202 {
203 modularDelphes->ProcessTask();
204
205 inputList.clear();
206 inputIterator->Reset();
207 while((candidate = static_cast<Candidate*>(inputIterator->Next())))
208 {
209 momentum = candidate->Momentum;
210 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
211 inputList.push_back(jet);
212 }
213 ClusterSequence sequence(inputList, *definition);
214 outputList.clear();
215 outputList = sorted_by_pt(sequence.inclusive_jets(10.0));
216 }
217
218 modularDelphes->Clear();
219 reader->Clear();
220 }
221 }
222
223 if(inputFile != stdin) fclose(inputFile);
224
225 ++i;
226 }
227 while(i < argc);
228
229 modularDelphes->FinishTask();
230
231 cout << "** Exiting..." << endl;
232
233 delete reader;
234 delete modularDelphes;
235 delete confReader;
236
237 return 0;
238 }
239 catch(runtime_error &e)
240 {
241 cerr << "** ERROR: " << e.what() << endl;
242 return 1;
243 }
244}
Note: See TracBrowser for help on using the repository browser.