CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
HGCalBestChoiceCodecImpl Class Reference

#include <HGCalBestChoiceCodecImpl.h>

Public Types

typedef HGCalBestChoiceDataPayload data_type
 

Public Member Functions

uint32_t adcnBits () const
 
double adcsaturation () const
 
void bestChoiceSelect (data_type &)
 
size_t dataLength () const
 
data_type decode (const std::vector< bool > &) const
 
std::vector< bool > encode (const data_type &) const
 
 HGCalBestChoiceCodecImpl (const edm::ParameterSet &conf)
 
void linearize (const std::vector< HGCDataFrame< DetId, HGCSample >> &, std::vector< std::pair< DetId, uint32_t > > &)
 
double linLSB () const
 
size_t nData () const
 
uint32_t tdcnBits () const
 
double tdcOnsetfC () const
 
double tdcsaturation () const
 
uint32_t triggerCellSaturationBits () const
 
void triggerCellSums (const HGCalTriggerGeometryBase &, const std::vector< std::pair< DetId, uint32_t > > &, data_type &)
 
uint32_t triggerCellTruncationBits () const
 

Private Attributes

double adcLSB_
 
uint32_t adcnBits_
 
double adcsaturation_
 
size_t dataLength_
 
double linLSB_
 
size_t nCellsInModule_
 
size_t nData_
 
double tdcLSB_
 
uint32_t tdcnBits_
 
double tdcOnsetfC_
 
double tdcsaturation_
 
uint32_t triggerCellSaturationBits_
 
uint32_t triggerCellTruncationBits_
 

Detailed Description

Definition at line 27 of file HGCalBestChoiceCodecImpl.h.

Member Typedef Documentation

Definition at line 30 of file HGCalBestChoiceCodecImpl.h.

Constructor & Destructor Documentation

HGCalBestChoiceCodecImpl::HGCalBestChoiceCodecImpl ( const edm::ParameterSet conf)

Definition at line 6 of file HGCalBestChoiceCodecImpl.cc.

References adcLSB_, adcnBits_, adcsaturation_, dataLength_, nCellsInModule_, nData_, funct::pow(), tdcLSB_, tdcnBits_, tdcsaturation_, triggerCellSaturationBits_, and triggerCellTruncationBits_.

6  :
7  nData_(conf.getParameter<uint32_t>("NData")),
8  dataLength_(conf.getParameter<uint32_t>("DataLength")),
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 {
19  // Cannot have more selected cells than the max number of cells
24 }
T getParameter(std::string const &) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40

Member Function Documentation

uint32_t HGCalBestChoiceCodecImpl::adcnBits ( ) const
inline
double HGCalBestChoiceCodecImpl::adcsaturation ( ) const
inline
void HGCalBestChoiceCodecImpl::bestChoiceSelect ( data_type data)

Definition at line 158 of file HGCalBestChoiceCodecImpl.cc.

References a, b, Exception, mps_fire::i, nCellsInModule_, nData_, and HGCalBestChoiceDataPayload::payload.

Referenced by HGCalBestChoiceCodec::setDataPayloadImpl().

160 {
161  // Store data payload in vector for energy sorting. Then refill the data payload after trigger
162  // cell selection.
163  // Probably not the most efficient way.
164  // Should check in the firmware how cells with the same energy are sorted
165 
166  // copy for sorting
167  std::vector< std::pair<uint32_t, uint32_t> > sortedtriggercells; // value, ID
168  sortedtriggercells.reserve(nCellsInModule_);
169  for(size_t i=0; i<nCellsInModule_; i++)
170  {
171  sortedtriggercells.push_back(std::make_pair(data.payload[i], i));
172  }
173  // sort, reverse order
174  sort(sortedtriggercells.begin(), sortedtriggercells.end(),
175  [](const std::pair<uint32_t, uint32_t>& a,
176  const std::pair<uint32_t, uint32_t>& b) -> bool
177  {
178  return a > b;
179  }
180  );
181  // keep only the first trigger cells
182  for(size_t i=nData_; i<nCellsInModule_; i++)
183  {
184  sortedtriggercells.at(i).first = 0;
185  }
186  for(const auto& value_id : sortedtriggercells)
187  {
188  if(value_id.second>=nCellsInModule_)
189  {
190  throw cms::Exception("BadGeometry")
191  << "Number of trigger cells in module too large for available data payload\n";
192  }
193  data.payload.at(value_id.second) = value_id.first;
194  }
195 }
double b
Definition: hdecay.h:120
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double a
Definition: hdecay.h:121
size_t HGCalBestChoiceCodecImpl::dataLength ( ) const
inline
HGCalBestChoiceCodecImpl::data_type HGCalBestChoiceCodecImpl::decode ( const std::vector< bool > &  data) const

