CMS 3D CMS Logo

HGCalTriggerCellBestChoiceCodecImpl.cc
Go to the documentation of this file.
1 
3 
4 
7  nData_(conf.getParameter<uint32_t>("NData")),
8  dataLength_(conf.getParameter<uint32_t>("DataLength")),
9  nCellsInModule_(conf.getParameter<uint32_t>("MaxCellsInModule")),
10  linLSB_(conf.getParameter<double>("linLSB")),
11  adcsaturation_(conf.getParameter<double>("adcsaturation")),
12  adcnBits_(conf.getParameter<uint32_t>("adcnBits")),
13  tdcsaturation_(conf.getParameter<double>("tdcsaturation")),
14  tdcnBits_(conf.getParameter<uint32_t>("tdcnBits")),
15  tdcOnsetfC_(conf.getParameter<double>("tdcOnsetfC")),
16  triggerCellTruncationBits_(conf.getParameter<uint32_t>("triggerCellTruncationBits"))
17 {
18  // Cannot have more selected cells than the max number of cells
23 }
24 
25 
26 
27 std::vector<bool>
30 {
31  // First nCellsInModule_ bits are encoding the map of selected trigger cells
32  // Followed by nData_ words of dataLength_ bits, corresponding to energy/transverse energy of
33  // the selected trigger cells
34  std::vector<bool> result(nCellsInModule_ + dataLength_*nData_, 0);
35  // No data: return vector of 0
36  if(data.payload.size()==0) return result;
37  // All trigger cells are in the same module
38  // Retrieve once the ordered list of trigger cells in this module
39  uint32_t module = geometry.getModuleFromTriggerCell(data.payload.begin()->detId());
40  HGCalTriggerGeometryBase::geom_ordered_set trigger_cells_in_module = geometry.getOrderedTriggerCellsFromModule(module);
41  // Convert payload into a map for later search
42  std::unordered_map<uint32_t, uint32_t> data_map; // (detid,energy)
43  for(const auto& triggercell : data.payload)
44  {
45  data_map.emplace(triggercell.detId(), triggercell.hwPt());
46  }
47  // Loop on trigger cell ids in module and check if energy in the cell
48  size_t index = 0; // index in module
49  size_t idata = 0; // counter for the number of non-zero energy values
50  for(const auto& triggercell_id : trigger_cells_in_module)
51  {
52  // Find if this trigger cell has data
53  const auto& data_itr = data_map.find(triggercell_id);
54  // if not data, increase index and skip
55  if(data_itr==data_map.end())
56  {
57  index++;
58  continue;
59  }
60  // else fill result vector with data
61  // (set the corresponding adress bit and fill energy if >0)
62  if(index>=nCellsInModule_)
63  {
64  throw cms::Exception("BadGeometry")
65  << "Number of trigger cells in module too large for available data payload\n";
66  }
67  uint32_t value = data_itr->second;
68  if(value>0)
69  {
70  if(idata>=nData_)
71  {
72  throw cms::Exception("BadData")
73  << "encode: Number of non-zero trigger cells larger than codec parameter\n"\
74  << " : Number of energy values = "<<nData_<<"\n";
75  }
76  // Set map bit to 1
77  result[index] = 1;
78  // Saturate and truncate energy values
79  if(value+1>(0x1u<<triggerCellSaturationBits_)) value = (0x1<<triggerCellSaturationBits_)-1;
80  for(size_t i=0; i<dataLength_; i++)
81  {
82  // remove the lowest bits (=triggerCellTruncationBits_)
83  result[nCellsInModule_ + idata*dataLength_ + i] = static_cast<bool>(value & (0x1<<(i+triggerCellTruncationBits_)));
84  }
85  idata++;
86  }
87  index++;
88  }
89  return result;
90 }
91 
94 decode(const std::vector<bool>& data, const uint32_t module, const HGCalTriggerGeometryBase& geometry) const
95 {
97  result.reset();
98  // TODO: could eventually reserve result memory to the max size of trigger cells
99  if(data.size()!=nCellsInModule_+dataLength_*nData_)
100  {
101  throw cms::Exception("BadData")
102  << "decode: data length ("<<data.size()<<") inconsistent with codec parameters:\n"\
103  << " : Map size = "<<nCellsInModule_<<"\n"\
104  << " : Number of energy values = "<<nData_<<"\n"\
105  << " : Energy value length = "<<dataLength_<<"\n";
106  }
107  HGCalTriggerGeometryBase::geom_ordered_set trigger_cells_in_module = geometry.getOrderedTriggerCellsFromModule(module);
108  size_t iselected = 0;
109  size_t index = 0;
110  for(const auto& triggercell : trigger_cells_in_module)
111  {
112  if(index>=nCellsInModule_)
113  {
114  throw cms::Exception("BadGeometry")
115  << "Number of trigger cells in module too large for available data payload\n";
116  }
117  if(data[index])
118  {
119  uint32_t value = 0;
120  for(size_t i=0;i<dataLength_;i++)
121  {
122  size_t ibit = nCellsInModule_+iselected*dataLength_+i;
123  if(data[ibit]) value |= (0x1<<i);
124  }
125  iselected++;
126  // Build trigger cell
127  if(value>0)
128  {
129  // Currently no hardware eta, phi and quality values
130  result.payload.emplace_back(reco::LeafCandidate::LorentzVector(),
131  value, 0, 0, 0, triggercell);
132  GlobalPoint point = geometry.getTriggerCellPosition(triggercell);
133  // 'value' is hardware, so p4 is meaningless, except for eta and phi
134  math::PtEtaPhiMLorentzVector p4((double)value/cosh(point.eta()), point.eta(), point.phi(), 0.);
135  result.payload.back().setP4(p4);
136  result.payload.back().setPosition(point);
137  }
138  }
139  index++;
140  }
141  return result;
142 }
143 
144 
145 void
147 linearize(const std::vector<HGCDataFrame<HGCalDetId,HGCSample>>& dataframes,
148  std::vector<std::pair<HGCalDetId, uint32_t > >& linearized_dataframes)
149 {
150  double amplitude; uint32_t amplitude_int;
151 
152 
153  for(const auto& frame : dataframes) {//loop on DIGI
154  if (frame[2].mode()) {//TOT mode
155  amplitude =( floor(tdcOnsetfC_/adcLSB_) + 1.0 )* adcLSB_ + double(frame[2].data()) * tdcLSB_;
156  }
157  else {//ADC mode
158  amplitude = double(frame[2].data()) * adcLSB_;
159  }
160 
161  amplitude_int = uint32_t (floor(amplitude/linLSB_+0.5));
162  if (amplitude_int>65535) amplitude_int = 65535;
163 
164  linearized_dataframes.push_back(std::make_pair (frame.id(), amplitude_int));
165  }
166 }
167 
168 
169 void
171 triggerCellSums(const HGCalTriggerGeometryBase& geometry, const std::vector<std::pair<HGCalDetId, uint32_t > >& linearized_dataframes, data_type& data)
172 {
173  if(linearized_dataframes.size()==0) return;
174  std::map<HGCalDetId, uint32_t> payload;
175  // sum energies in trigger cells
176  for(const auto& frame : linearized_dataframes)
177  {
178  HGCalDetId cellid(frame.first);
179  // find trigger cell associated to cell
180  uint32_t tcid = geometry.getTriggerCellFromCell(cellid);
181  HGCalDetId triggercellid( tcid );
182  payload.insert( std::make_pair(triggercellid, 0) ); // do nothing if key exists already
183  uint32_t value = frame.second;
184  payload[triggercellid] += value; // 32 bits integer should be largely enough
185 
186  }
187  uint32_t module = geometry.getModuleFromTriggerCell(payload.begin()->first);
188  HGCalTriggerGeometryBase::geom_ordered_set trigger_cells_in_module = geometry.getOrderedTriggerCellsFromModule(module);
189  // fill data payload
190  for(const auto& id_value : payload)
191  {
192  // Store only energy value and detid
193  // No need position here
194  data.payload.emplace_back(reco::LeafCandidate::LorentzVector(),
195  id_value.second, 0, 0, 0, id_value.first.rawId());
196  }
197 }
198 
199 void
202 {
203  // sort, reverse order
204  sort(data.payload.begin(), data.payload.end(),
205  [](const l1t::HGCalTriggerCell& a,
206  const l1t::HGCalTriggerCell& b) -> bool
207  {
208  return a.hwPt() > b.hwPt();
209  }
210  );
211  // keep only the first trigger cells
212  if(data.payload.size()>nData_) data.payload.resize(nData_);
213 }
214 
215 
void triggerCellSums(const HGCalTriggerGeometryBase &, const std::vector< std::pair< HGCalDetId, uint32_t > > &, data_type &)
virtual geom_ordered_set getOrderedTriggerCellsFromModule(const unsigned trigger_cell_det_id) const =0
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
void linearize(const std::vector< HGCDataFrame< HGCalDetId, HGCSample >> &, std::vector< std::pair< HGCalDetId, uint32_t > > &)
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
payload
payload postfix for testing
data_type decode(const std::vector< bool > &, const uint32_t, const HGCalTriggerGeometryBase &) const
double p4[4]
Definition: TauolaWrapper.h:92
Definition: value.py:1
int hwPt() const
Definition: L1Candidate.h:48
virtual unsigned getModuleFromTriggerCell(const unsigned trigger_cell_det_id) const =0
HGCalTriggerCellBestChoiceCodecImpl(const edm::ParameterSet &conf)
double b
Definition: hdecay.h:120
std::set< unsigned > geom_ordered_set
T eta() const
Definition: PV3DBase.h:76
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double a
Definition: hdecay.h:121
Definition: vlib.h:208
virtual unsigned getTriggerCellFromCell(const unsigned cell_det_id) const =0
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
std::vector< bool > encode(const data_type &, const HGCalTriggerGeometryBase &) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
virtual GlobalPoint getTriggerCellPosition(const unsigned trigger_cell_det_id) const =0