CMS 3D CMS Logo

EnergyScaleAnalyzer Class Reference

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

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

Inheritance diagram for EnergyScaleAnalyzer:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (edm::EventSetup const &)
virtual void endJob ()
 EnergyScaleAnalyzer (const edm::ParameterSet &)
 ~EnergyScaleAnalyzer ()

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_

Classes

struct  tree_structure_


Detailed Description

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

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

Definition at line 48 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_, and rootFile_.

00077 {
00078   
00079   hepMCLabel_ = ps.getParameter<std::string>("hepMCLabel");
00080 
00081   outputFile_   = ps.getParameter<std::string>("outputFile");
00082   rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE"); // open output file to store histograms
00083   
00084   evtN = 0;
00085 }

EnergyScaleAnalyzer::~EnergyScaleAnalyzer (  ) 

Definition at line 89 of file EnergyScaleAnalyzer.cc.

References rootFile_.

00091 {
00092   delete rootFile_;
00093 }


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 funct::abs(), kinem::delta_R(), dynamicHybridSuperClusters_cfi::dynamicHybridSuperClusters, fillTree(), MCTruth2::genEvent, edm::Event::getByLabel(), edm::Event::getManyByType(), hepMCLabel_, hybridSuperClusters_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, moduleLabel(), p, EnergyScaleAnalyzer::tree_structure_::parID, funct::sin(), std, tree_, xVtx_, yVtx_, and zVtx_.