Definition at line 62 of file HGCalBestChoiceCodecImpl.cc.

References b, EnergyCorrector::c, data, dataLength_, Exception, mps_fire::i, diffTreeTool::index, nCellsInModule_, nData_, HGCalBestChoiceDataPayload::payload, HGCalBestChoiceDataPayload::reset(), mps_fire::result, and relativeConstraints::value.

Referenced by HGCalBestChoiceCodec::decodeImpl().

64 {
66  result.reset();
68  {
69  throw cms::Exception("BadData")
70  << "decode: data length ("<<data.size()<<") inconsistent with codec parameters:\n"\
71  << " : Map size = "<<nCellsInModule_<<"\n"\
72  << " : Number of energy values = "<<nData_<<"\n"\
73  << " : Energy value length = "<<dataLength_<<"\n";
74  }
75  size_t c = 0;
76  for(size_t b=0; b<nCellsInModule_; b++)
77  {
78  if(data[b])
79  {
80  uint32_t value = 0;
81  for(size_t i=0;i<dataLength_;i++)
82  {
83  size_t index = nCellsInModule_+c*dataLength_+i;
84  if(data[index]) value |= (0x1<<i);
85  }
86  c++;
87  result.payload[b] = value;
88  }
89  }
90  return result;
91 }
Definition: value.py:1
HGCalBestChoiceDataPayload data_type
double b
Definition: hdecay.h:120
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
std::vector< bool > HGCalBestChoiceCodecImpl::encode ( const data_type data) const

Definition at line 29 of file HGCalBestChoiceCodecImpl.cc.

References data, dataLength_, Exception, mps_fire::i, nCellsInModule_, nData_, mps_fire::result, triggerCellSaturationBits_, and triggerCellTruncationBits_.

Referenced by HGCalBestChoiceCodec::encodeImpl().

31 {
32  // First nCellsInModule_ bits are encoding the map of selected trigger cells
33  // Followed by nData_ words of dataLength_ bits, corresponding to energy/transverse energy of
34  // the selected trigger cells
35  std::vector<bool> result(nCellsInModule_ + dataLength_*nData_, false);
36  size_t idata = 0;
37  for(size_t itc=0; itc<nCellsInModule_; itc++)
38  {
39  uint32_t value = data.payload.at(itc);
40  result[itc] = (value>0 ? 1 : 0);
41  if(value>0)
42  {
43  if(idata>=nData_)
44  {
45  throw cms::Exception("BadData")
46  << "encode: Number of non-zero trigger cells larger than codec parameter\n"\
47  << " : Number of energy values = "<<nData_<<"\n";
48  }
49  // Saturate and truncate energy values
50  if(value+1>(0x1u<<triggerCellSaturationBits_)) value = (0x1<<triggerCellSaturationBits_)-1;
51  for(size_t i=0; i<dataLength_; i++)
52  {
53  result[nCellsInModule_ + idata*dataLength_ + i] = static_cast<bool>(value & (0x1<<(i+triggerCellTruncationBits_)));// remove the lowest bits (=triggerCellTruncationBits_)
54  }
55  idata++;
56  }
57  }
58  return result;
59 }
Definition: value.py:1
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
void HGCalBestChoiceCodecImpl::linearize ( const std::vector< HGCDataFrame< DetId, HGCSample >> &  dataframes,
std::vector< std::pair< DetId, uint32_t > > &  linearized_dataframes 
)

