#include <SimpleConvertedPhotonAnalyzer.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
virtual void | endJob () |
SimpleConvertedPhotonAnalyzer (const edm::ParameterSet &) | |
virtual | ~SimpleConvertedPhotonAnalyzer () |
Private Member Functions | |
float | etaTransformation (float a, float b) |
Private Attributes | |
std::string | convertedPhotonCollection_ |
std::string | convertedPhotonCollectionProducer_ |
TFile * | fOutputFile_ |
std::string | fOutputFileName_ |
TH2F * | h2_tk_inPtVsR_ |
TH2F * | h2_tk_nHitsVsR_ |
TH1F * | h_deltaEta_ |
TH1F * | h_deltaPhi_ |
TH1F * | h_ErecoEMC_ |
TH1F * | h_MCConvE_ |
TH1F * | h_MCConvEta_ |
TH1F * | h_MCConvPt_ |
TH1F * | h_MCphoE_ |
TH1F * | h_MCphoEta_ |
TH1F * | h_MCphoPhi_ |
TH1F * | h_phoE_ |
TH1F * | h_phoEta_ |
TH1F * | h_phoPhi_ |
TH1F * | h_scE_ |
TH1F * | h_scEta_ |
TH1F * | h_scPhi_ |
std::string | HepMCLabel |
int | nEvt_ |
int | nMatched_ |
int | nMCPho_ |
std::string | SimHitLabel |
std::string | SimTkLabel |
std::string | SimVtxLabel |
PhotonMCTruthFinder * | thePhotonMCTruthFinder_ |
Definition at line 23 of file SimpleConvertedPhotonAnalyzer.h.
SimpleConvertedPhotonAnalyzer::SimpleConvertedPhotonAnalyzer | ( | const edm::ParameterSet & | pset | ) | [explicit] |
Definition at line 47 of file SimpleConvertedPhotonAnalyzer.cc.
References convertedPhotonCollection_, convertedPhotonCollectionProducer_, and edm::ParameterSet::getParameter().
: fOutputFileName_( pset.getUntrackedParameter<string>("HistOutFile",std::string("TestConversions.root")) ), fOutputFile_(0) { convertedPhotonCollectionProducer_ = pset.getParameter<std::string>("phoProducer"); convertedPhotonCollection_ = pset.getParameter<std::string>("convertedPhotonCollection"); // }
SimpleConvertedPhotonAnalyzer::~SimpleConvertedPhotonAnalyzer | ( | ) | [virtual] |
Definition at line 61 of file SimpleConvertedPhotonAnalyzer.cc.
References thePhotonMCTruthFinder_.
{ delete thePhotonMCTruthFinder_; }
void SimpleConvertedPhotonAnalyzer::analyze | ( | const edm::Event & | e, |
const edm::EventSetup & | |||
) | [virtual] |
Loop over recontructed photons
End loop over Reco particles
End loop over MC particles
Implements edm::EDAnalyzer.
Definition at line 156 of file SimpleConvertedPhotonAnalyzer.cc.
References convertedPhotonCollection_, convertedPhotonCollectionProducer_, gather_cfg::cout, delta, HLTFastRecoForTau_cff::deltaEta, SiPixelRawToDigiRegional_cfi::deltaPhi, etaTransformation(), HcalObjRepresent::Fill(), PhotonMCTruthFinder::find(), edm::Event::getByLabel(), h2_tk_inPtVsR_, h2_tk_nHitsVsR_, h_deltaEta_, h_deltaPhi_, h_ErecoEMC_, h_MCConvEta_, h_MCphoE_, h_MCphoEta_, h_MCphoPhi_, h_scE_, h_scEta_, h_scPhi_, i, edm::EventBase::id(), nEvt_, nMatched_, nMCPho_, Geom::pi(), funct::pow(), mathSSE::sqrt(), thePhotonMCTruthFinder_, and Geom::twoPi().
{ using namespace edm; const float etaPhiDistance=0.01; // Fiducial region //UNUSED const float TRK_BARL =0.9; const float BARL = 1.4442; // DAQ TDR p.290 const float END_LO = 1.566; const float END_HI = 2.5; // Electron mass //UNUSED const Float_t mElec= 0.000511; nEvt_++; LogInfo("ConvertedPhotonAnalyzer") << "ConvertedPhotonAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n"; // LogDebug("ConvertedPhotonAnalyzer") << "ConvertedPhotonAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n"; std::cout << "ConvertedPhotonAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n"; Handle<reco::ConversionCollection> convertedPhotonHandle; e.getByLabel(convertedPhotonCollectionProducer_, convertedPhotonCollection_ , convertedPhotonHandle); const reco::ConversionCollection phoCollection = *(convertedPhotonHandle.product()); std::cout << "ConvertedPhotonAnalyzer Converted photon collection size " << phoCollection.size() << "\n"; std::cout << " ConvertedPhotonAnalyzer Looking for MC truth " << "\n"; //get simtrack info std::vector<SimTrack> theSimTracks; std::vector<SimVertex> theSimVertices; edm::Handle<SimTrackContainer> SimTk; edm::Handle<SimVertexContainer> SimVtx; e.getByLabel("g4SimHits",SimTk); e.getByLabel("g4SimHits",SimVtx); theSimTracks.insert(theSimTracks.end(),SimTk->begin(),SimTk->end()); theSimVertices.insert(theSimVertices.end(),SimVtx->begin(),SimVtx->end()); std::vector<PhotonMCTruth> mcPhotons=thePhotonMCTruthFinder_->find (theSimTracks, theSimVertices); std::cout << " ConvertedPhotonAnalyzer mcPhotons size " << mcPhotons.size() << std::endl; // Loop over simulated photons //UNUSED int iDet=0; //UNUSED int iRadius=-1; //UNUSED int indPho=0; for ( std::vector<PhotonMCTruth>::const_iterator mcPho=mcPhotons.begin(); mcPho !=mcPhotons.end(); mcPho++) { float mcPhi= (*mcPho).fourMomentum().phi(); float mcEta= (*mcPho).fourMomentum().pseudoRapidity(); mcEta = etaTransformation(mcEta, (*mcPho).primaryVertex().z() ); if ( ! ( fabs(mcEta) <= BARL || ( fabs(mcEta) >= END_LO && fabs(mcEta) <=END_HI ) ) ) { continue; } // all ecal fiducial region std::cout << " ConvertedPhotonAnalyzer MC Photons before matching " << std::endl; std::cout << " ConvertedPhotonAnalyzer Photons isAconversion " << (*mcPho).isAConversion() << " mcMatchingPhoton energy " << (*mcPho).fourMomentum().e() << " conversion vertex R " << (*mcPho).vertex().perp() << " Z " << (*mcPho).vertex().z() << " x " << (*mcPho).vertex().x() << " y " <<(*mcPho).vertex().y() << " z " << (*mcPho).vertex().z() << std::endl; std::cout << " ConvertedPhotonAnalyzer mcEta " << mcEta<< " mcPhi " << mcPhi << std::endl; h_MCphoE_->Fill( (*mcPho).fourMomentum().e()); h_MCphoEta_->Fill( (*mcPho).fourMomentum().eta()); h_MCphoPhi_->Fill( (*mcPho).fourMomentum().phi()); // keep only visible conversions if ( (*mcPho).isAConversion() == 0 ) continue; nMCPho_++; h_MCConvEta_ ->Fill ( fabs( (*mcPho).fourMomentum().pseudoRapidity()) - 0.001); bool REJECTED; std::cout << " ConvertedPhotonAnalyzer Starting loop over photon candidates " << "\n"; for( reco::ConversionCollection::const_iterator iPho = phoCollection.begin(); iPho != phoCollection.end(); iPho++) { REJECTED=false; std::cout << " ConvertedPhotonAnalyzer Reco SC energy " << (*iPho).caloCluster()[0]->energy() << "\n"; float phiClu=(*iPho).caloCluster()[0]->phi(); float etaClu=(*iPho).caloCluster()[0]->eta(); float deltaPhi = phiClu-mcPhi; float deltaEta = etaClu-mcEta; if ( deltaPhi > Geom::pi() ) deltaPhi -= Geom::twoPi(); if ( deltaPhi < -Geom::pi() ) deltaPhi += Geom::twoPi(); deltaPhi=std::pow(deltaPhi,2); deltaEta=std::pow(deltaEta,2); float delta = deltaPhi+deltaEta ; if ( delta >= etaPhiDistance ) REJECTED=true; // if ( ! ( fabs(etaClu) <= BARL || ( fabs(etaClu) >= END_LO && fabs(etaClu) <=END_HI ) ) ) REJECTED=true; if ( REJECTED ) continue; std::cout << " MATCHED " << std::endl; nMatched_++; std::cout << " ConvertedPhotonAnalyzer Matching candidate " << std::endl; std::cout << " ConvertedPhotonAnalyzer Photons isAconversion " << (*mcPho).isAConversion() << " mcMatchingPhoton energy " << (*mcPho).fourMomentum().e() << " ConvertedPhotonAnalyzer conversion vertex R " << (*mcPho).vertex().perp() << " Z " << (*mcPho).vertex().z() << std::endl; h_ErecoEMC_->Fill( (*iPho).caloCluster()[0]->energy()/(*mcPho).fourMomentum().e()); h_deltaPhi_-> Fill ( (*iPho).caloCluster()[0]->position().phi()- mcPhi); h_deltaEta_-> Fill ( (*iPho).caloCluster()[0]->position().eta()- mcEta); h_scE_->Fill( (*iPho).caloCluster()[0]->energy() ); h_scEta_->Fill( (*iPho).caloCluster()[0]->position().eta() ); h_scPhi_->Fill( (*iPho).caloCluster()[0]->position().phi() ); for (unsigned int i=0; i<(*iPho).tracks().size(); i++) { std::cout << " ConvertedPhotonAnalyzer Reco Track charge " << (*iPho).tracks()[i]->charge() << " Num of RecHits " << (*iPho).tracks()[i]->recHitsSize() << " inner momentum " << sqrt ( (*iPho).tracks()[i]->innerMomentum().Mag2() ) << "\n"; h2_tk_nHitsVsR_ -> Fill ( (*mcPho).vertex().perp(), (*iPho).tracks()[i]->recHitsSize() ); h2_tk_inPtVsR_->Fill ( (*mcPho).vertex().perp(), sqrt( (*iPho).tracks()[i]->innerMomentum().Mag2() ) ); } } } }
void SimpleConvertedPhotonAnalyzer::beginJob | ( | void | ) | [virtual] |
Reco - MC
Reimplemented from edm::EDAnalyzer.
Definition at line 69 of file SimpleConvertedPhotonAnalyzer.cc.
References fOutputFile_, fOutputFileName_, h2_tk_inPtVsR_, h2_tk_nHitsVsR_, h_deltaEta_, h_deltaPhi_, h_ErecoEMC_, h_MCConvE_, h_MCConvEta_, h_MCConvPt_, h_MCphoE_, h_MCphoEta_, h_MCphoPhi_, h_phoE_, h_phoEta_, h_phoPhi_, h_scE_, h_scEta_, h_scPhi_, nEvt_, nMatched_, nMCPho_, and thePhotonMCTruthFinder_.
{ nEvt_=0; nMCPho_=0; nMatched_=0; thePhotonMCTruthFinder_ = new PhotonMCTruthFinder(); fOutputFile_ = new TFile( fOutputFileName_.c_str(), "RECREATE" ) ; h_ErecoEMC_ = new TH1F("deltaE"," delta(reco-mc) energy",100,0.,2.); h_deltaPhi_ = new TH1F("deltaPhi"," delta(reco-mc) phi",100,-0.1, 0.1); h_deltaEta_ = new TH1F("deltaEta"," delta(reco-mc) eta",100,-0.05, 0.05); h_MCphoE_ = new TH1F("MCphoE","MC photon energy",100,0.,100.); h_MCphoPhi_ = new TH1F("MCphoPhi","MC photon phi",40,-3.14, 3.14); h_MCphoEta_ = new TH1F("MCphoEta","MC photon eta",40,-3., 3.); h_MCConvE_ = new TH1F("MCConvE","MC converted photon energy",100,0.,100.); h_MCConvPt_ = new TH1F("MCConvPt","MC converted photon pt",100,0.,100.); h_MCConvEta_ = new TH1F("MCConvEta","MC converted photon eta",50, 0., 2.5); h_scE_ = new TH1F("scE","Uncorrected converted photons : SC Energy ",100,0., 200.); h_scEta_ = new TH1F("scEta","Uncorrected converted photons: SC Eta ",40,-3., 3.); h_scPhi_ = new TH1F("scPhi","Uncorrected converted photons: SC Phi ",40, -3.14, 3.14); // h_phoE_ = new TH1F("phoE","Uncorrected converted photons : Energy ",100,0., 200.); h_phoEta_ = new TH1F("phoEta","Uncorrected converted photons: Eta ",40,-3., 3.); h_phoPhi_ = new TH1F("phoPhi","Uncorrected converted photons: Phi ",40, -3.14, 3.14); // Recontructed tracks from converted photon candidates h2_tk_nHitsVsR_ = new TH2F("tknHitsVsR","Tracks Hits vs R ", 12,0.,120.,20,0.5, 20.5); h2_tk_inPtVsR_ = new TH2F("tkInPtvsR","Tracks inner Pt vs R ", 12,0.,120.,100,0., 100.); return ; }
void SimpleConvertedPhotonAnalyzer::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 314 of file SimpleConvertedPhotonAnalyzer.cc.
References gather_cfg::cout, fOutputFile_, and nEvt_.
{ fOutputFile_->Write() ; fOutputFile_->Close() ; edm::LogInfo("ConvertedPhotonAnalyzer") << "Analyzed " << nEvt_ << "\n"; // std::cout << "::endJob Analyzed " << nEvt_ << " events " << " with total " << nPho_ << " Photons " << "\n"; std::cout << "ConvertedPhotonAnalyzer::endJob Analyzed " << nEvt_ << " events " << "\n"; return ; }
float SimpleConvertedPhotonAnalyzer::etaTransformation | ( | float | a, |
float | b | ||
) | [private] |
Definition at line 116 of file SimpleConvertedPhotonAnalyzer.cc.
References ETA, etaBarrelEndcap, funct::log(), PI, R_ECAL, funct::tan(), and Z_Endcap.
Referenced by analyze().
{ //---Definitions const float PI = 3.1415927; //UNUSED const float TWOPI = 2.0*PI; //---Definitions for ECAL const float R_ECAL = 136.5; const float Z_Endcap = 328.0; const float etaBarrelEndcap = 1.479; //---ETA correction float Theta = 0.0 ; float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex; if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal); if(Theta<0.0) Theta = Theta+PI ; float ETA = - log(tan(0.5*Theta)); if( fabs(ETA) > etaBarrelEndcap ) { float Zend = Z_Endcap ; if(EtaParticle<0.0 ) Zend = -Zend ; float Zlen = Zend - Zvertex ; float RR = Zlen/sinh(EtaParticle); Theta = atan(RR/Zend); if(Theta<0.0) Theta = Theta+PI ; ETA = - log(tan(0.5*Theta)); } //---Return the result return ETA; //---end }
std::string SimpleConvertedPhotonAnalyzer::convertedPhotonCollection_ [private] |
Definition at line 60 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and SimpleConvertedPhotonAnalyzer().
std::string SimpleConvertedPhotonAnalyzer::convertedPhotonCollectionProducer_ [private] |
Definition at line 59 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and SimpleConvertedPhotonAnalyzer().
TFile* SimpleConvertedPhotonAnalyzer::fOutputFile_ [private] |
Definition at line 46 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by beginJob(), and endJob().
std::string SimpleConvertedPhotonAnalyzer::fOutputFileName_ [private] |
Definition at line 45 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by beginJob().
TH2F* SimpleConvertedPhotonAnalyzer::h2_tk_inPtVsR_ [private] |
Definition at line 92 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH2F* SimpleConvertedPhotonAnalyzer::h2_tk_nHitsVsR_ [private] |
Definition at line 90 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_deltaEta_ [private] |
Definition at line 64 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_deltaPhi_ [private] |
Definition at line 63 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_ErecoEMC_ [private] |
Definition at line 62 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_MCConvE_ [private] |
Definition at line 75 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_MCConvEta_ [private] |
Definition at line 77 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_MCConvPt_ [private] |
Definition at line 76 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_MCphoE_ [private] |
Definition at line 68 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_MCphoEta_ [private] |
Definition at line 70 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_MCphoPhi_ [private] |
Definition at line 69 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_phoE_ [private] |
Definition at line 85 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_phoEta_ [private] |
Definition at line 86 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_phoPhi_ [private] |
Definition at line 87 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_scE_ [private] |
Definition at line 81 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_scEta_ [private] |
Definition at line 82 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* SimpleConvertedPhotonAnalyzer::h_scPhi_ [private] |
Definition at line 83 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
std::string SimpleConvertedPhotonAnalyzer::HepMCLabel [private] |
Definition at line 53 of file SimpleConvertedPhotonAnalyzer.h.
int SimpleConvertedPhotonAnalyzer::nEvt_ [private] |
Definition at line 49 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
int SimpleConvertedPhotonAnalyzer::nMatched_ [private] |
Definition at line 51 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
int SimpleConvertedPhotonAnalyzer::nMCPho_ [private] |
Definition at line 50 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), and beginJob().
std::string SimpleConvertedPhotonAnalyzer::SimHitLabel [private] |
Definition at line 56 of file SimpleConvertedPhotonAnalyzer.h.
std::string SimpleConvertedPhotonAnalyzer::SimTkLabel [private] |
Definition at line 54 of file SimpleConvertedPhotonAnalyzer.h.
std::string SimpleConvertedPhotonAnalyzer::SimVtxLabel [private] |
Definition at line 55 of file SimpleConvertedPhotonAnalyzer.h.
Definition at line 43 of file SimpleConvertedPhotonAnalyzer.h.
Referenced by analyze(), beginJob(), and ~SimpleConvertedPhotonAnalyzer().