CMS 3D CMS Logo

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