Fork me on GitHub

source: git/examples/ExternalFastJetBasic.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: 7.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/*
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 ExternalFastJetBasic:
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 ExternalFastJetBasic.cpp -o ExternalFastJetBasic
41
42*/
43
44#include <stdexcept>
45#include <iostream>
46#include <sstream>
47
48#include <stdio.h>
49
50#include "TROOT.h"
51#include "TApplication.h"
52
53#include "TFile.h"
54#include "TObjArray.h"
55#include "TDatabasePDG.h"
56#include "TParticlePDG.h"
57#include "TLorentzVector.h"
58
59#include "modules/Delphes.h"
60#include "classes/DelphesClasses.h"
61#include "classes/DelphesFactory.h"
62
63#include "ExRootAnalysis/ExRootTreeWriter.h"
64#include "ExRootAnalysis/ExRootTreeBranch.h"
65#include "ExRootAnalysis/ExRootProgressBar.h"
66
67#include "fastjet/PseudoJet.hh"
68#include "fastjet/JetDefinition.hh"
69#include "fastjet/ClusterSequence.hh"
70
71using namespace std;
72using namespace fastjet;
73
74//---------------------------------------------------------------------------
75
76const int NEVENTS = 2;
77const int NPARTICLES = 4;
78
79float EVENTS[NEVENTS][NPARTICLES][11] =
80{
81 {
82 {-211, 1, -2.91e-01, 1.99e-01, -1.30e+00, 1.35e+00, 1.3957e-01, 0.0, 0.0, 0.0, 0.0},
83 { 211, 1, -2.62e-02, 1.90e-01, 1.42e+00, 1.44e+00, 1.3957e-01, 0.0, 0.0, 0.0, 0.0},
84 {-211, 1, -3.08e-02, -2.30e-01, 5.13e+00, 5.14e+00, 1.3957e-01, 0.0, 0.0, 0.0, 0.0},
85 { 211, 1, 1.21e-01, -3.52e-01, 2.12e+01, 2.12e+01, 1.3957e-01, 0.0, 0.0, 0.0, 0.0}
86 },
87 {
88 {22, 1, -7.80e-02, -9.70e-02, 2.03e+01, 2.03e+01, 0.0, 0.0, 0.0, 0.0, 0.0},
89 {22, 1, -1.16e-02, -2.08e-02, 4.89e-01, 4.90e-01, 0.0, 0.0, 0.0, 0.0, 0.0},
90 {22, 1, 5.64e-02, -1.34e-02, 5.81e+00, 5.82e+00, 0.0, 0.0, 0.0, 0.0, 0.0},
91 {22, 1, 2.35e-01, -1.04e-01, 9.12e+00, 9.12e+00, 0.0, 0.0, 0.0, 0.0, 0.0}
92 }
93};
94
95//---------------------------------------------------------------------------
96
97void ConvertInput(Int_t event, DelphesFactory *factory, TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray);
98
99int main(int argc, char *argv[])
100{
101 // Declaration of variables
102 ExRootConfReader *confReader = 0;
103 Delphes *modularDelphes = 0;
104 DelphesFactory *factory = 0;
105 TObjArray *allParticleOutputArray = 0;
106 TObjArray *stableParticleOutputArray = 0;
107 TObjArray *partonOutputArray = 0;
108
109 Int_t event;
110
111 TObjArray *inputArray = 0;
112 TIterator *inputIterator = 0;
113 Candidate *candidate = 0;
114 TLorentzVector momentum;
115
116 JetDefinition *definition = 0;
117 vector<PseudoJet> inputList, outputList;
118 PseudoJet jet;
119
120 gROOT->SetBatch();
121
122 int appargc = 1;
123 char appName[] = "ExternalFastJetBasic";
124 char *appargv[] = {appName};
125 TApplication app(appName, &appargc, appargv);
126
127 if(argc != 2)
128 {
129 cout << " Usage: " << appName << " config_file" << endl;
130 cout << " config_file - configuration file in Tcl format." << endl;
131 return 1;
132 }
133
134 try
135 {
136 // Initialization
137 confReader = new ExRootConfReader;
138 confReader->ReadFile(argv[1]);
139
140 modularDelphes = new Delphes("Delphes");
141 modularDelphes->SetConfReader(confReader);
142
143 factory = modularDelphes->GetFactory();
144
145 allParticleOutputArray = modularDelphes->ExportArray("allParticles");
146 stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
147 partonOutputArray = modularDelphes->ExportArray("partons");
148
149 modularDelphes->InitTask();
150
151 ClusterSequence::print_banner();
152 definition = new JetDefinition(antikt_algorithm, 0.5);
153
154 inputArray = modularDelphes->ImportArray("Calorimeter/towers");
155 inputIterator = inputArray->MakeIterator();
156
157 // Event loop
158 for(event = 0; event < NEVENTS; ++event)
159 {
160 modularDelphes->Clear();
161 ConvertInput(event, factory, allParticleOutputArray, stableParticleOutputArray, partonOutputArray);
162 modularDelphes->ProcessTask();
163
164 inputList.clear();
165 inputIterator->Reset();
166 while((candidate = static_cast<Candidate*>(inputIterator->Next())))
167 {
168 momentum = candidate->Momentum;
169 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
170 inputList.push_back(jet);
171 }
172 ClusterSequence sequence(inputList, *definition);
173 outputList.clear();
174 outputList = sorted_by_pt(sequence.inclusive_jets(0.0));
175
176 // tell the user what was done
177 // - the description of the algorithm used
178 // - show the inclusive jets as
179 // {index, rapidity, phi, pt}
180 //----------------------------------------------------------
181 cout << "Ran " << definition->description() << endl;
182
183 // label the columns
184 printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
185
186 // print out the details for each jet
187 for (unsigned int i = 0; i < outputList.size(); i++) {
188 printf("%5u %15.8f %15.8f %15.8f\n",
189 i, outputList[i].rap(), outputList[i].phi(),
190 outputList[i].perp());
191 }
192 }
193
194 // Finalization
195 modularDelphes->FinishTask();
196 delete modularDelphes;
197 delete confReader;
198 return 0;
199 }
200 catch(runtime_error &e)
201 {
202 cerr << "** ERROR: " << e.what() << endl;
203 return 1;
204 }
205}
206
207void ConvertInput(Int_t event, DelphesFactory *factory, TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray, TObjArray *partonOutputArray)
208{
209 Int_t particle;
210
211 Candidate *candidate;
212 TDatabasePDG *pdg;
213 TParticlePDG *pdgParticle;
214 Int_t pdgCode;
215
216 Int_t pid, status;
217 Double_t px, py, pz, e, mass;
218 Double_t x, y, z, t;
219
220 pdg = TDatabasePDG::Instance();
221
222 for(particle = 0; particle < NPARTICLES; ++particle)
223 {
224 pid = EVENTS[event][particle][0];
225 status = EVENTS[event][particle][1];
226 px = EVENTS[event][particle][2];
227 py = EVENTS[event][particle][3];
228 pz = EVENTS[event][particle][4];
229 e = EVENTS[event][particle][5];
230 mass = EVENTS[event][particle][6];
231 x = EVENTS[event][particle][7];
232 y = EVENTS[event][particle][8];
233 z = EVENTS[event][particle][9];
234 t = EVENTS[event][particle][10];
235
236 candidate = factory->NewCandidate();
237
238 candidate->PID = pid;
239 pdgCode = TMath::Abs(candidate->PID);
240
241 candidate->Status = status;
242
243 pdgParticle = pdg->GetParticle(pid);
244 candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
245 candidate->Mass = mass;
246
247 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
248
249 candidate->Position.SetXYZT(x, y, z, t);
250
251 allParticleOutputArray->Add(candidate);
252
253 if(!pdgParticle) return;
254
255 if(status == 1)
256 {
257 stableParticleOutputArray->Add(candidate);
258 }
259 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
260 {
261 partonOutputArray->Add(candidate);
262 }
263 }
264}
265
Note: See TracBrowser for help on using the repository browser.