/* ---- FastSim ---- 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 Smearing.cpp /// \brief executable for the FastSim #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 "H_BeamParticle.h" #include "H_BeamLine.h" #include "H_RomanPot.h" #include "interface/DataConverter.h" #include "interface/HEPEVTConverter.h" #include "interface/LHEFConverter.h" #include "interface/STDHEPConverter.h" #include "interface/SmearUtil.h" #include "interface/JetUtils.h" #include "interface/BFieldProp.h" #include "Utilities/Fastjet/include/fastjet/PseudoJet.hh" #include "Utilities/Fastjet/include/fastjet/ClusterSequence.hh" #include #include #include "interface/TreeClasses.h" using namespace std; //------------------------------------------------------------------------------ // //********************************** PYTHIA INFORMATION********************************* TSimpleArray TauHadr(const TClonesArray *GEN) { TIter it((TCollection*)GEN); it.Reset(); TRootGenParticle *gen1; TSimpleArray array,array2; while((gen1 = (TRootGenParticle*) it.Next())) { array.Add(gen1); } it.Reset(); bool tauhad; while((gen1 = (TRootGenParticle*) it.Next())) { tauhad=false; if(abs(gen1->PID)==15) { int d1=gen1->D1; int d2=gen1->D2; if((d1 < array.GetEntries()) && (d1 > 0) && (d2 < array.GetEntries()) && (d2 > 0)) { tauhad=true; for(int d=d1; d < d2+1; d++) { if(abs(array[d]->PID)== pE || abs(array[d]->PID)== pMU)tauhad=false; } } } if(tauhad)array2.Add(gen1); } return array2; } void PairingJet(TLorentzVector &JETSm, const TLorentzVector &JET, vector sorted_jetsS) { JETSm.SetPxPyPzE(0,0,0,0); float deltaRtest=5000; for (unsigned int i = 0; i < sorted_jetsS.size(); i++) { TLorentzVector Att; Att.SetPxPyPzE(sorted_jetsS[i].px(),sorted_jetsS[i].py(),sorted_jetsS[i].pz(),sorted_jetsS[i].E()); if(DeltaR(JET.Phi(),JET.Eta(),Att.Phi(),Att.Eta()) < deltaRtest) { deltaRtest = DeltaR(JET.Phi(),JET.Eta(),Att.Phi(),Att.Eta()); if(deltaRtest < 0.25) { JETSm = Att; } } } } int main(int argc, char *argv[]) { int appargc = 2; char *appName = "Smearing"; 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"; return -1; } 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 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("JetPTResol", RESOLJET::Class()); ExRootTreeBranch *branchelec = treeWriter->NewBranch("ElecEResol", RESOLELEC::Class()); ExRootTreeBranch *branchmuon = treeWriter->NewBranch("MuonPTResol", RESOLMUON::Class()); ExRootTreeBranch *branchtaujet = treeWriter->NewBranch("TauJetPTResol", TAUHAD::Class()); ExRootTreeBranch *branchetmis = treeWriter->NewBranch("ETmisResol",ETMIS::Class()); TRootGenParticle *particle; RESOLELEC * elementElec; RESOLMUON *elementMuon; RESOLJET *elementJet; TAUHAD *elementTaujet; ETMIS *elementEtmis; //read the datacard input file string DetDatacard("data/DataCardDet.dat"); if(argc==4) DetDatacard =argv[3]; RESOLution *DET = new RESOLution(); DET->ReadDataCard(DetDatacard); //Jet information JetsUtil *JETRUN = new JetsUtil(); //Propagation of tracks in the B field TrackPropagation *TRACP = new TrackPropagation(); TLorentzVector genMomentum(0,0,0,0);//TLorentzVector containing generator level information TLorentzVector recoMomentum(0,0,0,0);//TLorentzVector containing Reco level information LorentzVector jetMomentum; vector TrackCentral; vector input_particlesGEN;//for FastJet algorithm vector sorted_jetsGEN; vector input_particlesReco;//for FastJet algorithm vector sorted_jetsReco; vector towers; // Loop over all events Long64_t entry, allEntries = treeReader->GetEntries(); cout << "** Chain contains " << allEntries << " events" << endl; for(entry = 0; entry < allEntries; ++entry) { TLorentzVector PTmisReco(0,0,0,0); TLorentzVector PTmisGEN(0,0,0,0); treeReader->ReadEntry(entry); treeWriter->Clear(); if((entry % 100) == 0 && entry > 0 ) cout << "** Processing element # " << entry << endl; TSimpleArray bGen; itGen.Reset(); TrackCentral.clear(); TSimpleArray NFCentralQ; input_particlesReco.clear(); input_particlesGEN.clear(); towers.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(particle->Eta); //input generator level particle for jet algorithm if(particle->Status == 1 && eta < DET->MAX_CALO_FWD) { input_particlesGEN.push_back(fastjet::PseudoJet(genMomentum.Px(),genMomentum.Py(),genMomentum.Pz(), genMomentum.E())); } //Calculate ETMIS from generated particles if((pid == pNU1) || (pid == pNU2) || (pid == pNU3))PTmisGEN = PTmisGEN + genMomentum; if( (particle->Status == 1) && ((pid != pNU1) && (pid != pNU2) && (pid != pNU3)) && (fabs(particle->Eta) < DET->MAX_CALO_FWD) ) { recoMomentum = genMomentum; //use of the magnetic field propagation TRACP->Propagation(particle,recoMomentum); // cout<<"eta "<MAX_CALO_FWD DET->SmearElectron(recoMomentum); if(recoMomentum.E() !=0 && eta < DET->MAX_TRACKER){ elementElec=(RESOLELEC*) branchelec->NewEntry(); elementElec->E = genMomentum.E(); elementElec->SmearedE = recoMomentum.E();} break; // case pE case pGAMMA: // all photons with eta < DET->MAX_CALO_FWD DET->SmearElectron(recoMomentum); break; // case pGAMMA case pMU: // all muons with eta < DET->MAX_MU DET->SmearMu(recoMomentum); if(recoMomentum.E() !=0){ elementMuon = (RESOLMUON*) branchmuon->NewEntry(); elementMuon->OverPT = (1/genMomentum.Pt()); elementMuon->OverSmearedPT = (1/recoMomentum.Pt());} 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(recoMomentum, 0.7); break; // case hadron default: // all other final particles with eta < DET->MAX_CALO_FWD DET->SmearHadron(recoMomentum, 1.0); break; } // switch (pid) //information to reconstruct jets from reconstructed towers int charge=Charge(pid); if(recoMomentum.E() !=0 && pid != pMU) { if(charge == 0 || (charge !=0 && recoMomentum.Pt() >= DET->PT_TRACKS_MIN)){ towers.push_back(PhysicsTower(LorentzVector(recoMomentum.Px(),recoMomentum.Py(),recoMomentum.Pz(), recoMomentum.E()))); input_particlesReco.push_back(fastjet::PseudoJet(recoMomentum.Px(),recoMomentum.Py(),recoMomentum.Pz(), recoMomentum.E())); } } // all final charged particles if( (genMomentum.E()!=0) && (fabs(genMomentum.Eta()) < DET->MAX_TRACKER) && (genMomentum.Pt() > DET->PT_TRACKS_MIN ) && // pt too small to be taken into account ((rand()%100) < DET->TRACKING_EFF) && (charge!=0) ) {TrackCentral.push_back(genMomentum);} } // switch } // while //compute missing transverse energy from calo towers TLorentzVector Att(0.,0.,0.,0.); float ScalarEt=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); ScalarEt = ScalarEt + Att.Et(); PTmisReco = PTmisReco + Att; } elementEtmis= (ETMIS*) branchetmis->NewEntry(); elementEtmis->Et = (PTmisGEN).Pt(); elementEtmis->Ex = (-PTmisGEN).Px(); elementEtmis->SEt = ScalarEt; elementEtmis->EtSmeare = (PTmisReco).Pt()-(PTmisGEN).Pt(); elementEtmis->ExSmeare = (-PTmisReco).Px()-(PTmisGEN).Px(); //***************************** sorted_jetsGEN=JETRUN->RunJets(input_particlesGEN); sorted_jetsReco=JETRUN->RunJets(input_particlesReco); TSimpleArray TausHadr = TauHadr(branchGen); TLorentzVector JETreco(0,0,0,0); for (unsigned int i = 0; i < sorted_jetsGEN.size(); i++) { TLorentzVector JETgen(0,0,0,0); JETgen.SetPxPyPzE(sorted_jetsGEN[i].px(),sorted_jetsGEN[i].py(),sorted_jetsGEN[i].pz(),sorted_jetsGEN[i].E()); PairingJet(JETreco,JETgen,sorted_jetsReco); if(JETreco.Pt()>1) { elementJet= (RESOLJET*) branchjet->NewEntry(); elementJet->PT = JETgen.Et(); elementJet->SmearedPT = JETreco.Et()/JETgen.Et(); } } for (unsigned int i = 0; i < sorted_jetsReco.size(); i++) { TLorentzVector JETT(0,0,0,0); JETT.SetPxPyPzE(sorted_jetsReco[i].px(),sorted_jetsReco[i].py(),sorted_jetsReco[i].pz(),sorted_jetsReco[i].E()); if(fabs(JETT.Eta()) < (DET->MAX_TRACKER - DET->TAU_CONE_TRACKS)) { for(Int_t i=0; iPhi,TausHadr[i]->Eta,JETT.Phi(),JETT.Eta())<0.1) { elementTaujet= (TAUHAD*) branchtaujet->NewEntry(); elementTaujet->EnergieCen = (DET->EnergySmallCone(towers,JETT.Eta(),JETT.Phi())/JETT.E()); elementTaujet->NumTrack = DET->NumTracks(TrackCentral,DET->PT_TRACK_TAU,JETT.Eta(),JETT.Phi()); } } } } // for itJet : loop on all jets treeWriter->Fill(); } // Loop over all events treeWriter->Write(); cout << "** Exiting..." << endl; delete treeWriter; delete treeReader; delete DET; if(converter) delete converter; }