00117                                                                            {
00118   using namespace edm; // needed for all fwk related classes
00119   using namespace std;
00120 
00121   //  std::cout << "Proceccing event # " << ++evtN << std::endl;
00122   
00123   //Get containers for MC truth, SC etc. ===================================================
00124   // =======================================================================================
00125   // =======================================================================================
00126   Handle<HepMCProduct> hepMC;
00127   evt.getByLabel( hepMCLabel_, hepMC ) ;
00128   
00129   const HepMC::GenEvent* genEvent = hepMC->GetEvent();
00130   if ( !(hepMC.isValid())) {
00131     LogInfo("EnergyScaleAnalyzer") << "Could not get MC Product!";
00132     return;
00133   }
00134   
00135 
00136   //=======================For Vertex correction
00137   std::vector< Handle< HepMCProduct > > evtHandles ;
00138   evt.getManyByType( evtHandles ) ;
00139   
00140   for ( unsigned int i=0; i<evtHandles.size(); ++i) {
00141     if ( evtHandles[i].isValid() ) {
00142       const HepMC::GenEvent* evt = evtHandles[i]->GetEvent() ;
00143       
00144       // take only 1st vertex for now - it's been tested only of PGuns...
00145       //
00146       HepMC::GenEvent::vertex_const_iterator vtx = evt->vertices_begin() ;
00147       if ( evtHandles[i].provenance()->moduleLabel() == hepMCLabel_ ) {
00148         //Corrdinates of Vertex w.r.o. the point (0,0,0)
00149         xVtx_ = 0.1*(*vtx)->position().x();      
00150         yVtx_ = 0.1*(*vtx)->position().y();
00151         zVtx_ = 0.1*(*vtx)->position().z();  
00152       }
00153     }
00154   }
00155   //==============================================================================
00156   //Get handle to SC collections
00157 
00158   Handle<reco::SuperClusterCollection> hybridSuperClusters;
00159   try {
00160     evt.getByLabel("hybridSuperClusters","",hybridSuperClusters);
00161   }catch (cms::Exception& ex) {
00162     edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer hybridSuperClusters.";
00163   }
00164 
00165   Handle<reco::SuperClusterCollection> dynamicHybridSuperClusters;
00166   try {
00167     evt.getByLabel("dynamicHybridSuperClusters","",dynamicHybridSuperClusters);
00168   }catch (cms::Exception& ex) {
00169     edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer dynamicHybridSuperClusters.";
00170   }
00171 
00172   Handle<reco::SuperClusterCollection> fixedMatrixSuperClustersWithPS;
00173   try {
00174     evt.getByLabel("fixedMatrixSuperClustersWithPreshower","",fixedMatrixSuperClustersWithPS);
00175   }catch (cms::Exception& ex) {
00176     edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer fixedMatrixSuperClustersWithPreshower.";
00177   }
00178 
00179   //Corrected collections
00180   Handle<reco::SuperClusterCollection> correctedHybridSC;
00181   try {
00182     evt.getByLabel("correctedHybridSuperClusters","",correctedHybridSC);
00183   }catch (cms::Exception& ex) {
00184     edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedHybridSuperClusters.";
00185   }
00186 
00187   Handle<reco::SuperClusterCollection> correctedDynamicHybridSC;
00188   try{
00189     evt.getByLabel("correctedDynamicHybridSuperClusters","",correctedDynamicHybridSC);
00190   }catch (cms::Exception& ex) {
00191     edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedDynamicHybridSuperClusters.";
00192   }
00193   
00194   Handle<reco::SuperClusterCollection> correctedFixedMatrixSCWithPS;
00195   try {
00196     evt.getByLabel("correctedFixedMatrixSuperClustersWithPreshower","",correctedFixedMatrixSCWithPS);
00197   }catch (cms::Exception& ex ) {
00198     edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedFixedMatrixSuperClustersWithPreshower.";
00199   }
00200 
00201   const reco::SuperClusterCollection*  dSC  = dynamicHybridSuperClusters.product();
00202   const reco::SuperClusterCollection*  hSC  = hybridSuperClusters.product();
00203   const reco::SuperClusterCollection* fmSC  = fixedMatrixSuperClustersWithPS.product();
00204   const reco::SuperClusterCollection* chSC  = correctedHybridSC.product();
00205   const reco::SuperClusterCollection* cdSC  = correctedDynamicHybridSC.product();
00206   const reco::SuperClusterCollection* cfmSC = correctedFixedMatrixSCWithPS.product();
00207  
00208   // -----------------------  Print outs for debugging
00209   /*
00210   std::cout << "MC truth" << std::endl;
00211   int counterI = 0;
00212   for(  HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin();
00213         p != genEvent->particles_end(); ++p ) {
00214     if ( fabs((*p)->momentum().eta()) < 1.5 ) {
00215       std::cout << ++counterI << " " << (*p)->momentum().e() << " " 
00216                 << (*p)->momentum().phi() << " " << (*p)->momentum().eta() << std::endl;
00217     }
00218   }
00219 
00220   std::cout << "Standard clusters" << std::endl;
00221   counterI = 0;
00222   for(reco::SuperClusterCollection::const_iterator em = hSC->begin();
00223       em != hSC->end(); ++em ) 
00224     std::cout << ++counterI << " " << em->energy() << " " << em->position().phi() << " " << em->position().eta() << std::endl;
00225 
00226   std::cout << "Dynamic clusters" << std::endl;
00227   counterI = 0;
00228   for(reco::SuperClusterCollection::const_iterator em = dSC->begin();
00229       em != dSC->end(); ++em ) 
00230     std::cout << ++counterI << " " << em->energy() << " " << em->position().phi() << " " << em->position().eta() << std::endl;
00231 
00232   std::cout << "FixedMatrix clusters with PS" << std::endl;
00233   counterI = 0;
00234   for(reco::SuperClusterCollection::const_iterator em = fmSC->begin();
00235       em != fmSC->end(); ++em )
00236     std::cout << ++counterI << " " << em->energy() << " " << em->position().phi() << " " << em->position().eta() << std::endl;
00237   */
00238   // -----------------------------
00239   //=====================================================================  
00240   // All containers are loaded, perform the analysis 
00241   //====================================================================
00242 
00243   // --------------------------- Store MC particles
00244   HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin();
00245   
00246   // Search for MC electrons or photons that satisfy the criteria
00247   float min_eT = 5;
00248   float max_eta = 2.5;
00249   
00250   std::vector<HepMC::GenParticle* > mcParticles;
00251   //int counter = 0;
00252   for ( HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin(); 
00253         p != genEvent->particles_end(); 
00254         ++p ) {
00255     //LogInfo("EnergyScaleAnalyzer") << "Particle " << ++counter 
00256     //<< " PDG ID = " << (*p)->pdg_id() << " pT = " << (*p)->momentum().perp();
00257     // require photon or electron
00258     if ( (*p)->pdg_id() != 22 && abs((*p)->pdg_id()) != 11 ) continue;
00259     
00260     // require selection criteria
00261     bool satisfySelectionCriteria = 
00262       (*p)->momentum().perp() > min_eT &&
00263       fabs((*p)->momentum().eta()) < max_eta;
00264     
00265     if ( ! satisfySelectionCriteria ) continue;
00266     
00267     // EM MC particle is found, save it in the vector
00268     mcParticles.push_back(*p);
00269   }
00270   // separation in dR between 2 first MC particles
00271   // should not be used for MC samples with > 2 em objects generated!
00272   if ( mcParticles.size() == 2 ) {
00273       HepMC::GenParticle* mc1 = mcParticles[0];
00274       HepMC::GenParticle* mc2 = mcParticles[1];
00275       tree_.mc_sep = kinem::delta_R(mc1->momentum().eta(), mc1->momentum().phi(),
00276                                     mc2->momentum().eta(), mc2->momentum().phi());
00277   } else
00278     tree_.mc_sep = -100;
00279 
00280   
00281   // now loop over MC particles, find the match with SC and do everything we need
00282   // then save info in the tree for every MC particle
00283   for(std::vector<HepMC::GenParticle* >::const_iterator p = mcParticles.begin();
00284       p != mcParticles.end(); ++p) {
00285     HepMC::GenParticle* mc = *p;
00286     
00287     // Fill MC information
00288     tree_.mc_npar  = mcParticles.size();
00289     tree_.parID    = mc->pdg_id();
00290     tree_.mc_e     = mc->momentum().e();
00291     tree_.mc_et    = mc->momentum().e()*sin(mc->momentum().theta());
00292     tree_.mc_phi   = mc->momentum().phi();
00293     tree_.mc_eta   = mc->momentum().eta();
00294     tree_.mc_theta = mc->momentum().theta();
00295 
00296 
00297     //Call function to fill tree
00298     // scType coprreponds:
00299     // HybridSuperCluster                     -- 1
00300     // DynamicHybridSuperCluster              -- 2
00301     // FixedMatrixSuperClustersWithPreshower  -- 3
00302 
00303     fillTree( hSC, chSC, mc, tree_, xVtx_, yVtx_, zVtx_, 1);
00304     //    std::cout << " TYPE " << 1 << " : " << tree_.em_e << " : " << tree_.em_phi << " : " << tree_.em_eta << std::endl;
00305 
00306     fillTree( dSC, cdSC, mc, tree_, xVtx_, yVtx_, zVtx_, 2);
00307     //    std::cout << " TYPE " << 2 << " : " << tree_.em_e << " : " << tree_.em_phi << " : " << tree_.em_eta << std::endl;
00308 
00309     fillTree( fmSC, cfmSC, mc, tree_, xVtx_, yVtx_, zVtx_, 3);
00310     //    std::cout << " TYPE " << 3 << " : " << tree_.em_e << " : " << tree_.em_phi << " : " << tree_.em_eta << std::endl;
00311 
00312     //   mytree_->Fill();
00313   } // loop over particles  
00314 }

