CMS 3D CMS Logo

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