CMS 3D CMS Logo

JetMuonHitsIDHelper.cc
Go to the documentation of this file.
2 
5 
7 // #include "TrackingTools/TrackAssociator/interface/MuonDetIdAssociator.h"
25 
26 #include "TMath.h"
27 #include <vector>
28 #include <iostream>
29 
30 using namespace std;
31 
33  isRECO_ = true; // This will be "true" initially, then if the product isn't found, set to false once
34  numberOfHits1RPC_ = 0;
35  numberOfHits2RPC_ = 0;
36  numberOfHits3RPC_ = 0;
37  numberOfHits4RPC_ = 0;
38  numberOfHitsRPC_ = 0;
39  rpcRecHits_ = pset.getParameter<edm::InputTag>("rpcRecHits");
40 
41  input_rpchits_token_ = iC.consumes<RPCRecHitCollection>(rpcRecHits_);
42 }
43 
45  const edm::EventSetup& iSetup,
46  const reco::Jet& jet,
47  const int iDbg) {
48  // initialize
49  numberOfHits1RPC_ = 0;
50  numberOfHits2RPC_ = 0;
51  numberOfHits3RPC_ = 0;
52  numberOfHits4RPC_ = 0;
53  numberOfHitsRPC_ = 0;
54 
55  if (isRECO_) { // This will be "true" initially, then if the product isn't found, set to false once
56 
57  // Get tracking geometry
59  iSetup.get<GlobalTrackingGeometryRecord>().get(trackingGeometry);
60 
61  //####READ RPC RecHits Collection########
62  //#In config: RpcRecHits = cms.InputTag("rpcRecHits")
63  edm::Handle<RPCRecHitCollection> rpcRecHits_handle;
64  event.getByToken(input_rpchits_token_, rpcRecHits_handle);
65 
66  if (!rpcRecHits_handle.isValid()) {
67  // don't throw exception if not running on RECO
68  edm::LogWarning("DataNotAvailable") << "JetMuonHitsIDHelper will not be run at all, this is not a RECO file.";
69  isRECO_ = false;
70  return;
71  }
72 
73  //####calculate rpc variables for each jet########
74 
75  for (RPCRecHitCollection::const_iterator itRPC = rpcRecHits_handle->begin(), itRPCEnd = rpcRecHits_handle->end();
76  itRPC != itRPCEnd;
77  ++itRPC) {
78  RPCRecHit const& hit = *itRPC;
79  DetId detid = hit.geographicalId();
80  LocalPoint lp = hit.localPosition();
81  const GeomDet* gd = trackingGeometry->idToDet(detid);
82  GlobalPoint gp = gd->toGlobal(lp);
83  double dR2 = reco::deltaR(jet.eta(), jet.phi(), static_cast<double>(gp.eta()), static_cast<double>(gp.phi()));
84  if (dR2 < 0.5) {
85  RPCDetId rpcChamberId = (RPCDetId)detid;
86  numberOfHitsRPC_++;
87  if (rpcChamberId.station() == 1)
88  numberOfHits1RPC_++;
89  if (rpcChamberId.station() == 2)
90  numberOfHits2RPC_++;
91  if (rpcChamberId.station() == 3)
92  numberOfHits3RPC_++;
93  if (rpcChamberId.station() == 4)
94  numberOfHits4RPC_++;
95  }
96  }
97  }
98 }
T getParameter(std::string const &) const
double eta() const final
momentum pseudorapidity
Base class for all types of Jets.
Definition: Jet.h:20
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: RPCRecHit.h:37
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
bool isValid() const
Definition: HandleBase.h:70
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Definition: DetId.h:17
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:18
T get() const
Definition: EventSetup.h:73
const GeomDet * idToDet(DetId) const override
DetId geographicalId() const
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::LocalCoordinateSystemTag > LocalPoint
point in local coordinate system
Definition: Point3D.h:15
double phi() const final
momentum azimuthal angle
Definition: event.py:1
int station() const
Definition: RPCDetId.h:78