CMS 3D CMS Logo

GEMTPCollector.cc
Go to the documentation of this file.
3 
11 
12 using namespace emtf::phase2;
13 
15  : context_(context),
16  input_token_(i_consumes_collector.consumes<GEMTag::collection_type>(context.config_.gem_input_)) {}
17 
18 void GEMTPCollector::collect(const edm::Event& i_event, BXTPCMap& bx_tpc_map) const {
19  // Constants
20  static const int max_delta_roll = 1;
21  static const int max_delta_pad_ge11 = 4;
22  static const int max_delta_pad_ge21 = 4;
23 
24  // Read GEM digis
25  TPCollection tpc;
26 
28  i_event.getByToken(input_token_, gem_digis);
29 
30  auto chamber = gem_digis->begin();
31  auto chend = gem_digis->end();
32 
33  for (; chamber != chend; ++chamber) {
34  auto digi = (*chamber).second.first;
35  auto dend = (*chamber).second.second;
36 
37  for (; digi != dend; ++digi) {
38  // Short-Circuit: Ignore invalid digis
39  bool tp_valid = (*digi).isValid();
40 
41  if (!tp_valid) {
42  continue;
43  }
44 
45  // Append digi
46  tpc.emplace_back((*chamber).first, *digi);
47  }
48  }
49 
50  // Find Copads
51  std::map<std::pair<uint32_t, uint16_t>, std::vector<std::array<uint16_t, 3>>> chamber_copads_map;
52 
53  for (auto& tp_entry : tpc) {
54  const auto& tp_det_id = tp_entry.tp_.detId<GEMDetId>();
55  const GEMData& tp_data = tp_entry.tp_.getGEMData();
56 
57  const int tp_region = tp_det_id.region(); // 0: barrel, +/-1: endcap
58  const int tp_station = tp_det_id.station();
59  const int tp_ring = tp_det_id.ring();
60  const int tp_chamber = tp_det_id.chamber();
61  const int tp_layer = tp_det_id.layer();
62 
63  const uint16_t tp_roll = tp_det_id.roll();
64  const uint16_t tp_pad_lo = tp_data.pad_low;
65  const uint16_t tp_pad_hi = tp_data.pad_hi;
66 
67  const int tp_bx = tp_data.bx + this->context_.config_.gem_bx_shift_;
68 
69  GEMDetId tp_mod_det_id(tp_region, tp_ring, tp_station, 0, tp_chamber, 0);
70  auto key = std::make_pair(tp_mod_det_id.rawId(), tp_bx);
71 
72  if (tp_layer == 1) {
73  // Layer 1 is incidence
74  // If key does not exist, insert an empty vector. If key exists, do nothing.
75  decltype(chamber_copads_map)::mapped_type copads;
76  chamber_copads_map.insert({key, copads});
77  } else if (tp_layer == 2) {
78  // Layer 2 is coincidence
79  decltype(chamber_copads_map)::mapped_type::value_type copad{{tp_roll, tp_pad_lo, tp_pad_hi}};
80  chamber_copads_map[key].push_back(copad);
81  }
82  }
83 
84  // Map to BX
85  for (auto& tp_entry : tpc) {
86  const auto& tp_det_id = tp_entry.tp_.detId<GEMDetId>();
87  const GEMData& tp_data = tp_entry.tp_.getGEMData();
88 
89  const int tp_region = tp_det_id.region(); // 0: barrel, +/-1: endcap
90  const int tp_endcap = (tp_region == -1) ? 2 : tp_region; // 1: +endcap, 2: -endcap
91  const int tp_endcap_pm = (tp_endcap == 2) ? -1 : tp_endcap; // 1: +endcap, -1: -endcap
92  const int tp_station = tp_det_id.station();
93  const int tp_ring = tp_det_id.ring();
94  const int tp_layer = tp_det_id.layer();
95  const int tp_chamber = tp_det_id.chamber();
96 
97  const int tp_roll = tp_det_id.roll();
98  const int tp_pad = (tp_data.pad_low + tp_data.pad_hi) / 2;
99  const int tp_pad_lo = tp_data.pad_low;
100  const int tp_pad_hi = tp_data.pad_hi;
101 
102  const int tp_bx = tp_data.bx + this->context_.config_.gem_bx_shift_;
103 
104  // Get Copad Info
105  GEMDetId tp_mod_det_id(tp_region, tp_ring, tp_station, 0, tp_chamber, 0);
106  auto tp_copads_key = std::make_pair(tp_mod_det_id.rawId(), tp_bx);
107  auto tp_copads = chamber_copads_map.at(tp_copads_key);
108 
109  // Check Copads
110  bool tp_is_substitute = false;
111 
112  if (tp_layer == 1) {
113  // layer 1 is used as incidence
114  const bool is_ge21 = (tp_station == 2);
115 
116  auto match_fn = [&tp_roll, &tp_pad_lo, &tp_pad_hi, &is_ge21](const std::array<uint16_t, 3>& elem) {
117  // Unpack entry
118  // Compare roll and (pad_lo, pad_hi)-range with tolerance
119  const auto& [c_roll_tmp, c_pad_lo_tmp, c_pad_hi_tmp] = elem;
120  int c_roll_lo = static_cast<int>(c_roll_tmp) - max_delta_roll;
121  int c_roll_hi = static_cast<int>(c_roll_tmp) + max_delta_roll;
122  int c_pad_lo = static_cast<int>(c_pad_lo_tmp) - (is_ge21 ? max_delta_pad_ge21 : max_delta_pad_ge11);
123  int c_pad_hi = static_cast<int>(c_pad_hi_tmp) + (is_ge21 ? max_delta_pad_ge21 : max_delta_pad_ge11);
124 
125  // Two ranges overlap if (range_a_lo <= range_b_hi) and (range_a_hi >= range_b_lo)
126  return (tp_roll <= c_roll_hi) and (tp_roll >= c_roll_lo) and (tp_pad_lo <= c_pad_hi) and
127  (tp_pad_hi >= c_pad_lo);
128  };
129 
130  auto match = std::find_if(tp_copads.begin(), tp_copads.end(), match_fn);
131 
132  if (match != tp_copads.end()) {
133  // Has copad
134  tp_is_substitute = false;
135  } else if (tp_copads.empty()) {
136  // Kinda has copad
137  tp_is_substitute = true;
138  } else {
139  // Short-Circuit: Didn't find coincidence
140  continue;
141  }
142  } else if (tp_layer == 2) {
143  // Short-Circuit: layer 2 is used as coincidence
144  continue;
145  }
146 
147  // Calculate EMTF Info
148  const int tp_sector = csc::getTriggerSector(tp_station, tp_ring, tp_chamber);
149  const int tp_subsector = csc::getTriggerSubsector(tp_station, tp_chamber);
150  const int tp_csc_id = csc::getId(tp_station, tp_ring, tp_chamber);
151  const auto tp_csc_facing = csc::getFaceDirection(tp_station, tp_ring, tp_chamber);
152 
153  // Assertion checks
154  emtf_assert(kMinEndcap <= tp_endcap && tp_endcap <= kMaxEndcap);
155  emtf_assert(kMinTrigSector <= tp_sector && tp_sector <= kMaxTrigSector);
156  emtf_assert((0 <= tp_subsector) and (tp_subsector <= 2));
157  emtf_assert(1 <= tp_station && tp_station <= 2);
158  emtf_assert(tp_ring == 1);
159  emtf_assert((1 <= tp_chamber) and (tp_chamber <= 36));
160  emtf_assert(1 <= tp_csc_id && tp_csc_id <= 3);
161  emtf_assert(tp_station == 1 or (1 <= tp_roll && tp_roll <= 16));
162  emtf_assert(tp_station != 1 or (1 <= tp_roll && tp_roll <= 8));
163  emtf_assert(1 <= tp_layer && tp_layer <= 2);
164  emtf_assert((tp_station == 1 && 0 <= tp_pad && tp_pad <= 191) || (tp_station != 1));
165  emtf_assert((tp_station == 2 && 0 <= tp_pad && tp_pad <= 383) || (tp_station != 2));
166 
167  // Add info
168  tp_entry.info_.bx = tp_bx;
169 
170  tp_entry.info_.endcap = tp_endcap;
171  tp_entry.info_.endcap_pm = tp_endcap_pm;
172  tp_entry.info_.sector = tp_sector;
173  tp_entry.info_.subsector = tp_subsector;
174  tp_entry.info_.station = tp_station;
175  tp_entry.info_.ring = tp_ring;
176  tp_entry.info_.roll = tp_roll;
177  tp_entry.info_.layer = tp_layer;
178  tp_entry.info_.chamber = tp_chamber;
179 
180  tp_entry.info_.csc_id = tp_csc_id;
181  tp_entry.info_.csc_facing = tp_csc_facing;
182 
183  tp_entry.info_.flag_substitute = tp_is_substitute;
184 
185  bx_tpc_map[tp_bx].push_back(tp_entry);
186  }
187 }
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:528
constexpr int kMaxEndcap
Definition: EMTFConstants.h:7
const edm::EDGetToken input_token_
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
std::map< int, TPCollection > BXTPCMap
Definition: EMTFTypes.h:22
key
prepare the HTCondor submission files and eventually submit them
void collect(const edm::Event &, BXTPCMap &) const final
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
GEMTPCollector(const EMTFContext &, edm::ConsumesCollector &)
#define emtf_assert(expr)
Definition: DebugTools.h:18
int getTriggerSector(int ring, int station, int chamber)
Definition: CSCUtils.cc:70
Facing getFaceDirection(int station, int ring, int chamber)
Definition: CSCUtils.cc:100
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
const EMTFContext & context_
EMTFConfiguration config_
Definition: EMTFContext.h:39