CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
EnergyScaleAnalyzer Class Reference

#include <Validation/EcalClusters/src/EnergyScaleAnalyzer.cc>

Inheritance diagram for EnergyScaleAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Classes

struct  tree_structure_
 

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
virtual void endJob ()
 
 EnergyScaleAnalyzer (const edm::ParameterSet &)
 
 ~EnergyScaleAnalyzer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void fillTree (const reco::SuperClusterCollection *scColl, const reco::SuperClusterCollection *corrSCColl, HepMC::GenParticle *mc, tree_structure_ &tree_, float xV, float yV, float zV, int scType)
 

Private Attributes

float energyMax_
 
float etaMax_
 
float etaMaxVtx_
 
float eTMax_
 
float eTMaxVtx_
 
int evtN
 
std::string hepMCLabel_
 
TTree * mytree_
 
std::string outputFile_
 
float phiMax_
 
float phiMaxVtx_
 
float rClust_vtx_
 
TFile * rootFile_
 
float thetaMax_
 
float thetaMaxVtx_
 
tree_structure_ tree_
 
float xClust_vtx_
 
float xClust_zero_
 
float xVtx_
 
float yClust_vtx_
 
float yClust_zero_
 
float yVtx_
 
float zClust_vtx_
 
float zClust_zero_
 
float zVtx_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
CurrentProcessingContext const * currentContext () const
 
- 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

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 49 of file EnergyScaleAnalyzer.h.

Constructor & Destructor Documentation

EnergyScaleAnalyzer::EnergyScaleAnalyzer ( const edm::ParameterSet ps)
explicit

Definition at line 75 of file EnergyScaleAnalyzer.cc.

References evtN, edm::ParameterSet::getParameter(), hepMCLabel_, outputFile_, rootFile_, and AlCaHLTBitMon_QueryRunRegistry::string.

77 {
78 
79  hepMCLabel_ = ps.getParameter<std::string>("hepMCLabel");
80 
81  outputFile_ = ps.getParameter<std::string>("outputFile");
82  rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE"); // open output file to store histograms
83 
84  evtN = 0;
85 }
T getParameter(std::string const &) const
EnergyScaleAnalyzer::~EnergyScaleAnalyzer ( )

Definition at line 89 of file EnergyScaleAnalyzer.cc.

References rootFile_.

91 {
92  delete rootFile_;
93 }

Member Function Documentation

