CMS 3D CMS Logo

Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes

EBRecoSummary Class Reference

#include <EBRecoSummary.h>

Inheritance diagram for EBRecoSummary:
edm::EDAnalyzer

List of all members.

Public Member Functions

 EBRecoSummary (const edm::ParameterSet &)
 ~EBRecoSummary ()

Protected Attributes

edm::InputTag basicClusterCollection_EB_
double ethrEB_
edm::InputTag recHitCollection_EB_
edm::InputTag redRecHitCollection_EB_
double scEtThrEB_
edm::InputTag superClusterCollection_EB_

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
void convxtalid (int &, int &)
int diff_neta_s (int, int)
int diff_nphi_s (int, int)
virtual void endJob ()

Private Attributes

DQMStoredqmStore_
MonitorElementh_basicClusters_recHits_EB_recoFlag
MonitorElementh_recHits_EB_Chi2
MonitorElementh_recHits_EB_E1oE4
MonitorElementh_recHits_EB_energyMax
MonitorElementh_recHits_EB_recoFlag
MonitorElementh_recHits_EB_time
MonitorElementh_redRecHits_EB_recoFlag
MonitorElementh_superClusters_EB_E1oE4
MonitorElementh_superClusters_EB_nBC
MonitorElementh_superClusters_EB_phi
MonitorElementh_superClusters_eta
std::string prefixME_

Detailed Description

Definition at line 42 of file EBRecoSummary.h.


Constructor & Destructor Documentation

EBRecoSummary::EBRecoSummary ( const edm::ParameterSet ps) [explicit]

Definition at line 68 of file EBRecoSummary.cc.

References basicClusterCollection_EB_, dqmStore_, ethrEB_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), h_basicClusters_recHits_EB_recoFlag, h_recHits_EB_Chi2, h_recHits_EB_E1oE4, h_recHits_EB_energyMax, h_recHits_EB_recoFlag, h_recHits_EB_time, h_redRecHits_EB_recoFlag, h_superClusters_EB_E1oE4, h_superClusters_EB_nBC, h_superClusters_EB_phi, h_superClusters_eta, cppFunctionSkipper::operator, prefixME_, recHitCollection_EB_, redRecHitCollection_EB_, scEtThrEB_, and superClusterCollection_EB_.

{

  prefixME_ = ps.getUntrackedParameter<std::string>("prefixME", "");

  //now do what ever initialization is needed
  recHitCollection_EB_       = ps.getParameter<edm::InputTag>("recHitCollection_EB");
  redRecHitCollection_EB_    = ps.getParameter<edm::InputTag>("redRecHitCollection_EB");
  basicClusterCollection_EB_ = ps.getParameter<edm::InputTag>("basicClusterCollection_EB");
  superClusterCollection_EB_ = ps.getParameter<edm::InputTag>("superClusterCollection_EB");

  ethrEB_                    = ps.getParameter<double>("ethrEB");

  scEtThrEB_                 = ps.getParameter<double>("scEtThrEB");

  // DQM Store -------------------
  dqmStore_= edm::Service<DQMStore>().operator->();

  // Monitor Elements (ex THXD)
  dqmStore_->setCurrentFolder(prefixME_ + "/EBRecoSummary"); // to organise the histos in folders
     
  // ReducedRecHits ----------------------------------------------
  // ... barrel 
  h_redRecHits_EB_recoFlag = dqmStore_->book1D("redRecHits_EB_recoFlag","redRecHits_EB_recoFlag",16,-0.5,15.5);  

  // RecHits ---------------------------------------------- 
  // ... barrel
  h_recHits_EB_energyMax     = dqmStore_->book1D("recHits_EB_energyMax","recHitsEB_energyMax",110,-10,100);
  h_recHits_EB_Chi2          = dqmStore_->book1D("recHits_EB_Chi2","recHits_EB_Chi2",200,0,100);
  h_recHits_EB_time          = dqmStore_->book1D("recHits_EB_time","recHits_EB_time",200,-50,50);
  h_recHits_EB_E1oE4         = dqmStore_->book1D("recHits_EB_E1oE4","recHitsEB_E1oE4",150, 0, 1.5);
  h_recHits_EB_recoFlag      = dqmStore_->book1D("recHits_EB_recoFlag","recHits_EB_recoFlag",16,-0.5,15.5);  

  // Basic Clusters ----------------------------------------------    
  // ... associated barrel rec hits
  h_basicClusters_recHits_EB_recoFlag = dqmStore_->book1D("basicClusters_recHits_EB_recoFlag","basicClusters_recHits_EB_recoFlag",16,-0.5,15.5);  

  // Super Clusters ----------------------------------------------
  // ... barrel
  h_superClusters_EB_nBC     = dqmStore_->book1D("superClusters_EB_nBC","superClusters_EB_nBC",100,0.,100.);
  h_superClusters_EB_E1oE4   = dqmStore_->book1D("superClusters_EB_E1oE4","superClusters_EB_E1oE4",150,0,1.5);

  h_superClusters_eta        = dqmStore_->book1D("superClusters_eta","superClusters_eta",150,-3.,3.);
  h_superClusters_EB_phi     = dqmStore_->book1D("superClusters_EB_phi","superClusters_EB_phi",360,-3.1415927,3.1415927);
  
}
EBRecoSummary::~EBRecoSummary ( )