Definition at line 95 of file HGCalBestChoiceCodecImpl.cc.

References adcLSB_, CustomPhysics_cfi::amplitude, data, linLSB_, ALCARECOPromptCalibProdSiPixelAli0T_cff::mode, tdcLSB_, and tdcOnsetfC_.

Referenced by HGCalBestChoiceCodec::setDataPayloadImpl().

98 {
99  double amplitude; uint32_t amplitude_int;
100 
101 
102  for(const auto& frame : dataframes) {//loop on DIGI
103  if (frame[2].mode()) {//TOT mode
104  amplitude =( floor(tdcOnsetfC_/adcLSB_) + 1.0 )* adcLSB_ + double(frame[2].data()) * tdcLSB_;
105  }
106  else {//ADC mode
107  amplitude = double(frame[2].data()) * adcLSB_;
108  }
109 
110  amplitude_int = uint32_t (floor(amplitude/linLSB_+0.5));
111  if (amplitude_int>65535) amplitude_int = 65535;
112 
113  linearized_dataframes.push_back(std::make_pair (frame.id(), amplitude_int));
114  }
115 }
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double HGCalBestChoiceCodecImpl::linLSB ( ) const
inline
size_t HGCalBestChoiceCodecImpl::nData ( ) const
inline

Definition at line 46 of file HGCalBestChoiceCodecImpl.h.

uint32_t HGCalBestChoiceCodecImpl::tdcnBits ( ) const
inline
double HGCalBestChoiceCodecImpl::tdcOnsetfC ( ) const
inline
double HGCalBestChoiceCodecImpl::tdcsaturation ( ) const
inline
uint32_t HGCalBestChoiceCodecImpl::triggerCellSaturationBits ( ) const
inline

Definition at line 55 of file HGCalBestChoiceCodecImpl.h.

void HGCalBestChoiceCodecImpl::triggerCellSums ( const HGCalTriggerGeometryBase geometry,
const std::vector< std::pair< DetId, uint32_t > > &  linearized_dataframes,
data_type data 
)

Definition at line 119 of file HGCalBestChoiceCodecImpl.cc.

References SoftLeptonByDistance_cfi::distance, Exception, HGCalTriggerGeometryBase::getModuleFromTriggerCell(), HGCalTriggerGeometryBase::getOrderedTriggerCellsFromModule(), HGCalTriggerGeometryBase::getTriggerCellFromCell(), nCellsInModule_, HGCalBestChoiceDataPayload::payload, jets_cff::payload, and relativeConstraints::value.

Referenced by HGCalBestChoiceCodec::setDataPayloadImpl().

