CMS 3D CMS Logo

JetMuonHitsIDHelper.cc
Go to the documentation of this file.
2 
4 
20 
21 #include "TMath.h"
22 #include <vector>
23 #include <iostream>
24 
25 using namespace std;
26 
28  : trackingGeometryToken_(iC.esConsumes()) {
29  isRECO_ = true; // This will be "true" initially, then if the product isn't found, set to false once
34  numberOfHitsRPC_ = 0;
35  rpcRecHits_ = pset.getParameter<edm::InputTag>("rpcRecHits");
36 
38 }
39 
41  const edm::EventSetup& iSetup,
42  const reco::Jet& jet,
43  const int iDbg) {
44  // initialize
45  numberOfHits1RPC_ = 0;
46  numberOfHits2RPC_ = 0;
47  numberOfHits3RPC_ = 0;
48  numberOfHits4RPC_ = 0;
49  numberOfHitsRPC_ = 0;
50 
51  if (isRECO_) { // This will be "true" initially, then if the product isn't found, set to false once
52 
53  // Get tracking geometry
54  auto const& trackingGeometry = iSetup.getData(trackingGeometryToken_);
55 
56  //####READ RPC RecHits Collection########
57  //#In config: RpcRecHits = cms.InputTag("rpcRecHits")
58  edm::Handle<RPCRecHitCollection> rpcRecHits_handle;
59  event.getByToken(input_rpchits_token_, rpcRecHits_handle);
60 
61  if (!rpcRecHits_handle.isValid()) {
62  // don't throw exception if not running on RECO
63  edm::LogWarning("DataNotAvailable") << "JetMuonHitsIDHelper will not be run at all, this is not a RECO file.";
64  isRECO_ = false;
65  return;
66  }
67 
68  //####calculate rpc variables for each jet########
69 
70  for (RPCRecHitCollection::const_iterator itRPC = rpcRecHits_handle->begin(), itRPCEnd = rpcRecHits_handle->end();
71  itRPC != itRPCEnd;
72  ++itRPC) {
73  RPCRecHit const& hit = *itRPC;
74  DetId detid = hit.geographicalId();
75  LocalPoint lp = hit.localPosition();
76  const GeomDet* gd = trackingGeometry.idToDet(detid);
77  GlobalPoint gp = gd->toGlobal(lp);
78  double dR2 = reco::deltaR(jet.eta(), jet.phi(), static_cast<double>(gp.eta()), static_cast<double>(gp.phi()));
79  if (dR2 < 0.5) {
80  RPCDetId rpcChamberId = (RPCDetId)detid;
81  numberOfHitsRPC_++;
82  if (rpcChamberId.station() == 1)
83  numberOfHits1RPC_++;
84  if (rpcChamberId.station() == 2)
85  numberOfHits2RPC_++;
86  if (rpcChamberId.station() == 3)
87  numberOfHits3RPC_++;
88  if (rpcChamberId.station() == 4)
89  numberOfHits4RPC_++;
90  }
91  }
92  }
93 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
Base class for all types of Jets.
Definition: Jet.h:20
void calculate(const edm::Event &event, const edm::EventSetup &isetup, const reco::Jet &jet, const int iDbg=0)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::EDGetTokenT< RPCRecHitCollection > input_rpchits_token_
Definition: DetId.h:17
int station() const
Definition: RPCDetId.h:78
bool isValid() const
Definition: HandleBase.h:70
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:18
Log< level::Warning, false > LogWarning
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::LocalCoordinateSystemTag > LocalPoint
point in local coordinate system
Definition: Point3D.h:15
Definition: event.py:1