/* ---- 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 #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) { cout << " Usage: " << argv[0] << " input_file" << " output_file" << " data_card " << endl; cout << " input_list - list of files in Ntpl, StdHep of LHEF format," << endl; cout << " output_file - output file." << endl; cout << " data_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); //Trigger information Trigger *TRIG = new Trigger(); TRIG->TriggerReader("data/trigger.dat"); //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 elecPID; vector muon; vector muonPID; 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(); elecPID.clear(); muonPID.clear(); NFCentralQ.Clear(); itGen.Reset(); TrackCentral.clear(); towers.clear(); input_particles.clear(); // Loop over all particles in event while( (particle = (TRootGenParticle*) itGen.Next()) ) { genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E); int pid = abs(particle->PID); float eta = fabs(genMomentum.Eta()); //subarray of the quarks (i.e. not final particles, i.e status not equal to 1) // in the tracker (i.e. for those we know the tracks) // with enough p_T //// 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 ? eta < DET->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 == pMU && eta < DET->MAX_MU) || (pid != pMU && (pid != pNU1) && (pid != pNU2) && (pid != pNU3) && eta < DET->MAX_CALO_FWD) ) ) { // TRACP->Propagation(particle,genMomentum); eta = fabs(genMomentum.Eta()); switch(pid) { case pE: // all electrons with eta < DET->MAX_CALO_FWD DET->SmearElectron(genMomentum); electron.push_back(genMomentum); elecPID.push_back(particle->PID); break; // case pE case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD DET->SmearElectron(genMomentum); if(genMomentum.E()!=0 && eta < DET->MAX_TRACKER) { elementPhoton = (TRootPhoton*) branchPhoton->NewEntry(); elementPhoton->Set(genMomentum); } break; // case pGAMMA case pMU: // all muons with eta < DET->MAX_MU DET->SmearMu(genMomentum); muonPID.push_back(particle->PID); muon.push_back(genMomentum); 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 if(genMomentum.E() !=0) { if(pid !=pMU) { 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.SetPtEtaPhiE(CaloTower.Et(),CaloTower.iEta(),CaloTower.iPhi(),CaloTower.fourVector.E); genMomentumCalo.SetPxPyPzE(CaloTower.fourVector.px,CaloTower.fourVector.py,CaloTower.fourVector.pz,CaloTower.fourVector.E); elementCalo = (TRootCalo*) branchCalo->NewEntry(); elementCalo->Set(genMomentumCalo); } } // all final charged particles if( ((rand()%100) < DET->TRACKING_EFF) && (genMomentum.E()!=0) && (fabs(genMomentum.Eta()) < DET->MAX_TRACKER) && (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account (pid != pGAMMA) && (pid != pPI0) && (pid != pK0L) && (pid != pN) && (pid != pSIGMA0) && (pid != pDELTA0) && (pid != pK0S) // not charged particles : invisible by tracker ) { elementTracks = (TRootTracks*) branchTracks->NewEntry(); elementTracks->Set(genMomentum); TrackCentral.push_back(genMomentum); } } // switch VFD->ZDC(treeWriter,branchZDC,particle); VFD->RomanPots(treeWriter,branchRP220,branchFP420,particle); } // while for(unsigned int i=0; i < electron.size(); i++) { if(electron[i].E()!=0 && fabs(electron[i].Eta()) < DET->MAX_TRACKER && electron[i].Pt() > DET->ELEC_pt) { elementElec = (TRootElectron*) branchElectron->NewEntry(); elementElec->Set(electron[i]); elementElec->Charge = sign(elecPID[i]); elementElec->IsolFlag = DET->Isolation(electron[i].Phi(),electron[i].Eta(),TrackCentral,2.0); } } for(unsigned int i=0; i < muon.size(); i++) { if(muon[i].E()!=0 && fabs(muon[i].Eta()) < DET->MAX_MU && muon[i].Pt() > DET->MUON_pt) { elementMu = (TRootMuon*) branchMuon->NewEntry(); elementMu->Charge = sign(muonPID[i]); elementMu->Set(muon[i]); elementMu->IsolFlag = DET->Isolation(muon[i].Phi(),muon[i].Eta(),TrackCentral,2.0); } } // 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); treeWriter->Fill(); // Add here the trigger // Should test all the trigger table on the event, based on reconstructed objects // Trigger *TRIG = new Trigger(); // TRIG->TriggerReader("data/trigger.dat"); } // Loop over all events treeWriter->Write(); cout << "** Exiting..." << endl; delete treeWriter; delete treeReader; delete DET; if(converter) delete converter; todo("TODO"); }