CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
MCElectronAnalyzer Class Reference
Inheritance diagram for MCElectronAnalyzer:
edm::one::EDAnalyzer<> edm::one::EDAnalyzerBase 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::one::EDAnalyzer<>
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () 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
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex >
const & 
esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
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::array< std::vector< ModuleDescription const * > *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
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::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 30 of file MCElectronAnalyzer.cc.

Constructor & Destructor Documentation

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

Definition at line 79 of file MCElectronAnalyzer.cc.

80  : fOutputFileName_(pset.getUntrackedParameter<string>("HistOutFile", std::string("TestConversions.root"))),
81  fOutputFile_(nullptr) {}
T getUntrackedParameter(std::string const &, T const &) const
std::string fOutputFileName_
MCElectronAnalyzer::~MCElectronAnalyzer ( )
override

Definition at line 83 of file MCElectronAnalyzer.cc.

References theElectronMCTruthFinder_.

83 { delete theElectronMCTruthFinder_; }
ElectronMCTruthFinder * theElectronMCTruthFinder_

Member Function Documentation

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

Implements edm::one::EDAnalyzerBase.

Definition at line 158 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_.

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

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 85 of file MCElectronAnalyzer.cc.

References fOutputFile_, fOutputFileName_, h_BremEnergy_, h_BremFrac_, h_MCEleE_, h_MCEleEta_, h_MCElePhi_, nEvt_, p_BremVsEta_, p_BremVsR_, and theElectronMCTruthFinder_.

85  {
86  nEvt_ = 0;
87 
89 
90  fOutputFile_ = new TFile(fOutputFileName_.c_str(), "RECREATE");
91 
93  h_MCEleE_ = new TH1F("MCEleE", "MC ele energy", 100, 0., 200.);
94  h_MCElePhi_ = new TH1F("MCElePhi", "MC ele phi", 40, -3.14, 3.14);
95  h_MCEleEta_ = new TH1F("MCEleEta", "MC ele eta", 40, -3., 3.);
96  h_BremFrac_ = new TH1F("bremFrac", "brem frac ", 100, 0., 1.);
97  h_BremEnergy_ = new TH1F("BremE", "Brem energy", 100, 0., 200.);
98 
99  p_BremVsR_ = new TProfile("BremVsR", " Mean Brem energy vs R ", 48, 0., 120.);
100  p_BremVsEta_ = new TProfile("BremVsEta", " Mean Brem energy vs Eta ", 50, -2.5, 2.5);
101 
102  return;
103 }
ElectronMCTruthFinder * theElectronMCTruthFinder_
std::string fOutputFileName_
void MCElectronAnalyzer::endJob ( void  )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 221 of file MCElectronAnalyzer.cc.

References gather_cfg::cout, fOutputFile_, and nEvt_.

221  {
222  fOutputFile_->Write();
223  fOutputFile_->Close();
224 
225  edm::LogInfo("MCElectronAnalyzer") << "Analyzed " << nEvt_ << "\n";
226  std::cout << "MCElectronAnalyzer::endJob Analyzed " << nEvt_ << " events "
227  << "\n";
228 
229  return;
230 }
Log< level::Info, false > LogInfo
tuple cout
Definition: gather_cfg.py:144
float MCElectronAnalyzer::etaTransformation ( float  a,
float  b 
)
private

Definition at line 105 of file MCElectronAnalyzer.cc.

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

