CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
MCElectronAnalyzer Class Reference

#include <MCElectronAnalyzer.h>

Inheritance diagram for MCElectronAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 
 MCElectronAnalyzer (const edm::ParameterSet &)
 
 ~MCElectronAnalyzer () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

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 More...
 
int nEvt_
 
int nMatched_
 
TProfile * p_BremVsEta_
 
TProfile * p_BremVsR_
 
std::string SimHitLabel
 
std::string SimTkLabel
 
std::string SimVtxLabel
 
ElectronMCTruthFindertheElectronMCTruthFinder_
 
const TrackerGeometrytrackerGeom
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 22 of file MCElectronAnalyzer.h.

Constructor & Destructor Documentation

MCElectronAnalyzer::MCElectronAnalyzer ( const edm::ParameterSet pset)
explicit

Definition at line 44 of file MCElectronAnalyzer.cc.

45  : fOutputFileName_( pset.getUntrackedParameter<string>("HistOutFile",std::string("TestConversions.root")) ),
46  fOutputFile_(nullptr)
47 {
48 
49 
50 }
T getUntrackedParameter(std::string const &, T const &) const
std::string fOutputFileName_
MCElectronAnalyzer::~MCElectronAnalyzer ( )
override

Definition at line 54 of file MCElectronAnalyzer.cc.

References theElectronMCTruthFinder_.

54  {
55 
56 
58 
59 }
ElectronMCTruthFinder * theElectronMCTruthFinder_

Member Function Documentation

void MCElectronAnalyzer::analyze ( const edm::Event e,
const edm::EventSetup  
)
override

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_.

145 {
146 
147 
148  using namespace edm;
149  //const float etaPhiDistance=0.01;
150  // Fiducial region
151  //const float TRK_BARL =0.9;
152  //const float BARL = 1.4442; // DAQ TDR p.290
153  //const float END_LO = 1.566;
154  //const float END_HI = 2.5;
155  // Electron mass
156  //const Float_t mElec= 0.000511;
157 
158 
159  nEvt_++;
160  LogInfo("MCElectronAnalyzer") << "MCElectronAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
161  // LogDebug("MCElectronAnalyzer") << "MCElectronAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
162  std::cout << "MCElectronAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
163 
165  std::cout << " MCElectronAnalyzer Looking for MC truth " << "\n";
166 
167  //get simtrack info
168  std::vector<SimTrack> theSimTracks;
169  std::vector<SimVertex> theSimVertices;
170 
173  e.getByLabel("g4SimHits",SimTk);
174  e.getByLabel("g4SimHits",SimVtx);
175 
176  theSimTracks.insert(theSimTracks.end(),SimTk->begin(),SimTk->end());
177  theSimVertices.insert(theSimVertices.end(),SimVtx->begin(),SimVtx->end());
178  std::cout << " MCElectronAnalyzer This Event has " << theSimTracks.size() << " sim tracks " << std::endl;
179  std::cout << " MCElectronAnalyzer This Event has " << theSimVertices.size() << " sim vertices " << std::endl;
180  if ( theSimTracks.empty() ) std::cout << " Event number " << e.id() << " has NO sim tracks " << std::endl;
181 
182 
183  std::vector<ElectronMCTruth> MCElectronctrons=theElectronMCTruthFinder_->find (theSimTracks, theSimVertices);
184  std::cout << " MCElectronAnalyzer MCElectronctrons size " << MCElectronctrons.size() << std::endl;
185 
186  for ( std::vector<ElectronMCTruth>::const_iterator iEl=MCElectronctrons.begin(); iEl !=MCElectronctrons.end(); ++iEl ){
187 
188  h_MCEleE_->Fill ( (*iEl).fourMomentum().e() );
189  h_MCEleEta_->Fill ( (*iEl).fourMomentum().pseudoRapidity() );
190  h_MCElePhi_->Fill ( (*iEl).fourMomentum().phi() );
191 
192  float totBrem=0 ;
193  unsigned int iBrem ;
194  for ( iBrem=0; iBrem < (*iEl).bremVertices().size(); ++iBrem ) {
195  float rBrem= (*iEl).bremVertices()[iBrem].perp();
196  float etaBrem= (*iEl).bremVertices()[iBrem].eta();
197  if ( rBrem < 120 ) {
198  totBrem += (*iEl).bremMomentum()[iBrem].e();
199  p_BremVsR_ ->Fill ( rBrem, (*iEl).bremMomentum()[iBrem].e() );
200  p_BremVsEta_ ->Fill ( etaBrem, (*iEl).bremMomentum()[iBrem].e() );
201 
202  }
203  }
204 
205 
206  h_BremFrac_->Fill( totBrem/(*iEl).fourMomentum().e() );
207  h_BremEnergy_->Fill ( totBrem );
208 
209 
210 
211 
212  }
213 
214 
215 
216 
217 }
ElectronMCTruthFinder * theElectronMCTruthFinder_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:475
edm::EventID id() const
Definition: EventBase.h:60
HLT enums.
std::vector< ElectronMCTruth > find(const std::vector< SimTrack > &simTracks, const std::vector< SimVertex > &simVertices)
void MCElectronAnalyzer::beginJob ( void  )
overridevirtual

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_, reco::return(), and theElectronMCTruthFinder_.