void EnergyScaleAnalyzer::analyze ( const edm::Event evt,
const edm::EventSetup es 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 117 of file EnergyScaleAnalyzer.cc.

References abs, kinem::delta_R(), dynamicHybridSuperClusters_cfi::dynamicHybridSuperClusters, fillTree(), MCTruth::genEvent, configurableAnalysis::GenParticle, edm::Event::getByLabel(), edm::Event::getManyByType(), hepMCLabel_, unifiedSCCollection_cfi::hybridSuperClusters, i, EnergyScaleAnalyzer::tree_structure_::mc_e, EnergyScaleAnalyzer::tree_structure_::mc_et, EnergyScaleAnalyzer::tree_structure_::mc_eta, EnergyScaleAnalyzer::tree_structure_::mc_npar, EnergyScaleAnalyzer::tree_structure_::mc_phi, EnergyScaleAnalyzer::tree_structure_::mc_sep, EnergyScaleAnalyzer::tree_structure_::mc_theta, AlCaHLTBitMon_ParallelJobs::p, EnergyScaleAnalyzer::tree_structure_::parID, funct::sin(), tree_, xVtx_, yVtx_, and zVtx_.

117  {
118  using namespace edm; // needed for all fwk related classes
119  using namespace std;
120 
121  // std::cout << "Proceccing event # " << ++evtN << std::endl;
122 
123  //Get containers for MC truth, SC etc. ===================================================
124  // =======================================================================================
125  // =======================================================================================
126  Handle<HepMCProduct> hepMC;
127  evt.getByLabel( hepMCLabel_, hepMC ) ;
128 
129  const HepMC::GenEvent* genEvent = hepMC->GetEvent();
130  if ( !(hepMC.isValid())) {
131  LogInfo("EnergyScaleAnalyzer") << "Could not get MC Product!";
132  return;
133  }
134 
135 
136  //=======================For Vertex correction
137  std::vector< Handle< HepMCProduct > > evtHandles ;
138  evt.getManyByType( evtHandles ) ;
139 
140  for ( unsigned int i=0; i<evtHandles.size(); ++i) {
141  if ( evtHandles[i].isValid() ) {
142  const HepMC::GenEvent* evt = evtHandles[i]->GetEvent() ;
143 
144  // take only 1st vertex for now - it's been tested only of PGuns...
145  //
146  HepMC::GenEvent::vertex_const_iterator vtx = evt->vertices_begin() ;
147  if ( evtHandles[i].provenance()->moduleLabel() == hepMCLabel_ ) {
148  //Corrdinates of Vertex w.r.o. the point (0,0,0)
149  xVtx_ = 0.1*(*vtx)->position().x();
150  yVtx_ = 0.1*(*vtx)->position().y();
151  zVtx_ = 0.1*(*vtx)->position().z();
152  }
153  }
154  }
155  //==============================================================================
156  //Get handle to SC collections
157 
159  try {
160  evt.getByLabel("hybridSuperClusters","",hybridSuperClusters);
161  }catch (cms::Exception& ex) {
162  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer hybridSuperClusters.";
163  }
164 
166  try {
167  evt.getByLabel("dynamicHybridSuperClusters","",dynamicHybridSuperClusters);
168  }catch (cms::Exception& ex) {
169  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer dynamicHybridSuperClusters.";
170  }
171 
172  Handle<reco::SuperClusterCollection> fixedMatrixSuperClustersWithPS;
173  try {
174  evt.getByLabel("fixedMatrixSuperClustersWithPreshower","",fixedMatrixSuperClustersWithPS);
175  }catch (cms::Exception& ex) {
176  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer fixedMatrixSuperClustersWithPreshower.";
177  }
178 
179  //Corrected collections
180  Handle<reco::SuperClusterCollection> correctedHybridSC;
181  try {
182  evt.getByLabel("correctedHybridSuperClusters","",correctedHybridSC);
183  }catch (cms::Exception& ex) {
184  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedHybridSuperClusters.";
185  }
186 
187  Handle<reco::SuperClusterCollection> correctedDynamicHybridSC;
188  try{
189  evt.getByLabel("correctedDynamicHybridSuperClusters","",correctedDynamicHybridSC);
190  }catch (cms::Exception& ex) {
191  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedDynamicHybridSuperClusters.";
192  }
193 
194  Handle<reco::SuperClusterCollection> correctedFixedMatrixSCWithPS;
195  try {
196  evt.getByLabel("correctedFixedMatrixSuperClustersWithPreshower","",correctedFixedMatrixSCWithPS);
197  }catch (cms::Exception& ex ) {
198  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedFixedMatrixSuperClustersWithPreshower.";
199  }
200 
202  const reco::SuperClusterCollection* hSC = hybridSuperClusters.product();
203  const reco::SuperClusterCollection* fmSC = fixedMatrixSuperClustersWithPS.product();
204  const reco::SuperClusterCollection* chSC = correctedHybridSC.product();
205  const reco::SuperClusterCollection* cdSC = correctedDynamicHybridSC.product();
206  const reco::SuperClusterCollection* cfmSC = correctedFixedMatrixSCWithPS.product();
207 
208  // ----------------------- Print outs for debugging
209  /*
210  std::cout << "MC truth" << std::endl;
211  int counterI = 0;
212  for( HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin();
213  p != genEvent->particles_end(); ++p ) {
214  if ( fabs((*p)->momentum().eta()) < 1.5 ) {
215  std::cout << ++counterI << " " << (*p)->momentum().e() << " "
216  << (*p)->momentum().phi() << " " << (*p)->momentum().eta() << std::endl;
217  }
218  }
219 
220  std::cout << "Standard clusters" << std::endl;
221  counterI = 0;
222  for(reco::SuperClusterCollection::const_iterator em = hSC->begin();
223  em != hSC->end(); ++em )
224  std::cout << ++counterI << " " << em->energy() << " " << em->position().phi() << " " << em->position().eta() << std::endl;
225 
226  std::cout << "Dynamic clusters" << std::endl;
227  counterI = 0;
228  for(reco::SuperClusterCollection::const_iterator em = dSC->begin();
229  em != dSC->end(); ++em )
230  std::cout << ++counterI << " " << em->energy() << " " << em->position().phi() << " " << em->position().eta() << std::endl;
231 
232  std::cout << "FixedMatrix clusters with PS" << std::endl;
233  counterI = 0;
234  for(reco::SuperClusterCollection::const_iterator em = fmSC->begin();
235  em != fmSC->end(); ++em )
236  std::cout << ++counterI << " " << em->energy() << " " << em->position().phi() << " " << em->position().eta() << std::endl;
237  */
238  // -----------------------------
239  //=====================================================================
240  // All containers are loaded, perform the analysis
241  //====================================================================
242 
243  // --------------------------- Store MC particles
244  HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin();
245 
246  // Search for MC electrons or photons that satisfy the criteria
247  float min_eT = 5;
248  float max_eta = 2.5;
249 
250  std::vector<HepMC::GenParticle* > mcParticles;
251  //int counter = 0;
252  for ( HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin();
253  p != genEvent->particles_end();
254  ++p ) {
255  //LogInfo("EnergyScaleAnalyzer") << "Particle " << ++counter
256  //<< " PDG ID = " << (*p)->pdg_id() << " pT = " << (*p)->momentum().perp();
257  // require photon or electron
258  if ( (*p)->pdg_id() != 22 && abs((*p)->pdg_id()) != 11 ) continue;
259 
260  // require selection criteria
261  bool satisfySelectionCriteria =
262  (*p)->momentum().perp() > min_eT &&
263  fabs((*p)->momentum().eta()) < max_eta;
264 
265  if ( ! satisfySelectionCriteria ) continue;
266 
267  // EM MC particle is found, save it in the vector
268  mcParticles.push_back(*p);
269  }
270  // separation in dR between 2 first MC particles
271  // should not be used for MC samples with > 2 em objects generated!
272  if ( mcParticles.size() == 2 ) {
273  HepMC::GenParticle* mc1 = mcParticles[0];
274  HepMC::GenParticle* mc2 = mcParticles[1];
275  tree_.mc_sep = kinem::delta_R(mc1->momentum().eta(), mc1->momentum().phi(),
276  mc2->momentum().eta(), mc2->momentum().phi());
277  } else
278  tree_.mc_sep = -100;
279 
280 
281  // now loop over MC particles, find the match with SC and do everything we need
282  // then save info in the tree for every MC particle
283  for(std::vector<HepMC::GenParticle* >::const_iterator p = mcParticles.begin();
284  p != mcParticles.end(); ++p) {
285  HepMC::GenParticle* mc = *p;
286 
287  // Fill MC information
288  tree_.mc_npar = mcParticles.size();
289  tree_.parID = mc->pdg_id();
290  tree_.mc_e = mc->momentum().e();
291  tree_.mc_et = mc->momentum().e()*sin(mc->momentum().theta());
292  tree_.mc_phi = mc->momentum().phi();
293  tree_.mc_eta = mc->momentum().eta();
294  tree_.mc_theta = mc->momentum().theta();
295 
296 
297  //Call function to fill tree
298  // scType coprreponds:
299  // HybridSuperCluster -- 1
300  // DynamicHybridSuperCluster -- 2
301  // FixedMatrixSuperClustersWithPreshower -- 3
302 
303  fillTree( hSC, chSC, mc, tree_, xVtx_, yVtx_, zVtx_, 1);
304  // std::cout << " TYPE " << 1 << " : " << tree_.em_e << " : " << tree_.em_phi << " : " << tree_.em_eta << std::endl;
305 
306  fillTree( dSC, cdSC, mc, tree_, xVtx_, yVtx_, zVtx_, 2);
307  // std::cout << " TYPE " << 2 << " : " << tree_.em_e << " : " << tree_.em_phi << " : " << tree_.em_eta << std::endl;
308 
309  fillTree( fmSC, cfmSC, mc, tree_, xVtx_, yVtx_, zVtx_, 3);
310  // std::cout << " TYPE " << 3 << " : " << tree_.em_e << " : " << tree_.em_phi << " : " << tree_.em_eta << std::endl;
311 
312  // mytree_->Fill();
313  } // loop over particles
314 }
void getManyByType(std::vector< Handle< PROD > > &results) const
Definition: Event.h:395
int i
Definition: DBlmapReader.cc:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define abs(x)
Definition: mlp_lapack.h:159
void fillTree(const reco::SuperClusterCollection *scColl, const reco::SuperClusterCollection *corrSCColl, HepMC::GenParticle *mc, tree_structure_ &tree_, float xV, float yV, float zV, int scType)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
tuple genEvent
Definition: MCTruth.py:33
double delta_R(double eta1, double phi1, double eta2, double phi2)
Definition: AnglesUtil.h:116
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
void EnergyScaleAnalyzer::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 97 of file EnergyScaleAnalyzer.cc.

References EnergyScaleAnalyzer::tree_structure_::mc_npar, mytree_, and tree_.

97  {
98 //========================================================================
99 
100  mytree_ = new TTree("energyScale","");
101  TString treeVariables = "mc_npar/I:parID:mc_sep/F:mc_e:mc_et:mc_phi:mc_eta:mc_theta:"; // MC information
102  treeVariables += "em_dR/F:"; // MC <-> EM matching information
103  treeVariables += "em_isInCrack/I:em_scType:em_e/F:em_et:em_phi:em_eta:em_theta:em_nCell/I:em_nBC:"; // EM SC info
104  treeVariables += "em_pet/F:em_pe:em_peta:em_ptheta:" ; // EM SC physics (eta corrected information)
105 
106  treeVariables += "emCorr_e/F:emCorr_et:emCorr_eta:emCorr_phi:emCorr_theta:";// CMSSW standard corrections
107  treeVariables += "emCorr_pet/F:emCorr_peta:emCorr_ptheta:" ;// CMSSW standard physics
108 
109  treeVariables += "em_pw/F:em_ew:em_br" ; // EM widths pw -- phiWidth, ew -- etaWidth, ratios of pw/ew
110 
111  mytree_->Branch("energyScale",&(tree_.mc_npar),treeVariables);
112 
113 }
void EnergyScaleAnalyzer::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 443 of file EnergyScaleAnalyzer.cc.

References rootFile_.

443  {
444  //========================================================================
445  //Fill ROOT tree
446  rootFile_->Write();
447 }
void EnergyScaleAnalyzer::fillTree ( const reco::SuperClusterCollection scColl,
const reco::SuperClusterCollection corrSCColl,
HepMC::GenParticle *  mc,
tree_structure_ tree_,
float  xV,
float  yV,
float  zV,
int  scType 
)
private

fill tree with kinematic variables of corrected Super Cluster

Definition at line 317 of file EnergyScaleAnalyzer.cc.

References kinem::delta_R(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, EnergyScaleAnalyzer::tree_structure_::em_br, EnergyScaleAnalyzer::tree_structure_::em_dR, EnergyScaleAnalyzer::tree_structure_::em_e, EnergyScaleAnalyzer::tree_structure_::em_et, EnergyScaleAnalyzer::tree_structure_::em_eta, EnergyScaleAnalyzer::tree_structure_::em_ew, EnergyScaleAnalyzer::tree_structure_::em_isInCrack, EnergyScaleAnalyzer::tree_structure_::em_nBC, EnergyScaleAnalyzer::tree_structure_::em_nCell, EnergyScaleAnalyzer::tree_structure_::em_pe, EnergyScaleAnalyzer::tree_structure_::em_pet, EnergyScaleAnalyzer::tree_structure_::em_peta, EnergyScaleAnalyzer::tree_structure_::em_phi, EnergyScaleAnalyzer::tree_structure_::em_ptheta, EnergyScaleAnalyzer::tree_structure_::em_pw, EnergyScaleAnalyzer::tree_structure_::em_scType, EnergyScaleAnalyzer::tree_structure_::em_theta, EnergyScaleAnalyzer::tree_structure_::emCorr_e, EnergyScaleAnalyzer::tree_structure_::emCorr_et, EnergyScaleAnalyzer::tree_structure_::emCorr_eta, EnergyScaleAnalyzer::tree_structure_::emCorr_pet, EnergyScaleAnalyzer::tree_structure_::emCorr_peta, EnergyScaleAnalyzer::tree_structure_::emCorr_phi, EnergyScaleAnalyzer::tree_structure_::emCorr_ptheta, EnergyScaleAnalyzer::tree_structure_::emCorr_theta, relval_parameters_module::energy, energyMax_, eta(), etaMax_, etaMaxVtx_, eTMax_, eTMaxVtx_, create_public_lumi_plots::log, mytree_, phi, phiMax_, phiMaxVtx_, rClust_vtx_, funct::sin(), mathSSE::sqrt(), funct::tan(), thetaMax_, thetaMaxVtx_, xClust_vtx_, xClust_zero_, yClust_vtx_, yClust_zero_, zClust_vtx_, and zClust_zero_.

Referenced by analyze().

318  {
319 
320  // ----------------------------- SuperClusters before energy correction
321  reco::SuperClusterCollection::const_iterator em = scColl->end();
322  float energyMax = -100.0; // dummy energy of the matched SC
323  for(reco::SuperClusterCollection::const_iterator aClus = scColl->begin();
324  aClus != scColl->end(); ++aClus) {
325  // check the matching
326  float dR = kinem::delta_R(mc ->momentum().eta(), mc ->momentum().phi(),
327  aClus->position().eta(), aClus->position().phi());
328  if (dR < 0.4) { // a rather loose matching cut
329  float energy = aClus->energy();
330  if ( energy < energyMax ) continue;
331  energyMax = energy;
332  em = aClus;
333  }
334  }
335 
336  if (em == scColl->end() ) {
337  // std::cout << "No matching SC with type " << scType << " was found for MC particle! " << std::endl;
338  // std::cout << "Going to next type of SC. " << std::endl;
339  return;
340  }
341  // ------------
342 
343  tree_.em_scType = scType;
344 
345  tree_.em_isInCrack = 0;
346  double emAbsEta = fabs(em->position().eta());
347  // copied from RecoEgama/EgammaElectronAlgos/src/EgammaElectronClassification.cc
348  if ( emAbsEta < 0.018 ||
349  (emAbsEta > 0.423 && emAbsEta < 0.461) ||
350  (emAbsEta > 0.770 && emAbsEta < 0.806) ||
351  (emAbsEta > 1.127 && emAbsEta < 1.163) ||
352  (emAbsEta > 1.460 && emAbsEta < 1.558) )
353  tree_.em_isInCrack = 1;
354 
355  tree_.em_dR = kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(),
356  em->position().eta(), em->position().phi());
357  tree_.em_e = em->energy();
358  tree_.em_et = em->energy() * sin(em->position().theta());
359  tree_.em_phi = em->position().phi();
360  tree_.em_eta = em->position().eta();
361  tree_.em_theta = em->position().theta();
362  tree_.em_nCell = em->size();
363  tree_.em_nBC = em->clustersSize();
364 
365  //Get physics e, et etc:
366  //Coordinates of EM object with respect of the point (0,0,0)
367  xClust_zero_ = em->position().x();
368  yClust_zero_ = em->position().y();
369  zClust_zero_ = em->position().z();
370  //Coordinates of EM object w.r.o. the Vertex position
371  xClust_vtx_ = xClust_zero_ - xV;
372  yClust_vtx_ = yClust_zero_ - yV;
373  zClust_vtx_ = zClust_zero_ - zV;
374 
375  energyMax_ = em->energy();
376  thetaMax_ = em->position().theta();
377  etaMax_ = em->position().eta();
378  phiMax_ = em->position().phi();
380  if ( phiMax_ < 0) phiMax_ += 2*3.14159;
381 
384  etaMaxVtx_ = - log(tan(thetaMaxVtx_/2));
387  if ( phiMaxVtx_ < 0) phiMaxVtx_ += 2*3.14159;
388  //=============================
389  //parametres of EM object after vertex correction
394 
395 
396  //------------------------------- Get SC after energy correction
397  em = corrSCColl->end();
398  energyMax = -100.0; // dummy energy of the matched SC
399  for(reco::SuperClusterCollection::const_iterator aClus = corrSCColl->begin();
400  aClus != corrSCColl->end(); ++aClus) {
401  // check the matching
402  float dR = kinem::delta_R(mc ->momentum().eta(), mc ->momentum().phi(),
403  aClus->position().eta(), aClus->position().phi());
404  if (dR < 0.4) {
405  float energy = aClus->energy();
406  if ( energy < energyMax ) continue;
407  energyMax = energy;
408  em = aClus;
409  }
410  }
411 
412  if (em == corrSCColl->end() ) {
413  // std::cout << "No matching corrected SC with type " << scType << " was found for MC particle! " << std::endl;
414  // std::cout << "Going to next type of SC. " << std::endl;
415  return;
416  }
417  // ------------
418 
420  tree_.emCorr_e = em->energy();
421  tree_.emCorr_et = em->energy() * sin(em->position().theta());
422  tree_.emCorr_phi = em->position().phi();
423  tree_.emCorr_eta = em->position().eta();
424  tree_.emCorr_theta = em->position().theta();
425 
426  // =========== Eta and Theta wrt Vertex does not change after energy corrections are applied
427  // =========== So, no need to calculate them again
428 
432 
433  tree_.em_pw = em->phiWidth();
434  tree_.em_ew = em->etaWidth();
436 
437  mytree_->Fill();
438 
439 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T eta() const
T sqrt(T t)
Definition: SSEVec.h:48
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
double delta_R(double eta1, double phi1, double eta2, double phi2)
Definition: AnglesUtil.h:116
Definition: DDAxes.h:10

Member Data Documentation

float EnergyScaleAnalyzer::energyMax_
private

Definition at line 135 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::etaMax_
private

Definition at line 138 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::etaMaxVtx_
private

Definition at line 139 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::eTMax_
private

Definition at line 136 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::eTMaxVtx_
private

Definition at line 137 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

int EnergyScaleAnalyzer::evtN
private

Definition at line 146 of file EnergyScaleAnalyzer.h.

Referenced by EnergyScaleAnalyzer().

std::string EnergyScaleAnalyzer::hepMCLabel_
private

Definition at line 64 of file EnergyScaleAnalyzer.h.

Referenced by analyze(), and EnergyScaleAnalyzer().

TTree* EnergyScaleAnalyzer::mytree_
private

Definition at line 70 of file EnergyScaleAnalyzer.h.

Referenced by beginJob(), and fillTree().

std::string EnergyScaleAnalyzer::outputFile_
private

Definition at line 62 of file EnergyScaleAnalyzer.h.

Referenced by EnergyScaleAnalyzer().

float EnergyScaleAnalyzer::phiMax_
private

Definition at line 140 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::phiMaxVtx_
private

Definition at line 141 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::rClust_vtx_
private

Definition at line 133 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

TFile* EnergyScaleAnalyzer::rootFile_
private

Definition at line 67 of file EnergyScaleAnalyzer.h.

Referenced by endJob(), EnergyScaleAnalyzer(), and ~EnergyScaleAnalyzer().

float EnergyScaleAnalyzer::thetaMax_
private

Definition at line 142 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::thetaMaxVtx_
private

Definition at line 143 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

tree_structure_ EnergyScaleAnalyzer::tree_
private

Definition at line 118 of file EnergyScaleAnalyzer.h.

Referenced by analyze(), and beginJob().

float EnergyScaleAnalyzer::xClust_vtx_
private

Definition at line 129 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::xClust_zero_
private

Definition at line 125 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::xVtx_
private

Definition at line 121 of file EnergyScaleAnalyzer.h.

Referenced by analyze().

float EnergyScaleAnalyzer::yClust_vtx_
private

Definition at line 130 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::yClust_zero_
private

Definition at line 126 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::yVtx_
private

Definition at line 122 of file EnergyScaleAnalyzer.h.

Referenced by analyze().

float EnergyScaleAnalyzer::zClust_vtx_
private

Definition at line 131 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::zClust_zero_
private

Definition at line 127 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::zVtx_
private

Definition at line 123 of file EnergyScaleAnalyzer.h.

Referenced by analyze().