CMS 3D CMS Logo

RPCTPCollector.cc
Go to the documentation of this file.
3 
12 
13 using namespace emtf::phase2;
14 
16  : context_(context),
17  input_token_(i_consumes_collector.consumes<RPCTag::rechit_collection_type>(context.config_.rpc_input_)) {}
18 
19 void RPCTPCollector::collect(const edm::Event& i_event, BXTPCMap& bx_tpc_map) const {
20  // Constants
21  static const int clus_width_cut = 4;
22  static const int clus_width_cut_irpc = 6;
23 
24  // Read RPC digis
25  TPCollection tpc;
26 
28  i_event.getByToken(input_token_, rpc_digis);
29 
30  auto digi = rpc_digis->begin();
31  auto digi_end = rpc_digis->end();
32 
33  for (; digi != digi_end; ++digi) {
34  tpc.emplace_back(digi->rpcId(), *digi);
35  }
36 
37  // Map to BX
38  for (auto& tp_entry : tpc) {
39  const auto& tp_det_id = tp_entry.tp_.detId<RPCDetId>();
40  const RPCData& tp_data = tp_entry.tp_.getRPCData();
41 
42  const int tp_region = tp_det_id.region(); // 0: barrel, +/-1: endcap
43  const int tp_endcap = (tp_region == -1) ? 2 : tp_region; // 1: +endcap, 2: -endcap
44  const int tp_endcap_pm = (tp_endcap == 2) ? -1 : tp_endcap; // 1: +endcap, -1: -endcap
45 
46  // RPC sector is rotated by -20 deg relative to CSC sector.
47  // RPC sector 1 starts at -5 deg, CSC sector 1 starts at 15 deg.
48  const int tp_rpc_sector = tp_det_id.sector(); // 1 - 6 (60 degrees in phi, sector 1 begins at -5 deg)
49 
50  // RPC subsector is defined differently than CSC subsector.
51  // RPC subsector is used to label the chamber within a sector.
52  const int tp_rpc_subsector = tp_det_id.subsector();
53 
54  const int tp_station = tp_det_id.station(); // 1 - 4
55  const int tp_ring = tp_det_id.ring(); // 2 - 3 (increasing theta)
56  const int tp_roll =
57  tp_det_id.roll(); // 1 - 3 (decreasing theta; aka A - C; space between rolls is 9 - 15 in theta_fp)
58  const int tp_layer = tp_det_id.layer(); // Always 1 in the Endcap, 1 or 2 in the Barrel
59 
60  const int tp_strip = (tp_data.strip_low + tp_data.strip_hi) / 2; // in full-strip unit
61  const int tp_strip_lo = tp_data.strip_low;
62  const int tp_strip_hi = tp_data.strip_hi;
63  const int tp_clus_width = (tp_strip_hi - tp_strip_lo + 1);
64 
65  const bool tp_is_CPPF = tp_data.isCPPF;
66 
67  const int tp_bx = tp_data.bx + this->context_.config_.rpc_bx_shift_;
68 
69  // Check Ring
70  bool tp_is_substitute = (tp_ring == 3);
71 
72  // Calculate type
73  const bool tp_is_barrel = (tp_region == 0);
74 
75  rpc::Type tp_rpc_type;
76 
77  if ((!tp_is_barrel) && (tp_station >= 3) && (tp_ring == 1)) {
78  tp_rpc_type = rpc::Type::kiRPC;
79  } else {
80  tp_rpc_type = rpc::Type::kRPC;
81  }
82 
83  // Short-Circuit: Skip Barrel RPC (region = 0)
84  if (tp_region == 0) {
85  continue;
86  }
87 
88  // Short-Circuit: Skip Overlap region (RE1/3, RE2/3)
89  if (tp_station <= 2 && tp_ring == 3) {
90  continue;
91  }
92 
93  // Short-Circuit: Reject wide clusters
94  if (tp_rpc_type == rpc::Type::kiRPC) {
95  if (tp_clus_width > clus_width_cut_irpc) {
96  continue;
97  }
98  } else {
99  if (tp_clus_width > clus_width_cut) {
100  continue;
101  }
102  }
103 
104  // Calculate EMTF Info
105  int tp_chamber;
106 
107  if (tp_rpc_type == rpc::Type::kiRPC) {
108  tp_chamber = (tp_rpc_sector - 1) * 3 + tp_rpc_subsector;
109  } else {
110  tp_chamber = (tp_rpc_sector - 1) * 6 + tp_rpc_subsector;
111  }
112 
113  const int tp_sector = csc::getTriggerSector(tp_station, tp_ring, tp_chamber);
114  const int tp_subsector = csc::getTriggerSubsector(tp_station, tp_chamber);
115  const int tp_csc_id = csc::getId(tp_station, tp_ring, tp_chamber);
116  const auto tp_csc_facing = csc::getFaceDirection(tp_station, tp_ring, tp_chamber);
117 
118  // Assertion checks
119  emtf_assert(kMinEndcap <= tp_endcap && tp_endcap <= kMaxEndcap);
120  emtf_assert(kMinTrigSector <= tp_sector && tp_sector <= kMaxTrigSector);
121  emtf_assert(0 <= tp_subsector && tp_subsector <= 2);
122  emtf_assert(1 <= tp_station && tp_station <= 4);
123  emtf_assert(1 <= tp_chamber && tp_chamber <= 36);
124  emtf_assert((1 <= tp_csc_id) and (tp_csc_id <= 9));
125 
126  if (tp_rpc_type == rpc::Type::kiRPC) {
127  emtf_assert(tp_ring == 1);
128  emtf_assert(1 <= tp_roll && tp_roll <= 5);
129  emtf_assert(1 <= tp_strip && tp_strip <= 96);
130  } else {
131  emtf_assert(2 <= tp_ring && tp_ring <= 3);
132  emtf_assert(1 <= tp_roll && tp_roll <= 3);
133  emtf_assert(tp_is_CPPF || (1 <= tp_strip && tp_strip <= 32));
134  }
135 
136  emtf_assert(tp_data.valid);
137 
138  // Add info
139  tp_entry.info_.bx = tp_bx;
140 
141  tp_entry.info_.endcap = tp_endcap;
142  tp_entry.info_.endcap_pm = tp_endcap_pm;
143  tp_entry.info_.sector = tp_sector;
144  tp_entry.info_.subsector = tp_subsector;
145  tp_entry.info_.station = tp_station;
146  tp_entry.info_.ring = tp_ring;
147  tp_entry.info_.roll = tp_roll;
148  tp_entry.info_.layer = tp_layer;
149  tp_entry.info_.chamber = tp_chamber;
150 
151  tp_entry.info_.csc_id = tp_csc_id;
152  tp_entry.info_.csc_facing = tp_csc_facing;
153 
154  tp_entry.info_.rpc_type = tp_rpc_type;
155 
156  tp_entry.info_.flag_substitute = tp_is_substitute;
157 
158  bx_tpc_map[tp_bx].push_back(tp_entry);
159  }
160 }
int getId(int ring, int station, int chamber)
Definition: CSCUtils.cc:39
int getTriggerSubsector(int station, int chamber)
Definition: CSCUtils.cc:85
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:526
constexpr int kMaxEndcap
Definition: EMTFConstants.h:7
std::map< int, TPCollection > BXTPCMap
Definition: EMTFTypes.h:22
const EMTFContext & context_
#define emtf_assert(expr)
Definition: DebugTools.h:18
int getTriggerSector(int ring, int station, int chamber)
Definition: CSCUtils.cc:70
const edm::EDGetToken input_token_
Facing getFaceDirection(int station, int ring, int chamber)
Definition: CSCUtils.cc:100
RPCTPCollector(const EMTFContext &, edm::ConsumesCollector &)
constexpr int kMinTrigSector
Definition: EMTFConstants.h:20
constexpr int kMinEndcap
Definition: EMTFConstants.h:6
constexpr int kMaxTrigSector
Definition: EMTFConstants.h:21
std::vector< TPEntry > TPCollection
Definition: EMTFTypes.h:21
void collect(const edm::Event &, BXTPCMap &) const final
EMTFConfiguration config_
Definition: EMTFContext.h:39