63 {
64 
65 
66  nEvt_=0;
67 
69 
70 
71  fOutputFile_ = new TFile( fOutputFileName_.c_str(), "RECREATE" ) ;
72 
74  h_MCEleE_ = new TH1F("MCEleE","MC ele energy",100,0.,200.);
75  h_MCElePhi_ = new TH1F("MCElePhi","MC ele phi",40,-3.14, 3.14);
76  h_MCEleEta_ = new TH1F("MCEleEta","MC ele eta",40,-3., 3.);
77  h_BremFrac_ = new TH1F("bremFrac","brem frac ", 100, 0., 1.);
78  h_BremEnergy_ = new TH1F("BremE","Brem energy",100,0.,200.);
79 
80  p_BremVsR_ = new TProfile("BremVsR", " Mean Brem energy vs R ", 48, 0., 120.);
81  p_BremVsEta_ = new TProfile("BremVsEta", " Mean Brem energy vs Eta ", 50, -2.5, 2.5);
82 
83 
84 
85 
86  return ;
87 }
ElectronMCTruthFinder * theElectronMCTruthFinder_
return(e1-e2)*(e1-e2)+dp *dp
std::string fOutputFileName_
void MCElectronAnalyzer::endJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 222 of file MCElectronAnalyzer.cc.

References gather_cfg::cout, fOutputFile_, nEvt_, and reco::return().

223 {
224 
225 
226 
227 
228  fOutputFile_->Write() ;
229  fOutputFile_->Close() ;
230 
231  edm::LogInfo("MCElectronAnalyzer") << "Analyzed " << nEvt_ << "\n";
232  std::cout << "MCElectronAnalyzer::endJob Analyzed " << nEvt_ << " events " << "\n";
233 
234  return ;
235 }
return(e1-e2)*(e1-e2)+dp *dp
float MCElectronAnalyzer::etaTransformation ( float  a,
float  b 
)
private

Definition at line 90 of file MCElectronAnalyzer.cc.

References funct::abs(), ETA, etaBarrelEndcap, cmsBatch::log, PI, R_ECAL, funct::tan(), and Z_Endcap.