void EnergyScaleAnalyzer::beginJob ( edm::EventSetup const &   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 97 of file EnergyScaleAnalyzer.cc.

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

00097                                                   {
00098 //========================================================================
00099 
00100   mytree_ = new TTree("energyScale","");
00101   TString treeVariables = "mc_npar/I:parID:mc_sep/F:mc_e:mc_et:mc_phi:mc_eta:mc_theta:";    // MC information
00102   treeVariables += "em_dR/F:";                 // MC <-> EM matching information
00103   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    
00104   treeVariables += "em_pet/F:em_pe:em_peta:em_ptheta:"                     ;  // EM SC physics (eta corrected information)
00105 
00106   treeVariables += "emCorr_e/F:emCorr_et:emCorr_eta:emCorr_phi:emCorr_theta:";// CMSSW standard corrections
00107   treeVariables += "emCorr_pet/F:emCorr_peta:emCorr_ptheta:"                 ;// CMSSW standard physics  
00108 
00109   treeVariables += "em_pw/F:em_ew:em_br"                                     ;  // EM widths pw -- phiWidth, ew -- etaWidth, ratios of pw/ew
00110 
00111   mytree_->Branch("energyScale",&(tree_.mc_npar),treeVariables);
00112       
00113 }

void EnergyScaleAnalyzer::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 443 of file EnergyScaleAnalyzer.cc.

References rootFile_.

00443                             {
00444   //========================================================================
00445   //Fill ROOT tree
00446   rootFile_->Write();
00447 }

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(), em, 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_, funct::log(), mytree_, phi, phiMax_, phiMaxVtx_, rClust_vtx_, funct::sin(), size, funct::sqrt(), funct::tan(), thetaMax_, thetaMaxVtx_, xClust_vtx_, xClust_zero_, yClust_vtx_, yClust_zero_, zClust_vtx_, and zClust_zero_.

Referenced by analyze().

00318                                                                                                                              {
00319   
00320   // -----------------------------  SuperClusters before energy correction
00321   reco::SuperClusterCollection::const_iterator em = scColl->end();
00322   float energyMax = -100.0; // dummy energy of the matched SC
00323   for(reco::SuperClusterCollection::const_iterator aClus = scColl->begin();
00324       aClus != scColl->end(); ++aClus) {
00325     // check the matching
00326     float dR = kinem::delta_R(mc   ->momentum().eta(), mc   ->momentum().phi(), 
00327                               aClus->position().eta(), aClus->position().phi());
00328     if (dR <  0.4) { // a rather loose matching cut
00329       float energy = aClus->energy();
00330       if ( energy < energyMax ) continue;
00331       energyMax = energy;
00332       em = aClus;
00333     }
00334   }
00335   
00336   if (em == scColl->end() ) {
00337     //    std::cout << "No matching SC with type " << scType << " was found for MC particle! "  << std::endl;
00338     //    std::cout << "Going to next type of SC. " << std::endl;
00339     return; 
00340   }
00341   //  ------------  
00342 
00343   tree_.em_scType = scType;
00344 
00345   tree_.em_isInCrack = 0;
00346   double emAbsEta = fabs(em->position().eta());
00347   // copied from RecoEgama/EgammaElectronAlgos/src/EgammaElectronClassification.cc
00348   if ( emAbsEta < 0.018 ||
00349        (emAbsEta > 0.423 && emAbsEta < 0.461) || 
00350        (emAbsEta > 0.770 && emAbsEta < 0.806) || 
00351        (emAbsEta > 1.127 && emAbsEta < 1.163) || 
00352        (emAbsEta > 1.460 && emAbsEta < 1.558) )
00353     tree_.em_isInCrack = 1;
00354   
00355   tree_.em_dR = kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(), 
00356                               em->position().eta(), em->position().phi());
00357   tree_.em_e     = em->energy();
00358   tree_.em_et    = em->energy() * sin(em->position().theta());
00359   tree_.em_phi   = em->position().phi();
00360   tree_.em_eta   = em->position().eta();
00361   tree_.em_theta = em->position().theta();
00362   tree_.em_nCell = (em->getHitsByDetId()).size();
00363   tree_.em_nBC   = em->clustersSize();
00364   
00365   //Get physics e, et etc:
00366   //Coordinates of EM object with respect of the point (0,0,0)
00367   xClust_zero_ = em->position().x();
00368   yClust_zero_ = em->position().y();
00369   zClust_zero_ = em->position().z();
00370   //Coordinates of EM object w.r.o. the Vertex position
00371   xClust_vtx_ = xClust_zero_ - xV;
00372   yClust_vtx_ = yClust_zero_ - yV;
00373   zClust_vtx_ = zClust_zero_ - zV;
00374     
00375   energyMax_ = em->energy();
00376   thetaMax_ = em->position().theta();
00377   etaMax_ = em->position().eta();
00378   phiMax_ = em->position().phi();
00379   eTMax_ = energyMax_ * sin(thetaMax_);
00380   if ( phiMax_ < 0) phiMax_ += 2*3.14159;
00381   
00382   rClust_vtx_ = sqrt (xClust_vtx_*xClust_vtx_ + yClust_vtx_*yClust_vtx_ + zClust_vtx_*zClust_vtx_);
00383   thetaMaxVtx_ = acos(zClust_vtx_/rClust_vtx_);
00384   etaMaxVtx_   = - log(tan(thetaMaxVtx_/2));
00385   eTMaxVtx_    = energyMax_ * sin(thetaMaxVtx_); 
00386   phiMaxVtx_   = atan2(yClust_vtx_,xClust_vtx_); 
00387   if ( phiMaxVtx_ < 0) phiMaxVtx_ += 2*3.14159;
00388   //=============================
00389   //parametres of EM object after vertex correction
00390   tree_.em_pet    = eTMaxVtx_;
00391   tree_.em_pe     = tree_.em_pet/sin(thetaMaxVtx_);
00392   tree_.em_peta   = etaMaxVtx_;
00393   tree_.em_ptheta = thetaMaxVtx_;
00394   
00395   
00396   //-------------------------------   Get SC after energy correction 
00397   em = corrSCColl->end();
00398   energyMax = -100.0; // dummy energy of the matched SC
00399   for(reco::SuperClusterCollection::const_iterator aClus = corrSCColl->begin();
00400       aClus != corrSCColl->end(); ++aClus) {
00401     // check the matching
00402     float dR = kinem::delta_R(mc   ->momentum().eta(), mc   ->momentum().phi(), 
00403                               aClus->position().eta(), aClus->position().phi());
00404     if (dR <  0.4) {
00405       float energy = aClus->energy();
00406       if ( energy < energyMax ) continue;
00407       energyMax = energy;
00408       em = aClus;
00409     }
00410   }
00411   
00412   if (em == corrSCColl->end() ) {
00413     //    std::cout << "No matching corrected SC with type " << scType << " was found for MC particle! "  << std::endl;
00414     //    std::cout << "Going to next type of SC. " << std::endl;
00415     return; 
00416   }
00417   //  ------------  
00418   
00420     tree_.emCorr_e     = em->energy();
00421     tree_.emCorr_et    = em->energy() * sin(em->position().theta());
00422     tree_.emCorr_phi   = em->position().phi();
00423     tree_.emCorr_eta   = em->position().eta();
00424     tree_.emCorr_theta = em->position().theta();
00425       
00426     // =========== Eta and Theta wrt Vertex does not change after energy corrections are applied
00427     // =========== So, no need to calculate them again
00428     
00429     tree_.emCorr_peta   = tree_.em_peta;
00430     tree_.emCorr_ptheta = tree_.em_ptheta;
00431     tree_.emCorr_pet    = tree_.emCorr_e/cosh(tree_.emCorr_peta);
00432 
00433     tree_.em_pw = em->phiWidth();
00434     tree_.em_ew = em->etaWidth();
00435     tree_.em_br = tree_.em_pw/tree_.em_ew;
00436         
00437     mytree_->Fill();
00438     
00439 }


