Fork me on GitHub

source: git/converters/root2pileup.cpp@ 3874a9c

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

fix GPLv3 header

  • Property mode set to 100644
File size: 3.8 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#include <string>
23
24#include <signal.h>
25
26#include "TROOT.h"
27#include "TApplication.h"
28
29#include "TFile.h"
30#include "TClonesArray.h"
31
32#include "classes/DelphesClasses.h"
33#include "classes/DelphesPileUpWriter.h"
34
35#include "ExRootAnalysis/ExRootTreeReader.h"
36#include "ExRootAnalysis/ExRootProgressBar.h"
37
38using namespace std;
39
40//---------------------------------------------------------------------------
41
42static bool interrupted = false;
43
44void SignalHandler(int sig)
45{
46 interrupted = true;
47}
48
49//---------------------------------------------------------------------------
50
51int main(int argc, char *argv[])
52{
53 char appName[] = "root2pileup";
54 stringstream message;
55 TChain *inputChain = 0;
56 ExRootTreeReader *treeReader = 0;
57 TClonesArray *branchParticle = 0;
58 TIterator *itParticle = 0;
59 GenParticle *particle = 0;
60 DelphesPileUpWriter *writer = 0;
61 Long64_t entry, allEntries;
62 Int_t i;
63
64 if(argc < 3)
65 {
66 cout << " Usage: " << appName << " output_file" << " input_file(s)" << endl;
67 cout << " output_file - output binary pile-up file," << endl;
68 cout << " input_file(s) - input file(s) in ROOT format." << endl;
69 return 1;
70 }
71
72 signal(SIGINT, SignalHandler);
73
74 gROOT->SetBatch();
75
76 int appargc = 1;
77 char *appargv[] = {appName};
78 TApplication app(appName, &appargc, appargv);
79
80 try
81 {
82 inputChain = new TChain("Delphes");
83 for(i = 2; i < argc && !interrupted; ++i)
84 {
85 inputChain->Add(argv[i]);
86 }
87
88 treeReader = new ExRootTreeReader(inputChain);
89 branchParticle = treeReader->UseBranch("Particle");
90 itParticle = branchParticle->MakeIterator();
91
92 writer = new DelphesPileUpWriter(argv[1]);
93
94 allEntries = treeReader->GetEntries();
95 cout << "** Input file(s) contain(s) " << allEntries << " events" << endl;
96
97 if(allEntries > 0)
98 {
99 ExRootProgressBar progressBar(allEntries - 1);
100 // Loop over all events in the input file
101 for(entry = 0; entry < allEntries && !interrupted; ++entry)
102 {
103 if(!treeReader->ReadEntry(entry))
104 {
105 cerr << "** ERROR: cannot read event " << entry << endl;
106 break;
107 }
108
109 itParticle->Reset();
110 while((particle = static_cast<GenParticle*>(itParticle->Next())))
111 {
112 writer->WriteParticle(particle->PID,
113 particle->X, particle->Y, particle->Z, particle->T,
114 particle->Px, particle->Py, particle->Pz, particle->E);
115 }
116
117 writer->WriteEntry();
118
119 progressBar.Update(entry);
120 }
121 progressBar.Finish();
122
123 writer->WriteIndex();
124 }
125
126 cout << "** Exiting..." << endl;
127
128 delete writer;
129 delete itParticle;
130 delete treeReader;
131 delete inputChain;
132 return 0;
133 }
134 catch(runtime_error &e)
135 {
136 if(writer) delete writer;
137 if(itParticle) delete itParticle;
138 if(treeReader) delete treeReader;
139 if(inputChain) delete inputChain;
140 cerr << "** ERROR: " << e.what() << endl;
141 return 1;
142 }
143}
Note: See TracBrowser for help on using the repository browser.