#include <MCElectronAnalyzer.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
virtual void | endJob () |
MCElectronAnalyzer (const edm::ParameterSet &) | |
virtual | ~MCElectronAnalyzer () |
Private Member Functions | |
float | etaTransformation (float a, float b) |
float | phiNormalization (float &a) |
Private Attributes | |
TFile * | fOutputFile_ |
std::string | fOutputFileName_ |
TH1F * | h_BremEnergy_ |
TH1F * | h_BremFrac_ |
TH1F * | h_MCEleE_ |
TH1F * | h_MCEleEta_ |
TH1F * | h_MCElePhi_ |
std::string | HepMCLabel |
double | mcEta_ |
double | mcPhi_ |
global variable for the MC photon | |
int | nEvt_ |
int | nMatched_ |
TProfile * | p_BremVsEta_ |
TProfile * | p_BremVsR_ |
std::string | SimHitLabel |
std::string | SimTkLabel |
std::string | SimVtxLabel |
ElectronMCTruthFinder * | theElectronMCTruthFinder_ |
const TrackerGeometry * | trackerGeom |
Definition at line 23 of file MCElectronAnalyzer.h.
MCElectronAnalyzer::MCElectronAnalyzer | ( | const edm::ParameterSet & | pset | ) | [explicit] |
Definition at line 44 of file MCElectronAnalyzer.cc.
: fOutputFileName_( pset.getUntrackedParameter<string>("HistOutFile",std::string("TestConversions.root")) ), fOutputFile_(0) { }
MCElectronAnalyzer::~MCElectronAnalyzer | ( | ) | [virtual] |
Definition at line 54 of file MCElectronAnalyzer.cc.
References theElectronMCTruthFinder_.
{ delete theElectronMCTruthFinder_; }
void MCElectronAnalyzer::analyze | ( | const edm::Event & | e, |
const edm::EventSetup & | |||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 144 of file MCElectronAnalyzer.cc.
References gather_cfg::cout, ElectronMCTruthFinder::find(), edm::Event::getByLabel(), h_BremEnergy_, h_BremFrac_, h_MCEleE_, h_MCEleEta_, h_MCElePhi_, edm::EventBase::id(), nEvt_, p_BremVsEta_, p_BremVsR_, and theElectronMCTruthFinder_.
{ using namespace edm; //const float etaPhiDistance=0.01; // Fiducial region //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 //const Float_t mElec= 0.000511; nEvt_++; LogInfo("MCElectronAnalyzer") << "MCElectronAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n"; // LogDebug("MCElectronAnalyzer") << "MCElectronAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n"; std::cout << "MCElectronAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n"; std::cout << " MCElectronAnalyzer 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::cout << " MCElectronAnalyzer This Event has " << theSimTracks.size() << " sim tracks " << std::endl; std::cout << " MCElectronAnalyzer This Event has " << theSimVertices.size() << " sim vertices " << std::endl; if ( ! theSimTracks.size() ) std::cout << " Event number " << e.id() << " has NO sim tracks " << std::endl; std::vector<ElectronMCTruth> MCElectronctrons=theElectronMCTruthFinder_->find (theSimTracks, theSimVertices); std::cout << " MCElectronAnalyzer MCElectronctrons size " << MCElectronctrons.size() << std::endl; for ( std::vector<ElectronMCTruth>::const_iterator iEl=MCElectronctrons.begin(); iEl !=MCElectronctrons.end(); ++iEl ){ h_MCEleE_->Fill ( (*iEl).fourMomentum().e() ); h_MCEleEta_->Fill ( (*iEl).fourMomentum().pseudoRapidity() ); h_MCElePhi_->Fill ( (*iEl).fourMomentum().phi() ); float totBrem=0 ; unsigned int iBrem ; for ( iBrem=0; iBrem < (*iEl).bremVertices().size(); ++iBrem ) { float rBrem= (*iEl).bremVertices()[iBrem].perp(); float etaBrem= (*iEl).bremVertices()[iBrem].eta(); if ( rBrem < 120 ) { totBrem += (*iEl).bremMomentum()[iBrem].e(); p_BremVsR_ ->Fill ( rBrem, (*iEl).bremMomentum()[iBrem].e() ); p_BremVsEta_ ->Fill ( etaBrem, (*iEl).bremMomentum()[iBrem].e() ); } } h_BremFrac_->Fill( totBrem/(*iEl).fourMomentum().e() ); h_BremEnergy_->Fill ( totBrem ); } }
void MCElectronAnalyzer::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 62 of file MCElectronAnalyzer.cc.
References fOutputFile_, fOutputFileName_, h_BremEnergy_, h_BremFrac_, h_MCEleE_, h_MCEleEta_, h_MCElePhi_, nEvt_, p_BremVsEta_, p_BremVsR_, and theElectronMCTruthFinder_.
{ nEvt_=0; theElectronMCTruthFinder_ = new ElectronMCTruthFinder(); fOutputFile_ = new TFile( fOutputFileName_.c_str(), "RECREATE" ) ; h_MCEleE_ = new TH1F("MCEleE","MC ele energy",100,0.,200.); h_MCElePhi_ = new TH1F("MCElePhi","MC ele phi",40,-3.14, 3.14); h_MCEleEta_ = new TH1F("MCEleEta","MC ele eta",40,-3., 3.); h_BremFrac_ = new TH1F("bremFrac","brem frac ", 100, 0., 1.); h_BremEnergy_ = new TH1F("BremE","Brem energy",100,0.,200.); p_BremVsR_ = new TProfile("BremVsR", " Mean Brem energy vs R ", 48, 0., 120.); p_BremVsEta_ = new TProfile("BremVsEta", " Mean Brem energy vs Eta ", 50, -2.5, 2.5); return ; }
void MCElectronAnalyzer::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 222 of file MCElectronAnalyzer.cc.
References gather_cfg::cout, fOutputFile_, and nEvt_.
{ fOutputFile_->Write() ; fOutputFile_->Close() ; edm::LogInfo("MCElectronAnalyzer") << "Analyzed " << nEvt_ << "\n"; std::cout << "MCElectronAnalyzer::endJob Analyzed " << nEvt_ << " events " << "\n"; return ; }
float MCElectronAnalyzer::etaTransformation | ( | float | a, |
float | b | ||
) | [private] |
Definition at line 90 of file MCElectronAnalyzer.cc.
References abs, ETA, etaBarrelEndcap, create_public_lumi_plots::log, PI, R_ECAL, funct::tan(), and Z_Endcap.
{ //---Definitions const float PI = 3.1415927; //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( std::abs(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 }
float MCElectronAnalyzer::phiNormalization | ( | float & | a | ) | [private] |
Definition at line 125 of file MCElectronAnalyzer.cc.
TFile* MCElectronAnalyzer::fOutputFile_ [private] |
Definition at line 50 of file MCElectronAnalyzer.h.
Referenced by beginJob(), and endJob().
std::string MCElectronAnalyzer::fOutputFileName_ [private] |
Definition at line 49 of file MCElectronAnalyzer.h.
Referenced by beginJob().
TH1F* MCElectronAnalyzer::h_BremEnergy_ [private] |
Definition at line 72 of file MCElectronAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* MCElectronAnalyzer::h_BremFrac_ [private] |
Definition at line 71 of file MCElectronAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* MCElectronAnalyzer::h_MCEleE_ [private] |
Definition at line 68 of file MCElectronAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* MCElectronAnalyzer::h_MCEleEta_ [private] |
Definition at line 69 of file MCElectronAnalyzer.h.
Referenced by analyze(), and beginJob().
TH1F* MCElectronAnalyzer::h_MCElePhi_ [private] |
Definition at line 70 of file MCElectronAnalyzer.h.
Referenced by analyze(), and beginJob().
std::string MCElectronAnalyzer::HepMCLabel [private] |
Definition at line 62 of file MCElectronAnalyzer.h.
double MCElectronAnalyzer::mcEta_ [private] |
Definition at line 60 of file MCElectronAnalyzer.h.
double MCElectronAnalyzer::mcPhi_ [private] |
global variable for the MC photon
Definition at line 59 of file MCElectronAnalyzer.h.
int MCElectronAnalyzer::nEvt_ [private] |
Definition at line 55 of file MCElectronAnalyzer.h.
Referenced by analyze(), beginJob(), and endJob().
int MCElectronAnalyzer::nMatched_ [private] |
Definition at line 56 of file MCElectronAnalyzer.h.
TProfile* MCElectronAnalyzer::p_BremVsEta_ [private] |
Definition at line 75 of file MCElectronAnalyzer.h.
Referenced by analyze(), and beginJob().
TProfile* MCElectronAnalyzer::p_BremVsR_ [private] |
Definition at line 74 of file MCElectronAnalyzer.h.
Referenced by analyze(), and beginJob().
std::string MCElectronAnalyzer::SimHitLabel [private] |
Definition at line 65 of file MCElectronAnalyzer.h.
std::string MCElectronAnalyzer::SimTkLabel [private] |
Definition at line 63 of file MCElectronAnalyzer.h.
std::string MCElectronAnalyzer::SimVtxLabel [private] |
Definition at line 64 of file MCElectronAnalyzer.h.
Definition at line 45 of file MCElectronAnalyzer.h.
Referenced by analyze(), beginJob(), and ~MCElectronAnalyzer().
const TrackerGeometry* MCElectronAnalyzer::trackerGeom [private] |
Definition at line 47 of file MCElectronAnalyzer.h.