105  {
106  //---Definitions
107  const float PI = 3.1415927;
108  //const float TWOPI = 2.0*PI;
109 
110  //---Definitions for ECAL
111  const float R_ECAL = 136.5;
112  const float Z_Endcap = 328.0;
113  const float etaBarrelEndcap = 1.479;
114 
115  //---ETA correction
116 
117  float Theta = 0.0;
118  float ZEcal = R_ECAL * sinh(EtaParticle) + Zvertex;
119 
120  if (ZEcal != 0.0)
121  Theta = atan(R_ECAL / ZEcal);
122  if (Theta < 0.0)
123  Theta = Theta + PI;
124  float ETA = -log(tan(0.5 * Theta));
125 
126  if (std::abs(ETA) > etaBarrelEndcap) {
127  float Zend = Z_Endcap;
128  if (EtaParticle < 0.0)
129  Zend = -Zend;
130  float Zlen = Zend - Zvertex;
131  float RR = Zlen / sinh(EtaParticle);
132  Theta = atan(RR / Zend);
133  if (Theta < 0.0)
134  Theta = Theta + PI;
135  ETA = -log(tan(0.5 * Theta));
136  }
137  //---Return the result
138  return ETA;
139  //---end
140 }
static std::vector< std::string > checklist log
static constexpr float R_ECAL
#define ETA
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr float etaBarrelEndcap
#define PI
Definition: QcdUeDQM.h:37
static constexpr float ZEcal
static constexpr float Z_Endcap
float MCElectronAnalyzer::phiNormalization ( float &  a)
private

Definition at line 142 of file MCElectronAnalyzer.cc.

References phi, PI, and TWOPI.

142  {
143  //---Definitions
144  const float PI = 3.1415927;
145  const float TWOPI = 2.0 * PI;
146 
147  if (phi > PI) {
148  phi = phi - TWOPI;
149  }
150  if (phi < -PI) {
151  phi = phi + TWOPI;
152  }
153 
154  // cout << " Float_t PHInormalization out " << PHI << endl;
155  return phi;
156 }
#define PI
Definition: QcdUeDQM.h:37
const float TWOPI
2.0 * PI
Definition: IceTypes.h:23

Member Data Documentation

TFile* MCElectronAnalyzer::fOutputFile_
private

Definition at line 50 of file MCElectronAnalyzer.cc.

Referenced by beginJob(), and endJob().

std::string MCElectronAnalyzer::fOutputFileName_
private

Definition at line 49 of file MCElectronAnalyzer.cc.

Referenced by beginJob().

TH1F* MCElectronAnalyzer::h_BremEnergy_
private

Definition at line 68 of file MCElectronAnalyzer.cc.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_BremFrac_
private

Definition at line 67 of file MCElectronAnalyzer.cc.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_MCEleE_
private

Definition at line 64 of file MCElectronAnalyzer.cc.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_MCEleEta_
private

Definition at line 65 of file MCElectronAnalyzer.cc.

Referenced by analyze(), and beginJob().

TH1F* MCElectronAnalyzer::h_MCElePhi_
private

Definition at line 66 of file MCElectronAnalyzer.cc.

Referenced by analyze(), and beginJob().

std::string MCElectronAnalyzer::HepMCLabel
private

Definition at line 59 of file MCElectronAnalyzer.cc.

double MCElectronAnalyzer::mcEta_
private

Definition at line 57 of file MCElectronAnalyzer.cc.

double MCElectronAnalyzer::mcPhi_
private

global variable for the MC photon

Definition at line 56 of file MCElectronAnalyzer.cc.

int MCElectronAnalyzer::nEvt_
private

Definition at line 52 of file MCElectronAnalyzer.cc.

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

int MCElectronAnalyzer::nMatched_
private

Definition at line 53 of file MCElectronAnalyzer.cc.

TProfile* MCElectronAnalyzer::p_BremVsEta_
private

Definition at line 71 of file MCElectronAnalyzer.cc.

Referenced by analyze(), and beginJob().

TProfile* MCElectronAnalyzer::p_BremVsR_
private

Definition at line 70 of file MCElectronAnalyzer.cc.

Referenced by analyze(), and beginJob().

std::string MCElectronAnalyzer::SimHitLabel
private

Definition at line 62 of file MCElectronAnalyzer.cc.

std::string MCElectronAnalyzer::SimTkLabel
private

Definition at line 60 of file MCElectronAnalyzer.cc.

std::string MCElectronAnalyzer::SimVtxLabel
private

Definition at line 61 of file MCElectronAnalyzer.cc.

ElectronMCTruthFinder* MCElectronAnalyzer::theElectronMCTruthFinder_
private

Definition at line 45 of file MCElectronAnalyzer.cc.

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

const TrackerGeometry* MCElectronAnalyzer::trackerGeom
private

Definition at line 47 of file MCElectronAnalyzer.cc.