CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
reco::helper::JetMuonHitsIDHelper Class Reference

#include <JetMuonHitsIDHelper.h>

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 ()
 
 JetMuonHitsIDHelper (edm::ParameterSet const &pset, edm::ConsumesCollector &&iC)
 
int numberOfHits1RPC () const
 
int numberOfHits2RPC () const
 
int numberOfHits3RPC () const
 
int numberOfHits4RPC () const
 
int numberOfHitsRPC () const
 
 ~JetMuonHitsIDHelper ()
 

Private Attributes

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

Detailed Description

Definition at line 19 of file JetMuonHitsIDHelper.h.

Constructor & Destructor Documentation

reco::helper::JetMuonHitsIDHelper::JetMuonHitsIDHelper ( )
inline

Definition at line 23 of file JetMuonHitsIDHelper.h.

References muonDTDigis_cfi::pset.

23 {}
reco::helper::JetMuonHitsIDHelper::JetMuonHitsIDHelper ( edm::ParameterSet const &  pset,
edm::ConsumesCollector &&  iC 
)

Definition at line 32 of file JetMuonHitsIDHelper.cc.

References edm::ParameterSet::getParameter().

33 {
34  isRECO_ = true; // This will be "true" initially, then if the product isn't found, set to false once
39  numberOfHitsRPC_ = 0;
40  rpcRecHits_ = pset.getParameter<edm::InputTag>("rpcRecHits");
41 
43 
44 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< RPCRecHitCollection > input_rpchits_token_
reco::helper::JetMuonHitsIDHelper::~JetMuonHitsIDHelper ( )
inline

Definition at line 25 of file JetMuonHitsIDHelper.h.

References calculate(), fillDescription(), initValues(), and metsig::jet.

25 {}

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 49 of file JetMuonHitsIDHelper.cc.

References reco::deltaR(), reco::LeafCandidate::eta(), TrackingRecHit::geographicalId(), edm::EventSetup::get(), runTauDisplay::gp, GlobalTrackingGeometry::idToDet(), edm::HandleBase::isValid(), RPCRecHit::localPosition(), reco::LeafCandidate::phi(), RPCDetId::station(), and GeomDet::toGlobal().

Referenced by JetIDProducer::produce(), and ~JetMuonHitsIDHelper().

51 {
52 
53  // initialize
58  numberOfHitsRPC_ = 0;
59 
60 
61  if ( isRECO_ ) { // This will be "true" initially, then if the product isn't found, set to false once
62 
63  // Get tracking geometry
65  iSetup.get<GlobalTrackingGeometryRecord> ().get(trackingGeometry);
66 
67  //####READ RPC RecHits Collection########
68  //#In config: RpcRecHits = cms.InputTag("rpcRecHits")
69  edm::Handle<RPCRecHitCollection> rpcRecHits_handle;
70  event.getByToken(input_rpchits_token_, rpcRecHits_handle);
71 
72 
73  if ( ! rpcRecHits_handle.isValid() ) {
74  // don't throw exception if not running on RECO
75  edm::LogWarning("DataNotAvailable") << "JetMuonHitsIDHelper will not be run at all, this is not a RECO file.";
76  isRECO_ = false;
77  return;
78  }
79 
80  //####calculate rpc variables for each jet########
81 
82  for ( RPCRecHitCollection::const_iterator itRPC = rpcRecHits_handle->begin(),
83  itRPCEnd = rpcRecHits_handle->end();
84  itRPC != itRPCEnd; ++itRPC) {
85  RPCRecHit const & hit = *itRPC;
86  DetId detid = hit.geographicalId();
87  LocalPoint lp = hit.localPosition();
88  const GeomDet* gd = trackingGeometry->idToDet(detid);
89  GlobalPoint gp = gd->toGlobal(lp);
90  double dR2 = reco::deltaR(jet.eta(), jet.phi(),
91  static_cast<double>( gp.eta() ), static_cast<double>(gp.phi()) );
92  if (dR2 < 0.5) {
93  RPCDetId rpcChamberId = (RPCDetId) detid;
95  if (rpcChamberId.station() == 1)
97  if (rpcChamberId.station() == 2)
99  if (rpcChamberId.station() == 3)
101  if (rpcChamberId.station() == 4)
103  }
104  }
105  }
106 }
double eta() const final
momentum pseudorapidity
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: RPCRecHit.h:53
bool isValid() const
Definition: HandleBase.h:74
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
edm::EDGetTokenT< RPCRecHitCollection > input_rpchits_token_
Definition: DetId.h:18
T eta() const
Definition: PV3DBase.h:76
const GeomDet * idToDet(DetId) const override
DetId geographicalId() const
double phi() const final
momentum azimuthal angle
int station() const
Definition: RPCDetId.h:96
void reco::helper::JetMuonHitsIDHelper::fillDescription ( edm::ParameterSetDescription iDesc)

Referenced by ~JetMuonHitsIDHelper().

void reco::helper::JetMuonHitsIDHelper::initValues ( )

Referenced by ~JetMuonHitsIDHelper().

int reco::helper::JetMuonHitsIDHelper::numberOfHits1RPC ( ) const
inline

Definition at line 36 of file JetMuonHitsIDHelper.h.

References numberOfHits1RPC_.

int reco::helper::JetMuonHitsIDHelper::numberOfHits2RPC ( ) const
inline

Definition at line 37 of file JetMuonHitsIDHelper.h.

References numberOfHits2RPC_.

Referenced by JetIDProducer::produce().

int reco::helper::JetMuonHitsIDHelper::numberOfHits3RPC ( ) const
inline

Definition at line 38 of file JetMuonHitsIDHelper.h.

References numberOfHits3RPC_.

Referenced by JetIDProducer::produce().

int reco::helper::JetMuonHitsIDHelper::numberOfHits4RPC ( ) const
inline

Definition at line 39 of file JetMuonHitsIDHelper.h.

References numberOfHits4RPC_.

int reco::helper::JetMuonHitsIDHelper::numberOfHitsRPC ( ) const
inline

Definition at line 40 of file JetMuonHitsIDHelper.h.

References numberOfHitsRPC_.

Referenced by JetIDProducer::produce().

Member Data Documentation

edm::EDGetTokenT<RPCRecHitCollection> reco::helper::JetMuonHitsIDHelper::input_rpchits_token_
private

Definition at line 53 of file JetMuonHitsIDHelper.h.

bool reco::helper::JetMuonHitsIDHelper::isRECO_
private

Definition at line 45 of file JetMuonHitsIDHelper.h.

int reco::helper::JetMuonHitsIDHelper::numberOfHits1RPC_
private

Definition at line 47 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits1RPC().

int reco::helper::JetMuonHitsIDHelper::numberOfHits2RPC_
private

Definition at line 48 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits2RPC().

int reco::helper::JetMuonHitsIDHelper::numberOfHits3RPC_
private

Definition at line 49 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits3RPC().

int reco::helper::JetMuonHitsIDHelper::numberOfHits4RPC_
private

Definition at line 50 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHits4RPC().

int reco::helper::JetMuonHitsIDHelper::numberOfHitsRPC_
private

Definition at line 51 of file JetMuonHitsIDHelper.h.

Referenced by numberOfHitsRPC().

edm::InputTag reco::helper::JetMuonHitsIDHelper::rpcRecHits_
private

Definition at line 44 of file JetMuonHitsIDHelper.h.