CMS 3D CMS Logo

Public Member Functions | Private Attributes

reco::helper::JetMuonHitsIDHelper Class Reference

#include <JetMuonHitsIDHelper.h>

List of all members.

Public Member Functions

void calculate (const edm::Event &event, const edm::EventSetup &isetup, const reco::Jet &jet, const int iDbg=0)
void fillDescription (edm::ParameterSetDescription &iDesc)
void initValues ()
 JetMuonHitsIDHelper (edm::ParameterSet const &pset)
 JetMuonHitsIDHelper ()
int numberOfHits1RPC () const
int numberOfHits2RPC () const
int numberOfHits3RPC () const
int numberOfHits4RPC () const
int numberOfHitsRPC () const
 ~JetMuonHitsIDHelper ()

Private Attributes

bool isRECO_
int numberOfHits1RPC_
int numberOfHits2RPC_
int numberOfHits3RPC_
int numberOfHits4RPC_
int numberOfHitsRPC_
edm::InputTag rpcRecHits_

Detailed Description

Definition at line 17 of file JetMuonHitsIDHelper.h.


Constructor & Destructor Documentation

reco::helper::JetMuonHitsIDHelper::JetMuonHitsIDHelper ( ) [inline]

Definition at line 21 of file JetMuonHitsIDHelper.h.

{}
reco::helper::JetMuonHitsIDHelper::JetMuonHitsIDHelper ( edm::ParameterSet const &  pset)

Definition at line 33 of file JetMuonHitsIDHelper.cc.

References edm::ParameterSet::getParameter().

{
  isRECO_ = true; // This will be "true" initially, then if the product isn't found, set to false once
  numberOfHits1RPC_ = 0;
  numberOfHits2RPC_ = 0;
  numberOfHits3RPC_ = 0;
  numberOfHits4RPC_ = 0;
  numberOfHitsRPC_ = 0;
  rpcRecHits_ = pset.getParameter<edm::InputTag>("rpcRecHits");
 
}
reco::helper::JetMuonHitsIDHelper::~JetMuonHitsIDHelper ( ) [inline]

Definition at line 23 of file JetMuonHitsIDHelper.h.

{} 

Member Function Documentation

void reco::helper::JetMuonHitsIDHelper::calculate ( const edm::Event event,
const edm::EventSetup isetup,
const reco::Jet jet,
const int  iDbg = 0 
)

Definition at line 48 of file JetMuonHitsIDHelper.cc.

References deltaR(), cond::rpcobgas::detid, reco::LeafCandidate::eta(), TrackingRecHit::geographicalId(), edm::EventSetup::get(), edm::HandleBase::isValid(), RPCRecHit::localPosition(), reco::LeafCandidate::phi(), RPCDetId, RPCDetId::station(), and GeomDet::toGlobal().

Referenced by JetIDProducer::produce().

{

  // initialize
  numberOfHits1RPC_ = 0;
  numberOfHits2RPC_ = 0;
  numberOfHits3RPC_ = 0;
  numberOfHits4RPC_ = 0;
  numberOfHitsRPC_ = 0;


  if ( isRECO_ ) { // This will be "true" initially, then if the product isn't found, set to false once

    // Get tracking geometry
    edm::ESHandle<GlobalTrackingGeometry> trackingGeometry;
    iSetup.get<GlobalTrackingGeometryRecord> ().get(trackingGeometry);
    
    //####READ RPC RecHits Collection########
    //#In config: RpcRecHits     = cms.InputTag("rpcRecHits")
    edm::Handle<RPCRecHitCollection> rpcRecHits_handle;
    event.getByLabel(rpcRecHits_, rpcRecHits_handle);


    if ( ! rpcRecHits_handle.isValid()  ) {
      // don't throw exception if not running on RECO
      edm::LogWarning("DataNotAvailable") << "JetMuonHitsIDHelper will not be run at all, this is not a RECO file.";
      isRECO_ = false;
      return;
    }

    //####calculate rpc  variables for each jet########
    
    for (  RPCRecHitCollection::const_iterator itRPC = rpcRecHits_handle->begin(),
             itRPCEnd = rpcRecHits_handle->end();
           itRPC != itRPCEnd; ++itRPC) {
      RPCRecHit const & hit = *itRPC;
      DetId detid = hit.geographicalId();
      LocalPoint lp = hit.localPosition();
      const GeomDet* gd = trackingGeometry->idToDet(detid);
      GlobalPoint gp = gd->toGlobal(lp);
      double dR2 = reco::deltaR(jet.eta(), jet.phi(), 
                                static_cast<double>( gp.eta() ), static_cast<double>(gp.phi()) );
      if (dR2 < 0.5) {
        RPCDetId rpcChamberId = (RPCDetId) detid;
        numberOfHitsRPC_++;
        if (rpcChamberId.station() == 1)
          numberOfHits1RPC_++;
        if (rpcChamberId.station() == 2)
          numberOfHits2RPC_++;
        if (rpcChamberId.station() == 3)
          numberOfHits3RPC_++;
        if (rpcChamberId.station() == 4)
          numberOfHits4RPC_++;
      }
    }
  }
}
void reco::helper::JetMuonHitsIDHelper::fillDescription ( edm::ParameterSetDescription iDesc)
void reco::helper::JetMuonHitsIDHelper::initValues ( )
int reco::helper::JetMuonHitsIDHelper::numberOfHits1RPC ( ) const [inline]

Definition at line 34 of file JetMuonHitsIDHelper.h.

References numberOfHits1RPC_.

{ return numberOfHits1RPC_;}
int reco::helper::JetMuonHitsIDHelper::numberOfHits2RPC ( ) const [inline]

Definition at line 35 of file JetMuonHitsIDHelper.h.

References numberOfHits2RPC_.

Referenced by JetIDProducer::produce().

{ return numberOfHits2RPC_;}
int reco::helper::JetMuonHitsIDHelper::numberOfHits3RPC ( ) const [inline]

Definition at line 36 of file JetMuonHitsIDHelper.h.

References numberOfHits3RPC_.

Referenced by JetIDProducer::produce().

{ return numberOfHits3RPC_;}
int reco::helper::JetMuonHitsIDHelper::numberOfHits4RPC ( ) const [inline]

Definition at line 37 of file JetMuonHitsIDHelper.h.

References numberOfHits4RPC_.

{ return numberOfHits4RPC_;}
int reco::helper::JetMuonHitsIDHelper::numberOfHitsRPC ( ) const [inline]

Definition at line 38 of file JetMuonHitsIDHelper.h.

References numberOfHitsRPC_.

Referenced by JetIDProducer::produce().

{ return numberOfHitsRPC_ ;}

Member Data Documentation

Definition at line 43 of file JetMuonHitsIDHelper.h.

Definition at line 45 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits1RPC().

Definition at line 46 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits2RPC().

Definition at line 47 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits3RPC().

Definition at line 48 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits4RPC().

Definition at line 49 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHitsRPC().

Definition at line 42 of file JetMuonHitsIDHelper.h.