CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
RPCIntegrator.h
Go to the documentation of this file.
1 #ifndef Phase2L1Trigger_DTTrigger_RPCIntegrator_h
2 #define Phase2L1Trigger_DTTrigger_RPCIntegrator_h
3 
10 
20 
24 
25 //DT geometry
31 
33 
38  int rpcFlag;
39  int rpc_bx;
40  double rpc_t0;
41  RPCMetaprimitive(RPCDetId rpc_id_construct,
42  const RPCRecHit* rpc_cluster_construct,
43  GlobalPoint global_position_construct,
44  int rpcFlag_construct,
45  int rpc_bx_construct,
46  double rpc_t0_construct)
47  : rpc_id(rpc_id_construct),
48  rpc_cluster(rpc_cluster_construct),
49  global_position(global_position_construct),
50  rpcFlag(rpcFlag_construct),
51  rpc_bx(rpc_bx_construct),
52  rpc_t0(rpc_t0_construct) {}
53 };
54 
56 public:
59 
60  void initialise(const edm::EventSetup& iEventSetup, double shift_back_fromDT);
61  void finish();
62 
64  void matchWithDTAndUseRPCTime(std::vector<cmsdt::metaPrimitive>& dt_metaprimitives);
65  void makeRPCOnlySegments();
66  void storeRPCSingleHits();
67  void removeRPCHitsUsed();
68 
71  RPCDetId rpcDetId, int rpc_bx, double rpc_time, double rpc_global_phi, double phiB, int rpc_flag);
72 
73  double phiBending(RPCMetaprimitive* rpc_hit_1, RPCMetaprimitive* rpc_hit_2);
74  int phiInDTTPFormat(double rpc_global_phi, int rpcSector);
75  GlobalPoint RPCGlobalPosition(RPCDetId rpcId, const RPCRecHit& rpcIt) const;
76  double phi_DT_MP_conv(double rpc_global_phi, int rpcSector);
77  bool hasPosRF_rpc(int wh, int sec) const;
78 
79  std::vector<L1Phase2MuDTPhDigi> rpcRecHits_translated_;
80  std::vector<RPCMetaprimitive> RPCMetaprimitives_;
81 
82 private:
83  const bool m_debug_;
86  double m_phi_window_;
90 
93 
94  // Constant geometry values
95  //R[stat][layer] - radius of rpc station/layer from center of CMS
96  static constexpr double R_[2][2] = {{410.0, 444.8}, {492.7, 527.3}};
97  static constexpr double distance_between_two_rpc_layers_ = 35; // in cm
98 
99  double shift_back_;
100 };
101 #endif
double m_phi_window_
Definition: RPCIntegrator.h:86
void storeRPCSingleHits()
int m_max_quality_to_overwrite_t0_
Definition: RPCIntegrator.h:84
void matchWithDTAndUseRPCTime(std::vector< cmsdt::metaPrimitive > &dt_metaprimitives)
L1Phase2MuDTPhDigi createL1Phase2MuDTPhDigi(RPCDetId rpcDetId, int rpc_bx, double rpc_time, double rpc_global_phi, double phiB, int rpc_flag)
void removeRPCHitsUsed()
double phi_DT_MP_conv(double rpc_global_phi, int rpcSector)
std::vector< RPCMetaprimitive > RPCMetaprimitives_
Definition: RPCIntegrator.h:80
GlobalPoint global_position
Definition: RPCIntegrator.h:37
RPCGeometry const * rpcGeo_
Definition: RPCIntegrator.h:92
GlobalPoint RPCGlobalPosition(RPCDetId rpcId, const RPCRecHit &rpcIt) const
double shift_back_
Definition: RPCIntegrator.h:99
void makeRPCOnlySegments()
RPCIntegrator(const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
RPCMetaprimitive(RPCDetId rpc_id_construct, const RPCRecHit *rpc_cluster_construct, GlobalPoint global_position_construct, int rpcFlag_construct, int rpc_bx_construct, double rpc_t0_construct)
Definition: RPCIntegrator.h:41
static constexpr double distance_between_two_rpc_layers_
Definition: RPCIntegrator.h:97
std::vector< L1Phase2MuDTPhDigi > rpcRecHits_translated_
Definition: RPCIntegrator.h:79
const bool m_debug_
Definition: RPCIntegrator.h:83
bool hasPosRF_rpc(int wh, int sec) const
const RPCRecHit * rpc_cluster
Definition: RPCIntegrator.h:36
void prepareMetaPrimitives(edm::Handle< RPCRecHitCollection > rpcRecHits)
double phiBending(RPCMetaprimitive *rpc_hit_1, RPCMetaprimitive *rpc_hit_2)
bool m_storeAllRPCHits_
Definition: RPCIntegrator.h:87
int phiInDTTPFormat(double rpc_global_phi, int rpcSector)
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomH_
Definition: RPCIntegrator.h:88
void initialise(const edm::EventSetup &iEventSetup, double shift_back_fromDT)
static constexpr double R_[2][2]
Definition: RPCIntegrator.h:96
RPCMetaprimitive * matchDTwithRPC(cmsdt::metaPrimitive *dt_metaprimitive)
edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeomH_
Definition: RPCIntegrator.h:89
DTGeometry const * dtGeo_
Definition: RPCIntegrator.h:91