90  {
91 
92 //---Definitions
93  const float PI = 3.1415927;
94  //const float TWOPI = 2.0*PI;
95 
96 //---Definitions for ECAL
97  const float R_ECAL = 136.5;
98  const float Z_Endcap = 328.0;
99  const float etaBarrelEndcap = 1.479;
100 
101 //---ETA correction
102 
103  float Theta = 0.0 ;
104  float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex;
105 
106  if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
107  if(Theta<0.0) Theta = Theta+PI ;
108  float ETA = - log(tan(0.5*Theta));
109 
110  if( std::abs(ETA) > etaBarrelEndcap )
111  {
112  float Zend = Z_Endcap ;
113  if(EtaParticle<0.0 ) Zend = -Zend ;
114  float Zlen = Zend - Zvertex ;
115  float RR = Zlen/sinh(EtaParticle);
116  Theta = atan(RR/Zend);
117  if(Theta<0.0) Theta = Theta+PI ;
118  ETA = - log(tan(0.5*Theta));
119  }
120 //---Return the result
121  return ETA;
122 //---end
123 }
#define ETA
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define PI
Definition: QcdUeDQM.h:36
static const float etaBarrelEndcap
static const float Z_Endcap
static const float R_ECAL
float MCElectronAnalyzer::phiNormalization ( float &  a)
private

Definition at line 125 of file MCElectronAnalyzer.cc.

References phi, PI, and TWOPI.

126 {
127 //---Definitions
128  const float PI = 3.1415927;
129  const float TWOPI = 2.0*PI;
130 
131 
132  if(phi > PI) {phi = phi - TWOPI;}
133  if(phi < -PI) {phi = phi + TWOPI;}
134 
135  // cout << " Float_t PHInormalization out " << PHI << endl;
136  return phi;
137 
138 }
#define TWOPI
Definition: DQMSourcePi0.cc:37
#define PI
Definition: QcdUeDQM.h:36

Member Data Documentation

TFile* MCElectronAnalyzer::fOutputFile_
private

Definition at line 49 of file MCElectronAnalyzer.h.

Referenced by beginJob(), and endJob().

std::string MCElectronAnalyzer::fOutputFileName_
private

Definition at line 48 of file MCElectronAnalyzer.h.

Referenced by beginJob().

TH1F* MCElectronAnalyzer::h_BremEnergy_
private

Definition at line 71 of file MCElectronAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_BremFrac_
private

Definition at line 70 of file MCElectronAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_MCEleE_
private

Definition at line 67 of file MCElectronAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_MCEleEta_
private

Definition at line 68 of file MCElectronAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_MCElePhi_
private

Definition at line 69 of file MCElectronAnalyzer.h.

Referenced by analyze(), and beginJob().

std::string MCElectronAnalyzer::HepMCLabel
private

Definition at line 61 of file MCElectronAnalyzer.h.

double MCElectronAnalyzer::mcEta_
private

Definition at line 59 of file MCElectronAnalyzer.h.

double MCElectronAnalyzer::mcPhi_
private

global variable for the MC photon

Definition at line 58 of file MCElectronAnalyzer.h.

int MCElectronAnalyzer::nEvt_
private

Definition at line 54 of file MCElectronAnalyzer.h.

Referenced by analyze(), beginJob(), and endJob().

int MCElectronAnalyzer::nMatched_
private

Definition at line 55 of file MCElectronAnalyzer.h.

TProfile* MCElectronAnalyzer::p_BremVsEta_
private

Definition at line 74 of file MCElectronAnalyzer.h.

Referenced by analyze(), and beginJob().

TProfile* MCElectronAnalyzer::p_BremVsR_
private

Definition at line 73 of file MCElectronAnalyzer.h.

Referenced by analyze(), and beginJob().

std::string MCElectronAnalyzer::SimHitLabel
private

Definition at line 64 of file MCElectronAnalyzer.h.

std::string MCElectronAnalyzer::SimTkLabel
private

Definition at line 62 of file MCElectronAnalyzer.h.

std::string MCElectronAnalyzer::SimVtxLabel
private

Definition at line 63 of file MCElectronAnalyzer.h.

ElectronMCTruthFinder* MCElectronAnalyzer::theElectronMCTruthFinder_
private

Definition at line 44 of file MCElectronAnalyzer.h.

Referenced by analyze(), beginJob(), and ~MCElectronAnalyzer().

const TrackerGeometry* MCElectronAnalyzer::trackerGeom
private

Definition at line 46 of file MCElectronAnalyzer.h.