CMS 3D CMS Logo

JetMuonHitsIDHelper.cc
Go to the documentation of this file.
2 
5 
7 // #include "TrackingTools/TrackAssociator/interface/MuonDetIdAssociator.h"
24 
25 #include "TMath.h"
26 #include <vector>
27 #include <iostream>
28 
29 using namespace std;
30 
31 
33 {
34  isRECO_ = true; // This will be "true" initially, then if the product isn't found, set to false once
35  numberOfHits1RPC_ = 0;
36  numberOfHits2RPC_ = 0;
37  numberOfHits3RPC_ = 0;
38  numberOfHits4RPC_ = 0;
39  numberOfHitsRPC_ = 0;
40  rpcRecHits_ = pset.getParameter<edm::InputTag>("rpcRecHits");
41 
42  input_rpchits_token_ = iC.consumes<RPCRecHitCollection>(rpcRecHits_);
43 
44 }
45 
46 
47 
48 
50  const reco::Jet &jet, const int iDbg )
51 {
52 
53  // initialize
54  numberOfHits1RPC_ = 0;
55  numberOfHits2RPC_ = 0;
56  numberOfHits3RPC_ = 0;
57  numberOfHits4RPC_ = 0;
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;
94  numberOfHitsRPC_++;
95  if (rpcChamberId.station() == 1)
96  numberOfHits1RPC_++;
97  if (rpcChamberId.station() == 2)
98  numberOfHits2RPC_++;
99  if (rpcChamberId.station() == 3)
100  numberOfHits3RPC_++;
101  if (rpcChamberId.station() == 4)
102  numberOfHits4RPC_++;
103  }
104  }
105  }
106 }
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:54
LocalPoint localPosition() const override
Return the 3-dimensional local position.
Definition: RPCRecHit.h:53
void calculate(const edm::Event &event, const edm::EventSetup &isetup, const reco::Jet &jet, const int iDbg=0)
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
bool isValid() const
Definition: HandleBase.h:74
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:59
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:17
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:96