[01f457a] | 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 |
|
---|
[975405a] | 19 | #include <algorithm>
|
---|
[d7d2da3] | 20 | #include <stdexcept>
|
---|
| 21 | #include <iostream>
|
---|
| 22 | #include <sstream>
|
---|
| 23 | #include <memory>
|
---|
| 24 |
|
---|
| 25 | #include <map>
|
---|
[975405a] | 26 | #include <vector>
|
---|
[d7d2da3] | 27 |
|
---|
| 28 | #include <stdlib.h>
|
---|
| 29 | #include <signal.h>
|
---|
| 30 | #include <stdio.h>
|
---|
| 31 |
|
---|
| 32 | #include "TROOT.h"
|
---|
| 33 | #include "TApplication.h"
|
---|
| 34 |
|
---|
| 35 | #include "TFile.h"
|
---|
| 36 | #include "TObjArray.h"
|
---|
| 37 | #include "TStopwatch.h"
|
---|
| 38 | #include "TDatabasePDG.h"
|
---|
| 39 | #include "TParticlePDG.h"
|
---|
| 40 | #include "TLorentzVector.h"
|
---|
| 41 |
|
---|
| 42 | #include "modules/Delphes.h"
|
---|
| 43 | #include "classes/DelphesStream.h"
|
---|
| 44 | #include "classes/DelphesClasses.h"
|
---|
| 45 | #include "classes/DelphesFactory.h"
|
---|
| 46 |
|
---|
| 47 | #include "ExRootAnalysis/ExRootTreeWriter.h"
|
---|
| 48 | #include "ExRootAnalysis/ExRootTreeBranch.h"
|
---|
| 49 | #include "ExRootAnalysis/ExRootProgressBar.h"
|
---|
| 50 |
|
---|
| 51 | #include "FWCore/FWLite/interface/AutoLibraryLoader.h"
|
---|
| 52 |
|
---|
| 53 | #include "DataFormats/FWLite/interface/Event.h"
|
---|
| 54 | #include "DataFormats/FWLite/interface/Handle.h"
|
---|
| 55 | #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
|
---|
[f29758e] | 56 | #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
|
---|
| 57 | #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
|
---|
[d7d2da3] | 58 |
|
---|
| 59 | using namespace std;
|
---|
| 60 |
|
---|
| 61 | //---------------------------------------------------------------------------
|
---|
| 62 |
|
---|
[f29758e] | 63 | void ConvertInput(fwlite::Event &event, Long64_t eventCounter,
|
---|
| 64 | ExRootTreeBranch *branchEvent, DelphesFactory *factory,
|
---|
| 65 | TObjArray *allParticleOutputArray, TObjArray *stableParticleOutputArray,
|
---|
| 66 | TObjArray *partonOutputArray)
|
---|
[d7d2da3] | 67 | {
|
---|
[f29758e] | 68 | fwlite::Handle< GenEventInfoProduct > handleGenEventInfo;
|
---|
| 69 |
|
---|
[d7d2da3] | 70 | fwlite::Handle< vector< reco::GenParticle > > handleParticle;
|
---|
[975405a] | 71 | vector< reco::GenParticle >::const_iterator itParticle;
|
---|
| 72 |
|
---|
| 73 | vector< const reco::Candidate * > vectorCandidate;
|
---|
| 74 | vector< const reco::Candidate * >::iterator itCandidate;
|
---|
| 75 |
|
---|
[f29758e] | 76 | handleGenEventInfo.getByLabel(event, "generator");
|
---|
[d7d2da3] | 77 | handleParticle.getByLabel(event, "genParticles");
|
---|
| 78 |
|
---|
[f29758e] | 79 | HepMCEvent *element;
|
---|
[d7d2da3] | 80 | Candidate *candidate;
|
---|
| 81 | TDatabasePDG *pdg;
|
---|
| 82 | TParticlePDG *pdgParticle;
|
---|
| 83 | Int_t pdgCode;
|
---|
| 84 |
|
---|
| 85 | Int_t pid, status;
|
---|
[80d4a34] | 86 | Double_t px, py, pz, e, mass;
|
---|
[d7d2da3] | 87 | Double_t x, y, z;
|
---|
| 88 |
|
---|
[f29758e] | 89 | element = static_cast<HepMCEvent *>(branchEvent->NewEntry());
|
---|
| 90 |
|
---|
| 91 | element->Number = eventCounter;
|
---|
| 92 |
|
---|
| 93 | element->ProcessID = handleGenEventInfo->signalProcessID();
|
---|
| 94 | element->MPI = 1;
|
---|
[529fe78] | 95 | element->Weight = handleGenEventInfo->weight();
|
---|
[f29758e] | 96 | element->Scale = handleGenEventInfo->qScale();
|
---|
| 97 | element->AlphaQED = handleGenEventInfo->alphaQED();
|
---|
| 98 | element->AlphaQCD = handleGenEventInfo->alphaQCD();
|
---|
| 99 |
|
---|
| 100 | element->ID1 = 0;
|
---|
| 101 | element->ID2 = 0;
|
---|
| 102 | element->X1 = 0.0;
|
---|
| 103 | element->X2 = 0.0;
|
---|
| 104 | element->ScalePDF = 0.0;
|
---|
| 105 | element->PDF1 = 0.0;
|
---|
| 106 | element->PDF2 = 0.0;
|
---|
| 107 |
|
---|
| 108 | element->ReadTime = 0.0;
|
---|
| 109 | element->ProcTime = 0.0;
|
---|
| 110 |
|
---|
[d7d2da3] | 111 | pdg = TDatabasePDG::Instance();
|
---|
| 112 |
|
---|
[975405a] | 113 | for(itParticle = handleParticle->begin(); itParticle != handleParticle->end(); ++itParticle)
|
---|
[d7d2da3] | 114 | {
|
---|
[975405a] | 115 | vectorCandidate.push_back(&*itParticle);
|
---|
| 116 | }
|
---|
| 117 |
|
---|
| 118 | for(itParticle = handleParticle->begin(); itParticle != handleParticle->end(); ++itParticle)
|
---|
| 119 | {
|
---|
| 120 | const reco::GenParticle &particle = *itParticle;
|
---|
[d7d2da3] | 121 |
|
---|
[1e1f73f] | 122 | pid = particle.pdgId();
|
---|
| 123 | status = particle.status();
|
---|
[80d4a34] | 124 | px = particle.px(); py = particle.py(); pz = particle.pz(); e = particle.energy(); mass = particle.mass();
|
---|
[1e1f73f] | 125 | x = particle.vx(); y = particle.vy(); z = particle.vz();
|
---|
[d7d2da3] | 126 |
|
---|
[1e1f73f] | 127 | candidate = factory->NewCandidate();
|
---|
[d7d2da3] | 128 |
|
---|
[1e1f73f] | 129 | candidate->PID = pid;
|
---|
| 130 | pdgCode = TMath::Abs(candidate->PID);
|
---|
[d7d2da3] | 131 |
|
---|
[1e1f73f] | 132 | candidate->Status = status;
|
---|
[d7d2da3] | 133 |
|
---|
[f29758e] | 134 | if(particle.mother())
|
---|
| 135 | {
|
---|
| 136 | itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.mother());
|
---|
| 137 | if(itCandidate != vectorCandidate.end()) candidate->M1 = distance(vectorCandidate.begin(), itCandidate);
|
---|
| 138 | }
|
---|
| 139 |
|
---|
[975405a] | 140 | itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(0));
|
---|
| 141 | if(itCandidate != vectorCandidate.end()) candidate->D1 = distance(vectorCandidate.begin(), itCandidate);
|
---|
| 142 |
|
---|
| 143 | itCandidate = find(vectorCandidate.begin(), vectorCandidate.end(), particle.daughter(particle.numberOfDaughters() - 1));
|
---|
| 144 | if(itCandidate != vectorCandidate.end()) candidate->D2 = distance(vectorCandidate.begin(), itCandidate);
|
---|
| 145 |
|
---|
[1e1f73f] | 146 | pdgParticle = pdg->GetParticle(pid);
|
---|
| 147 | candidate->Charge = pdgParticle ? Int_t(pdgParticle->Charge()/3.0) : -999;
|
---|
[bba4646] | 148 | candidate->Mass = mass;
|
---|
[d7d2da3] | 149 |
|
---|
[1e1f73f] | 150 | candidate->Momentum.SetPxPyPzE(px, py, pz, e);
|
---|
[d7d2da3] | 151 |
|
---|
[1e1f73f] | 152 | candidate->Position.SetXYZT(x, y, z, 0.0);
|
---|
[d7d2da3] | 153 |
|
---|
[1e1f73f] | 154 | allParticleOutputArray->Add(candidate);
|
---|
[d7d2da3] | 155 |
|
---|
[2f82259] | 156 | if(!pdgParticle) continue;
|
---|
[d7d2da3] | 157 |
|
---|
[1e1f73f] | 158 | if(status == 1)
|
---|
| 159 | {
|
---|
| 160 | stableParticleOutputArray->Add(candidate);
|
---|
| 161 | }
|
---|
| 162 | else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
|
---|
| 163 | {
|
---|
| 164 | partonOutputArray->Add(candidate);
|
---|
[d7d2da3] | 165 | }
|
---|
| 166 | }
|
---|
| 167 | }
|
---|
| 168 |
|
---|
| 169 | //---------------------------------------------------------------------------
|
---|
| 170 |
|
---|
| 171 | static bool interrupted = false;
|
---|
| 172 |
|
---|
| 173 | void SignalHandler(int sig)
|
---|
| 174 | {
|
---|
| 175 | interrupted = true;
|
---|
| 176 | }
|
---|
| 177 |
|
---|
| 178 | //---------------------------------------------------------------------------
|
---|
| 179 |
|
---|
| 180 | int main(int argc, char *argv[])
|
---|
| 181 | {
|
---|
| 182 | char appName[] = "DelphesCMSFWLite";
|
---|
| 183 | stringstream message;
|
---|
| 184 | TFile *inputFile = 0;
|
---|
| 185 | TFile *outputFile = 0;
|
---|
| 186 | TStopwatch eventStopWatch;
|
---|
| 187 | ExRootTreeWriter *treeWriter = 0;
|
---|
[f29758e] | 188 | ExRootTreeBranch *branchEvent = 0;
|
---|
[d7d2da3] | 189 | ExRootConfReader *confReader = 0;
|
---|
| 190 | Delphes *modularDelphes = 0;
|
---|
| 191 | DelphesFactory *factory = 0;
|
---|
| 192 | TObjArray *allParticleOutputArray = 0, *stableParticleOutputArray = 0, *partonOutputArray = 0;
|
---|
| 193 | Int_t i;
|
---|
[5ca3d52] | 194 | Long64_t eventCounter, numberOfEvents;
|
---|
[d7d2da3] | 195 |
|
---|
| 196 | if(argc < 4)
|
---|
| 197 | {
|
---|
| 198 | cout << " Usage: " << appName << " config_file" << " output_file" << " input_file(s)" << endl;
|
---|
| 199 | cout << " config_file - configuration file in Tcl format," << endl;
|
---|
| 200 | cout << " output_file - output file in ROOT format," << endl;
|
---|
| 201 | cout << " input_file(s) - input file(s) in ROOT format." << endl;
|
---|
| 202 | return 1;
|
---|
| 203 | }
|
---|
| 204 |
|
---|
| 205 | signal(SIGINT, SignalHandler);
|
---|
| 206 |
|
---|
| 207 | gROOT->SetBatch();
|
---|
| 208 |
|
---|
| 209 | int appargc = 1;
|
---|
| 210 | char *appargv[] = {appName};
|
---|
| 211 | TApplication app(appName, &appargc, appargv);
|
---|
| 212 |
|
---|
| 213 | AutoLibraryLoader::enable();
|
---|
| 214 |
|
---|
| 215 | try
|
---|
| 216 | {
|
---|
| 217 | outputFile = TFile::Open(argv[2], "CREATE");
|
---|
| 218 |
|
---|
| 219 | if(outputFile == NULL)
|
---|
| 220 | {
|
---|
| 221 | message << "can't open " << argv[2] << endl;
|
---|
| 222 | throw runtime_error(message.str());
|
---|
| 223 | }
|
---|
| 224 |
|
---|
| 225 | treeWriter = new ExRootTreeWriter(outputFile, "Delphes");
|
---|
| 226 |
|
---|
[f29758e] | 227 | branchEvent = treeWriter->NewBranch("Event", HepMCEvent::Class());
|
---|
| 228 |
|
---|
[d7d2da3] | 229 | confReader = new ExRootConfReader;
|
---|
| 230 | confReader->ReadFile(argv[1]);
|
---|
| 231 |
|
---|
| 232 | modularDelphes = new Delphes("Delphes");
|
---|
| 233 | modularDelphes->SetConfReader(confReader);
|
---|
| 234 | modularDelphes->SetTreeWriter(treeWriter);
|
---|
| 235 |
|
---|
| 236 | factory = modularDelphes->GetFactory();
|
---|
| 237 | allParticleOutputArray = modularDelphes->ExportArray("allParticles");
|
---|
| 238 | stableParticleOutputArray = modularDelphes->ExportArray("stableParticles");
|
---|
| 239 | partonOutputArray = modularDelphes->ExportArray("partons");
|
---|
| 240 |
|
---|
| 241 | modularDelphes->InitTask();
|
---|
| 242 |
|
---|
| 243 | for(i = 3; i < argc && !interrupted; ++i)
|
---|
| 244 | {
|
---|
| 245 | cout << "** Reading " << argv[i] << endl;
|
---|
| 246 |
|
---|
| 247 | inputFile = TFile::Open(argv[i]);
|
---|
| 248 |
|
---|
| 249 | if(inputFile == NULL)
|
---|
| 250 | {
|
---|
| 251 | message << "can't open " << argv[i] << endl;
|
---|
| 252 | throw runtime_error(message.str());
|
---|
| 253 | }
|
---|
| 254 |
|
---|
| 255 | fwlite::Event event(inputFile);
|
---|
| 256 |
|
---|
[5ca3d52] | 257 | numberOfEvents = event.size();
|
---|
[d7d2da3] | 258 |
|
---|
[5ca3d52] | 259 | if(numberOfEvents <= 0) continue;
|
---|
[d7d2da3] | 260 |
|
---|
[a0538b9] | 261 | // ExRootProgressBar progressBar(numberOfEvents - 1);
|
---|
| 262 | ExRootProgressBar progressBar(-1);
|
---|
[d7d2da3] | 263 |
|
---|
| 264 | // Loop over all objects
|
---|
[5ca3d52] | 265 | eventCounter = 0;
|
---|
[d7d2da3] | 266 | modularDelphes->Clear();
|
---|
| 267 | treeWriter->Clear();
|
---|
| 268 | for(event.toBegin(); !event.atEnd() && !interrupted; ++event)
|
---|
| 269 | {
|
---|
[f29758e] | 270 | ConvertInput(event, eventCounter, branchEvent, factory,
|
---|
| 271 | allParticleOutputArray, stableParticleOutputArray, partonOutputArray);
|
---|
[d7d2da3] | 272 | modularDelphes->ProcessTask();
|
---|
| 273 |
|
---|
| 274 | treeWriter->Fill();
|
---|
| 275 |
|
---|
| 276 | modularDelphes->Clear();
|
---|
| 277 | treeWriter->Clear();
|
---|
| 278 |
|
---|
[a0538b9] | 279 | progressBar.Update(eventCounter, eventCounter);
|
---|
[5ca3d52] | 280 | ++eventCounter;
|
---|
[d7d2da3] | 281 | }
|
---|
[a0538b9] | 282 |
|
---|
| 283 | progressBar.Update(eventCounter, eventCounter, kTRUE);
|
---|
[d7d2da3] | 284 | progressBar.Finish();
|
---|
| 285 |
|
---|
| 286 | inputFile->Close();
|
---|
| 287 | }
|
---|
| 288 |
|
---|
| 289 | modularDelphes->FinishTask();
|
---|
| 290 | treeWriter->Write();
|
---|
| 291 |
|
---|
| 292 | cout << "** Exiting..." << endl;
|
---|
| 293 |
|
---|
| 294 | delete modularDelphes;
|
---|
| 295 | delete confReader;
|
---|
| 296 | delete treeWriter;
|
---|
| 297 | delete outputFile;
|
---|
| 298 |
|
---|
| 299 | return 0;
|
---|
| 300 | }
|
---|
| 301 | catch(runtime_error &e)
|
---|
| 302 | {
|
---|
| 303 | if(treeWriter) delete treeWriter;
|
---|
| 304 | if(outputFile) delete outputFile;
|
---|
| 305 | cerr << "** ERROR: " << e.what() << endl;
|
---|
| 306 | return 1;
|
---|
| 307 | }
|
---|
| 308 | }
|
---|