Definition at line 117 of file EBRecoSummary.cc.

{
        // do anything here that needs to be done at desctruction time
        // (e.g. close files, deallocate resources etc.)
}

Member Function Documentation

void EBRecoSummary::analyze ( const edm::Event ev,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 129 of file EBRecoSummary.cc.

References abs, basicClusterCollection_EB_, edm::SortedCollection< T, SORT >::begin(), EcalClusterTools::eBottom(), EcalClusterTools::eLeft(), jptDQMConfig_cff::eMax, edm::SortedCollection< T, SORT >::end(), relval_parameters_module::energy, EcalClusterTools::eRight(), eta(), ethrEB_, EcalClusterTools::eTop(), create_public_lumi_plots::exp, MonitorElement::Fill(), HcalObjRepresent::Fill(), edm::SortedCollection< T, SORT >::find(), edm::EventSetup::get(), edm::Event::getByLabel(), h_basicClusters_recHits_EB_recoFlag, h_recHits_EB_Chi2, h_recHits_EB_E1oE4, h_recHits_EB_energyMax, h_recHits_EB_recoFlag, h_recHits_EB_time, h_redRecHits_EB_recoFlag, h_superClusters_EB_E1oE4, h_superClusters_EB_nBC, h_superClusters_EB_phi, h_superClusters_eta, EBDetId::ieta(), edm::HandleBase::isValid(), phi, edm::ESHandle< T >::product(), edm::Handle< T >::product(), recHitCollection_EB_, redRecHitCollection_EB_, scEtThrEB_, funct::sin(), superClusterCollection_EB_, EcalTools::swissCross(), and cond::rpcobgas::time.

{
  
  //Get the magnetic field
  edm::ESHandle<MagneticField> theMagField;
  iSetup.get<IdealMagneticFieldRecord>().get(theMagField);

  // calo topology
  edm::ESHandle<CaloTopology> pTopology;
  iSetup.get<CaloTopologyRecord>().get(pTopology);
  const CaloTopology *topology = pTopology.product();

  // --- REDUCED REC HITS ------------------------------------------------------------------------------------- 
  edm::Handle<EcalRecHitCollection> redRecHitsEB;
  ev.getByLabel( redRecHitCollection_EB_, redRecHitsEB );
  const EcalRecHitCollection* theBarrelEcalredRecHits = redRecHitsEB.product () ;
  if ( ! redRecHitsEB.isValid() ) {
    edm::LogWarning("EBRecoSummary") << "redRecHitsEB not found"; 
  }
  
  for ( EcalRecHitCollection::const_iterator itr = theBarrelEcalredRecHits->begin () ;
        itr != theBarrelEcalredRecHits->end () ;++itr)
  {
      
    h_redRecHits_EB_recoFlag->Fill( itr -> recoFlag() );
  
  }

  // --- REC HITS ------------------------------------------------------------------------------------- 
  
  // ... barrel
  edm::Handle<EcalRecHitCollection> recHitsEB;
  ev.getByLabel( recHitCollection_EB_, recHitsEB );
  const EcalRecHitCollection* theBarrelEcalRecHits = recHitsEB.product () ;
  if ( ! recHitsEB.isValid() ) {
    edm::LogWarning("EBRecoSummary") << "recHitsEB not found"; 
  }

  float maxRecHitEnergyEB = -999.;
  
  EBDetId ebid_MrecHitEB;

  for ( EcalRecHitCollection::const_iterator itr = theBarrelEcalRecHits->begin () ;
        itr != theBarrelEcalRecHits->end () ;++itr)
    {

      EBDetId ebid( itr -> id() );
      
      h_recHits_EB_recoFlag      -> Fill( itr -> recoFlag() );

      // max E rec hit
      if (itr -> energy() > maxRecHitEnergyEB ){
        maxRecHitEnergyEB = itr -> energy() ;
      }       

      if ( itr -> energy() > ethrEB_ ){
        h_recHits_EB_Chi2          -> Fill( itr -> chi2() );
        h_recHits_EB_time          -> Fill( itr -> time() );
      }

      float R4 = EcalTools::swissCross( ebid, *theBarrelEcalRecHits, 0. );
      
      if ( itr -> energy() > 3. && abs(ebid.ieta())!=85 )  h_recHits_EB_E1oE4-> Fill( R4 );
      
    }
  
  h_recHits_EB_energyMax         -> Fill( maxRecHitEnergyEB );
  
  //--- BASIC CLUSTERS --------------------------------------------------------------

  // ... barrel
  edm::Handle<reco::BasicClusterCollection> basicClusters_EB_h;
  ev.getByLabel( basicClusterCollection_EB_, basicClusters_EB_h );
  const reco::BasicClusterCollection* theBarrelBasicClusters = basicClusters_EB_h.product () ;
  if ( ! basicClusters_EB_h.isValid() ) {
    edm::LogWarning("EBRecoSummary") << "basicClusters_EB_h not found"; 
  }

  for (reco::BasicClusterCollection::const_iterator itBC = theBarrelBasicClusters->begin(); 
       itBC != theBarrelBasicClusters->end(); ++itBC ) {
         
    //Get the associated RecHits
    const std::vector<std::pair<DetId,float> > & hits= itBC->hitsAndFractions();
    for (std::vector<std::pair<DetId,float> > ::const_iterator rh = hits.begin(); rh!=hits.end(); ++rh){
      
      EBRecHitCollection::const_iterator itrechit = theBarrelEcalRecHits->find((*rh).first);
      if (itrechit==theBarrelEcalRecHits->end()) continue;
      h_basicClusters_recHits_EB_recoFlag -> Fill ( itrechit -> recoFlag() );
    
    }
  
  }
 
  // Super Clusters
  // ... barrel
  edm::Handle<reco::SuperClusterCollection> superClusters_EB_h;
  ev.getByLabel( superClusterCollection_EB_, superClusters_EB_h );
  const reco::SuperClusterCollection* theBarrelSuperClusters = superClusters_EB_h.product () ;
  if ( ! superClusters_EB_h.isValid() ) {
    edm::LogWarning("EBRecoSummary") << "superClusters_EB_h not found"; 
  }

  for (reco::SuperClusterCollection::const_iterator itSC = theBarrelSuperClusters->begin(); 
       itSC != theBarrelSuperClusters->end(); ++itSC ) {
    
    double scEt = itSC -> energy() * sin(2.*atan( exp(- itSC->position().eta() )));
    
    if (scEt < scEtThrEB_ ) continue;

    h_superClusters_EB_nBC    -> Fill( itSC -> clustersSize());
    h_superClusters_eta       -> Fill( itSC -> eta() );
    h_superClusters_EB_phi    -> Fill( itSC -> phi() );
 
    float E1 = EcalClusterTools::eMax   ( *itSC, theBarrelEcalRecHits);
    float E4 = EcalClusterTools::eTop   ( *itSC, theBarrelEcalRecHits, topology )+
               EcalClusterTools::eRight ( *itSC, theBarrelEcalRecHits, topology )+
               EcalClusterTools::eBottom( *itSC, theBarrelEcalRecHits, topology )+
               EcalClusterTools::eLeft  ( *itSC, theBarrelEcalRecHits, topology );

    if ( E1 > 3. ) h_superClusters_EB_E1oE4  -> Fill( 1.- E4/E1);
    
  }

}
void EBRecoSummary::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 257 of file EBRecoSummary.cc.

{
}
void EBRecoSummary::convxtalid ( int &  ,
int &   
) [private]
int EBRecoSummary::diff_neta_s ( int  ,
int   
) [private]
int EBRecoSummary::diff_nphi_s ( int  ,
int   
) [private]
void EBRecoSummary::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 263 of file EBRecoSummary.cc.

{}

Member Data Documentation

Definition at line 95 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 61 of file EBRecoSummary.h.

Referenced by EBRecoSummary().

double EBRecoSummary::ethrEB_ [protected]

Definition at line 98 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 80 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 74 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 76 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 73 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 77 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 75 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 69 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 85 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 84 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 88 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 87 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

std::string EBRecoSummary::prefixME_ [private]

Definition at line 63 of file EBRecoSummary.h.

Referenced by EBRecoSummary().

Definition at line 93 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 94 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

double EBRecoSummary::scEtThrEB_ [protected]

Definition at line 100 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().

Definition at line 96 of file EBRecoSummary.h.

Referenced by analyze(), and EBRecoSummary().