Member Data Documentation

float EnergyScaleAnalyzer::energyMax_ [private]

Definition at line 134 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::etaMax_ [private]

Definition at line 137 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::etaMaxVtx_ [private]

Definition at line 138 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::eTMax_ [private]

Definition at line 135 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::eTMaxVtx_ [private]

Definition at line 136 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

int EnergyScaleAnalyzer::evtN [private]

Definition at line 145 of file EnergyScaleAnalyzer.h.

Referenced by EnergyScaleAnalyzer().

std::string EnergyScaleAnalyzer::hepMCLabel_ [private]

Definition at line 63 of file EnergyScaleAnalyzer.h.

Referenced by analyze(), and EnergyScaleAnalyzer().

TTree* EnergyScaleAnalyzer::mytree_ [private]

Definition at line 69 of file EnergyScaleAnalyzer.h.

Referenced by beginJob(), and fillTree().

std::string EnergyScaleAnalyzer::outputFile_ [private]

Definition at line 61 of file EnergyScaleAnalyzer.h.

Referenced by EnergyScaleAnalyzer().

float EnergyScaleAnalyzer::phiMax_ [private]

Definition at line 139 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::phiMaxVtx_ [private]

Definition at line 140 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::rClust_vtx_ [private]

Definition at line 132 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