121 {
122  if(linearized_dataframes.empty()) return;
123  std::map<HGCalDetId, uint32_t> payload;
124  // sum energies in trigger cells
125  for(const auto& frame : linearized_dataframes)
126  {
127  DetId cellid(frame.first);
128  // find trigger cell associated to cell
129  uint32_t tcid = geometry.getTriggerCellFromCell(cellid);
130  HGCalDetId triggercellid( tcid );
131  payload.insert( std::make_pair(triggercellid, 0) ); // do nothing if key exists already
132  // FIXME: need to transform ADC and TDC to the same linear scale on 12 bits
133  uint32_t value = frame.second; // 'value' has to be a 12 bit word
134  payload[triggercellid] += value; // 32 bits integer should be largely enough (maximum 7 12-bits sums are done)
135 
136  }
137  uint32_t module = geometry.getModuleFromTriggerCell(payload.begin()->first);
138  HGCalTriggerGeometryBase::geom_ordered_set trigger_cells_in_module = geometry.getOrderedTriggerCellsFromModule(module);
139  // fill data payload
140  for(const auto& id_value : payload)
141  {
142  // find the index of the trigger cell in the module (not necessarily equal to .cell())
143  // FIXME: std::distance is linear with size for sets (no random access). In order to have constant
144  // access would require to convert the set into a vector.
145  uint32_t id = std::distance(trigger_cells_in_module.begin(),trigger_cells_in_module.find(id_value.first));
146  //uint32_t id = id_value.first.cell();
147  //std::cerr<<"cell id in trigger cell sum: "<<id<<"("<<id_value.first.wafer()<<","<<id_value.first.cell()<<")\n";
148  if(id>=nCellsInModule_)
149  {
150  throw cms::Exception("BadGeometry")
151  << "Number of trigger cells in module too large for available data payload\n";
152  }
153  data.payload.at(id) = id_value.second;
154  }
155 }
virtual geom_ordered_set getOrderedTriggerCellsFromModule(const unsigned trigger_cell_det_id) const =0
Definition: value.py:1
Definition: DetId.h:18
virtual unsigned getModuleFromTriggerCell(const unsigned trigger_cell_det_id) const =0
std::set< unsigned > geom_ordered_set
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
Definition: vlib.h:208
virtual unsigned getTriggerCellFromCell(const unsigned cell_det_id) const =0
uint32_t HGCalBestChoiceCodecImpl::triggerCellTruncationBits ( ) const
inline

Member Data Documentation

double HGCalBestChoiceCodecImpl::adcLSB_
private

Definition at line 68 of file HGCalBestChoiceCodecImpl.h.

Referenced by HGCalBestChoiceCodecImpl(), and linearize().

uint32_t HGCalBestChoiceCodecImpl::adcnBits_
private

Definition at line 64 of file HGCalBestChoiceCodecImpl.h.

Referenced by HGCalBestChoiceCodecImpl().

double HGCalBestChoiceCodecImpl::adcsaturation_
private

Definition at line 63 of file HGCalBestChoiceCodecImpl.h.

Referenced by HGCalBestChoiceCodecImpl().

size_t HGCalBestChoiceCodecImpl::dataLength_
private

Definition at line 60 of file HGCalBestChoiceCodecImpl.h.

Referenced by decode(), encode(), and HGCalBestChoiceCodecImpl().

double HGCalBestChoiceCodecImpl::linLSB_
private

Definition at line 62 of file HGCalBestChoiceCodecImpl.h.

Referenced by linearize().

size_t HGCalBestChoiceCodecImpl::nCellsInModule_
private
size_t HGCalBestChoiceCodecImpl::nData_
private
double HGCalBestChoiceCodecImpl::tdcLSB_
private

Definition at line 69 of file HGCalBestChoiceCodecImpl.h.

Referenced by HGCalBestChoiceCodecImpl(), and linearize().

uint32_t HGCalBestChoiceCodecImpl::tdcnBits_
private

Definition at line 66 of file HGCalBestChoiceCodecImpl.h.

Referenced by HGCalBestChoiceCodecImpl().

double HGCalBestChoiceCodecImpl::tdcOnsetfC_
private

Definition at line 67 of file HGCalBestChoiceCodecImpl.h.

Referenced by linearize().

double HGCalBestChoiceCodecImpl::tdcsaturation_
private

Definition at line 65 of file HGCalBestChoiceCodecImpl.h.

Referenced by HGCalBestChoiceCodecImpl().

uint32_t HGCalBestChoiceCodecImpl::triggerCellSaturationBits_
private

Definition at line 71 of file HGCalBestChoiceCodecImpl.h.

Referenced by encode(), and HGCalBestChoiceCodecImpl().

uint32_t HGCalBestChoiceCodecImpl::triggerCellTruncationBits_
private

Definition at line 70 of file HGCalBestChoiceCodecImpl.h.

Referenced by encode(), and HGCalBestChoiceCodecImpl().