/* ---- Delphes ---- A Fast Simulator for general purpose LHC detector S. Ovyn ~~~~ severine.ovyn@uclouvain.be Center for Particle Physics and Phenomenology (CP3) Universite Catholique de Louvain (UCL) Louvain-la-Neuve, Belgium */ /// \file Delphes.cpp /// \brief executable for the Delphes #include "TChain.h" #include "TApplication.h" #include "Utilities/ExRootAnalysis/interface/ExRootTreeReader.h" #include "Utilities/ExRootAnalysis/interface/ExRootTreeWriter.h" #include "Utilities/ExRootAnalysis/interface/ExRootTreeBranch.h" #include "interface/DataConverter.h" #include "interface/HEPEVTConverter.h" #include "interface/LHEFConverter.h" #include "interface/STDHEPConverter.h" #include "interface/SmearUtil.h" #include "interface/BFieldProp.h" #include "interface/TriggerUtil.h" #include "interface/VeryForward.h" #include "interface/JetUtils.h" #include "interface/FrogUtil.h" #include #include using namespace std; //------------------------------------------------------------------------------ void todo(string filename) { ifstream infile(filename.c_str()); cout << "** TODO list ..." << endl; while(infile.good()) { string temp; getline(infile,temp); cout << "*" << temp << endl; } cout << "** done...\n"; } //------------------------------------------------------------------------------ int main(int argc, char *argv[]) { int appargc = 2; char *appName = "Delphes"; char *appargv[] = {appName, "-b"}; TApplication app(appName, &appargc, appargv); if(argc != 4 && argc != 3 && argc != 5) { cout << " Usage: " << argv[0] << " input_file output_file [detector_card] [trigger_card] " << endl; cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl; cout << " output_file - output file." << endl; cout << " detector_card - Datacard containing resolution variables for the detector simulation (optional) "< outputfilename.length() ) { cout << "output_file should be a .root file!\n"; exit(1); } //create output log-file name string forLog = outputfilename; string LogName = forLog.erase(forLog.find(".root")); LogName = LogName+"_run.log"; TFile *outputFile = TFile::Open(outputfilename.c_str(), "RECREATE"); // Creates the file, but should be closed just after outputFile->Close(); string line; ifstream infile(inputFileList.c_str()); infile >> line; // the first line determines the type of input files //read the datacard input file string DetDatacard(""); if(argc>=4) DetDatacard =argv[3]; //Smearing information RESOLution *DET = new RESOLution(); DET->ReadDataCard(DetDatacard); DET->Logfile(LogName); //read the trigger input file string TrigDatacard("data/trigger.dat"); if(argc==5) TrigDatacard =argv[4]; //Trigger information TriggerTable *TRIGT = new TriggerTable(); TRIGT->TriggerCardReader(TrigDatacard.c_str()); TRIGT->PrintTriggerTable(LogName); //Propagation of tracks in the B field TrackPropagation *TRACP = new TrackPropagation(); //Jet information JetsUtil *JETRUN = new JetsUtil(); //VFD information VeryForward * VFD = new VeryForward(); //todo(LogName.c_str()); DataConverter *converter=0; if(strstr(line.c_str(),".hep")) { cout<<"#**********************************************************************"<UseBranch("Particle"); TIter itGen((TCollection*)branchGen); //write the output root file ExRootTreeWriter *treeWriter = new ExRootTreeWriter(outputfilename, "Analysis"); ExRootTreeBranch *branchJet = treeWriter->NewBranch("Jet", TRootJet::Class()); ExRootTreeBranch *branchTauJet = treeWriter->NewBranch("TauJet", TRootTauJet::Class()); ExRootTreeBranch *branchElectron = treeWriter->NewBranch("Electron", TRootElectron::Class()); ExRootTreeBranch *branchMuon = treeWriter->NewBranch("Muon", TRootMuon::Class()); ExRootTreeBranch *branchPhoton = treeWriter->NewBranch("Photon", TRootPhoton::Class()); ExRootTreeBranch *branchTracks = treeWriter->NewBranch("Tracks", TRootTracks::Class()); ExRootTreeBranch *branchETmis = treeWriter->NewBranch("ETmis", TRootETmis::Class()); ExRootTreeBranch *branchCalo = treeWriter->NewBranch("CaloTower", TRootCalo::Class()); ExRootTreeBranch *branchZDC = treeWriter->NewBranch("ZDChits", TRootZdcHits::Class()); ExRootTreeBranch *branchRP220 = treeWriter->NewBranch("RP220hits", TRootRomanPotHits::Class()); ExRootTreeBranch *branchFP420 = treeWriter->NewBranch("FP420hits", TRootRomanPotHits::Class()); TRootGenParticle *particle; TRootETmis *elementEtmis; TRootElectron *elementElec; TRootMuon *elementMu; TRootPhoton *elementPhoton; TRootTracks *elementTracks; TRootCalo *elementCalo; TLorentzVector genMomentum(0,0,0,0); TLorentzVector genMomentumCalo(0,0,0,0); LorentzVector jetMomentum; vector input_particles;//for FastJet algorithm vector sorted_jets; vector TrackCentral; vector towers; vector electron; vector muon; vector gamma; TSimpleArray NFCentralQ; // Loop over all events Long64_t entry, allEntries = treeReader->GetEntries(); cout << "** Chain contains " << allEntries << " events" << endl; for(entry = 0; entry < allEntries; ++entry) { TLorentzVector PTmis(0,0,0,0); treeReader->ReadEntry(entry); treeWriter->Clear(); if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl; electron.clear(); muon.clear(); gamma.clear(); NFCentralQ.Clear(); TrackCentral.clear(); towers.clear(); input_particles.clear(); // Loop over all particles in event itGen.Reset(); while( (particle = (TRootGenParticle*) itGen.Next()) ) { int pid = abs(particle->PID); //// This subarray is needed for the B-jet algorithm // optimization for speed : put first PID condition, then ETA condition, then either pt or status if( (pid <= pB || pid == pGLUON) &&// is it a light quark or a gluon, i.e. is it one of these : u,d,c,s,b,g ? fabs(particle->Eta) < DET->CEN_max_tracker && particle->Status != 1 && particle->PT > DET->PT_QUARKS_MIN ) { NFCentralQ.Add(particle); } // keeps only final particles, visible by the central detector, including the fiducial volume // the ordering of conditions have been optimised for speed : put first the STATUS condition // // if( (particle->Status == 1) && ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) && (fabs(particle->Eta) < DET->CEN_max_calo_fwd) ) { genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E); if(DET->FLAG_bfield==1)TRACP->Propagation(particle,genMomentum); float eta=fabs(genMomentum.Eta()); switch(pid) { case pE: // all electrons with eta < DET->MAX_CALO_FWD DET->SmearElectron(genMomentum); if(genMomentum.E()!=0 && eta < DET->CEN_max_tracker && genMomentum.Pt() > DET->PTCUT_elec){ electron.push_back(ParticleUtil(genMomentum,particle->PID)); } break; // case pE case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD DET->SmearElectron(genMomentum); if(genMomentum.E()!=0 && eta < DET->CEN_max_tracker && genMomentum.Pt() > DET->PTCUT_gamma) { gamma.push_back(ParticleUtil(genMomentum,particle->PID)); } break; // case pGAMMA case pMU: // all muons with eta < DET->MAX_MU DET->SmearMu(genMomentum); if(genMomentum.E()!=0 && eta < DET->CEN_max_mu && genMomentum.Pt() > DET->PTCUT_muon){ muon.push_back(ParticleUtil(genMomentum,particle->PID)); } break; // case pMU case pLAMBDA: // all lambdas with eta < DET->MAX_CALO_FWD case pK0S: // all K0s with eta < DET->MAX_CALO_FWD DET->SmearHadron(genMomentum, 0.7); break; // case hadron default: // all other final particles with eta < DET->MAX_CALO_FWD DET->SmearHadron(genMomentum, 1.0); break; } // switch (pid) // all final particles but muons and neutrinos // for calorimetric towers and mission PT int charge=Charge(pid); if(genMomentum.E() !=0 && pid != pMU) { if(charge == 0 || (charge !=0 && genMomentum.Pt() >= DET->TRACK_ptmin)){ PhysicsTower CaloTower = PhysicsTower(LorentzVector(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E())); towers.push_back(CaloTower); // create a fastjet::PseudoJet with these components and put it onto // back of the input_particles vector input_particles.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E())); genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E); elementCalo = (TRootCalo*) branchCalo->NewEntry(); elementCalo->Set(genMomentumCalo); DET->BinEtaPhi(genMomentumCalo.Phi(), genMomentumCalo.Eta(), elementCalo->Phi, elementCalo->Eta); } } // all final charged particles if( (genMomentum.E()!=0) && (fabs(genMomentum.Eta()) < DET->CEN_max_tracker) && (genMomentum.Pt() > DET->TRACK_ptmin ) && // pt too small to be taken into account ((rand()%100) < DET->TRACK_eff) && (charge!=0) ) { elementTracks = (TRootTracks*) branchTracks->NewEntry(); elementTracks->Set(genMomentum); TrackCentral.push_back(genMomentum); } } // switch if(DET->FLAG_vfd==1) { VFD->ZDC(treeWriter,branchZDC,particle); VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle); } } // while DET->SortedVector(electron); for(unsigned int i=0; i < electron.size(); i++) { elementElec = (TRootElectron*) branchElectron->NewEntry(); elementElec->Set(electron[i].Px(),electron[i].Py(),electron[i].Pz(),electron[i].E()); elementElec->Charge = sign(electron[i].PID()); elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0); } DET->SortedVector(muon); for(unsigned int i=0; i < muon.size(); i++) { elementMu = (TRootMuon*) branchMuon->NewEntry(); elementMu->Charge = sign(muon[i].PID()); elementMu->Set(muon[i].Px(),muon[i].Py(),muon[i].Pz(),muon[i].E()); elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0); } DET->SortedVector(gamma); for(unsigned int i=0; i < gamma.size(); i++) { elementPhoton = (TRootPhoton*) branchPhoton->NewEntry(); elementPhoton->Set(gamma[i].Px(),gamma[i].Py(),gamma[i].Pz(),gamma[i].E()); } // computes the Missing Transverse Momentum TLorentzVector Att(0.,0.,0.,0.); for(unsigned int i=0; i < towers.size(); i++) { Att.SetPxPyPzE(towers[i].fourVector.px,towers[i].fourVector.py,towers[i].fourVector.pz,towers[i].fourVector.E); PTmis = PTmis + Att; } elementEtmis = (TRootETmis*) branchETmis->NewEntry(); elementEtmis->ET = (PTmis).Pt(); elementEtmis->Phi = (-PTmis).Phi(); elementEtmis->Px = (-PTmis).Px(); elementEtmis->Py = (-PTmis).Py(); //***************************** sorted_jets=JETRUN->RunJets(input_particles); JETRUN->RunJetBtagging(treeWriter, branchJet,sorted_jets,NFCentralQ); JETRUN->RunTauJets(treeWriter,branchTauJet,sorted_jets,towers, TrackCentral); // Add here the trigger // Should test all the trigger table on the event, based on reconstructed objects treeWriter->Fill(); } // Loop over all events treeWriter->Write(); delete treeWriter; //running the trigger in case the FLAG trigger is put to 1 in the datacard if(DET->FLAG_trigger == 1) { TChain chainT("Analysis"); chainT.Add(outputfilename.c_str()); ExRootTreeReader *treeReaderT = new ExRootTreeReader(&chainT); TClonesArray *branchElecTrig = treeReaderT->UseBranch("Electron"); TClonesArray *branchMuonTrig = treeReaderT->UseBranch("Muon"); TClonesArray *branchJetTrig = treeReaderT->UseBranch("Jet"); TClonesArray *branchTauJetTrig = treeReaderT->UseBranch("TauJet"); TClonesArray *branchPhotonTrig = treeReaderT->UseBranch("Photon"); TClonesArray *branchETmisTrig = treeReaderT->UseBranch("ETmis"); ExRootTreeWriter *treeWriterT = new ExRootTreeWriter(outputfilename, "Trigger"); ExRootTreeBranch *branchTrigger = treeWriterT->NewBranch("TrigResult", TRootTrigger::Class()); Long64_t entryT, allEntriesT = treeReaderT->GetEntries(); cout << "** Chain contains " << allEntriesT << " events" << endl; for(entryT = 0; entryT < allEntriesT; ++entryT) { treeWriterT->Clear(); treeReaderT->ReadEntry(entryT); TRIGT->GetGlobalResult(branchElecTrig, branchMuonTrig,branchJetTrig, branchTauJetTrig,branchPhotonTrig, branchETmisTrig,branchTrigger); treeWriterT->Fill(); } treeWriterT->Write(); delete treeWriterT; } //FROG display if(DET->FLAG_frog == 1) { FrogDisplay *FROG = new FrogDisplay(); FROG->BuidEvents(outputfilename,DET->NEvents_Frog); FROG->BuildGeom(DetDatacard); } cout << "** Exiting..." << endl; delete treeReader; delete DET; delete TRIGT; delete TRACP; delete JETRUN; delete VFD; if(converter) delete converter; todo("TODO"); }