TFile* EnergyScaleAnalyzer::rootFile_ [private]

Definition at line 66 of file EnergyScaleAnalyzer.h.

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

float EnergyScaleAnalyzer::thetaMax_ [private]

Definition at line 141 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::thetaMaxVtx_ [private]

Definition at line 142 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

tree_structure_ EnergyScaleAnalyzer::tree_ [private]

Definition at line 117 of file EnergyScaleAnalyzer.h.

Referenced by analyze(), and beginJob().

float EnergyScaleAnalyzer::xClust_vtx_ [private]

Definition at line 128 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::xClust_zero_ [private]

Definition at line 124 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::xVtx_ [private]

Definition at line 120 of file EnergyScaleAnalyzer.h.

Referenced by analyze().

float EnergyScaleAnalyzer::yClust_vtx_ [private]

Definition at line 129 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::yClust_zero_ [private]

Definition at line 125 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::yVtx_ [private]

Definition at line 121 of file EnergyScaleAnalyzer.h.

Referenced by analyze().

float EnergyScaleAnalyzer::zClust_vtx_ [private]

Definition at line 130 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::zClust_zero_ [private]

Definition at line 126 of file EnergyScaleAnalyzer.h.

Referenced by fillTree().

float EnergyScaleAnalyzer::zVtx_ [private]

Definition at line 122 of file EnergyScaleAnalyzer.h.

Referenced by analyze().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:20:23 2009 for CMSSW by  doxygen 1.5.4