CMS 3D CMS Logo

HGCalConcentratorAutoEncoderImpl.cc
Go to the documentation of this file.
3 #include <iomanip>
4 
5 // Following example of implementing graphloading from here:
6 // https://gitlab.cern.ch/mrieger/CMSSW-TensorFlowExamples/-/blob/master/GraphLoading/
7 
9  : cellRemap_(conf.getParameter<std::vector<int>>("cellRemap")),
10  cellRemapNoDuplicates_(conf.getParameter<std::vector<int>>("cellRemapNoDuplicates")),
11  encoderShape_(conf.getParameter<std::vector<uint>>("encoderShape")),
12  decoderShape_(conf.getParameter<std::vector<uint>>("decoderShape")),
13  bitsPerInput_(conf.getParameter<int>("nBitsPerInput")),
14  maxBitsPerOutput_(conf.getParameter<int>("maxBitsPerOutput")),
15  outputBitsPerLink_(conf.getParameter<std::vector<int>>("bitsPerLink")),
16  modelFilePaths_(conf.getParameter<std::vector<edm::ParameterSet>>("modelFiles")),
17  linkToGraphMap_(conf.getParameter<std::vector<unsigned int>>("linkToGraphMap")),
18  zeroSuppresionThreshold_(conf.getParameter<double>("zeroSuppresionThreshold")),
19  bitShiftNormalization_(conf.getParameter<bool>("bitShiftNormalization")),
20  saveEncodedValues_(conf.getParameter<bool>("saveEncodedValues")),
21  preserveModuleSum_(conf.getParameter<bool>("preserveModuleSum")) {
22  // find total size of the expected input shape
23  // used for checking the maximum size used in cell Remap
24  nInputs_ = 1;
25  for (const auto& i : encoderShape_) {
26  nInputs_ *= i;
27  }
28 
29  // check the size of the inputs shapes
30  if (encoderShape_.size() != encoderTensorDims_) {
31  throw cms::Exception("BadInitialization")
32  << "Encoder input shapes are currently expected to be " << encoderTensorDims_ << " values";
33  }
34 
35  if (decoderShape_.size() != decoderTensorDims_) {
36  throw cms::Exception("BadInitialization")
37  << "Encoder input shapes are currently expected to be " << decoderTensorDims_ << " values long";
38  }
39 
40  if (cellRemap_.size() != nInputs_) {
41  throw cms::Exception("BadInitialization")
42  << "Size of cellRemap (" << cellRemap_.size()
43  << ") does not agree with the total size specified for the encoder inputs based on the encoderShape variable ("
44  << nInputs_ << ")";
45  }
46 
47  if (cellRemap_.size() != cellRemapNoDuplicates_.size()) {
48  throw cms::Exception("BadInitialization")
49  << "Size of cellRemap (" << cellRemap_.size() << ") does not agree with size of cellRemapNoDuplicates ("
50  << cellRemapNoDuplicates_.size() << ")";
51  }
52 
53  for (unsigned i = 0; i < cellRemap_.size(); i++) {
54  if (cellRemap_[i] > nTriggerCells_ - 1) {
55  throw cms::Exception("BadInitialization")
56  << "cellRemap value " << cellRemap_[i] << " is larger than the number of trigger cells " << nTriggerCells_;
57  }
59  throw cms::Exception("BadInitialization") << "cellRemapNoDuplicates value " << cellRemapNoDuplicates_[i]
60  << " is larger than the number of trigger cells " << nTriggerCells_;
61  }
62  }
63 
65 
66  for (const auto& modelFilePset : modelFilePaths_) {
67  std::string encoderPath = modelFilePset.getParameter<edm::FileInPath>("encoderModelFile").fullPath();
68  std::string decoderPath = modelFilePset.getParameter<edm::FileInPath>("decoderModelFile").fullPath();
69 
70  graphDef_encoder_ = std::unique_ptr<tensorflow::GraphDef>{tensorflow::loadGraphDef(encoderPath)};
71 
72  // create a new session and add the graphDef
73  session_encoder_.push_back(
74  std::unique_ptr<tensorflow::Session>{tensorflow::createSession(graphDef_encoder_.get())});
75 
76  graphDef_decoder_ = std::unique_ptr<tensorflow::GraphDef>{tensorflow::loadGraphDef(decoderPath)};
77 
78  // create a new session and add the graphDef
79  session_decoder_.push_back(
80  std::unique_ptr<tensorflow::Session>{tensorflow::createSession(graphDef_decoder_.get())});
81 
82  //extract encoder tenser names from first graph, check that rest of the names are consistent
83  if (modelFilePset == modelFilePaths_.front()) {
84  inputTensorName_encoder_ = graphDef_encoder_.get()->node(0).name();
85  outputTensorName_encoder_ = graphDef_encoder_.get()->node(graphDef_encoder_.get()->node_size() - 1).name();
86  inputTensorName_decoder_ = graphDef_decoder_.get()->node(0).name();
87  outputTensorName_decoder_ = graphDef_decoder_.get()->node(graphDef_decoder_.get()->node_size() - 1).name();
88  } else {
89  if (inputTensorName_encoder_ != graphDef_encoder_.get()->node(0).name()) {
90  throw cms::Exception("BadInitialization") << "provided list of encoder graphs have different input nodes";
91  }
92  if (outputTensorName_encoder_ != graphDef_encoder_.get()->node(graphDef_encoder_.get()->node_size() - 1).name()) {
93  throw cms::Exception("BadInitialization") << "provided list of encoder graphs have different output nodes";
94  }
95  if (inputTensorName_decoder_ != graphDef_decoder_.get()->node(0).name()) {
96  throw cms::Exception("BadInitialization") << "provided list of decoder graphs have different input nodes";
97  }
98  if (outputTensorName_decoder_ != graphDef_decoder_.get()->node(graphDef_decoder_.get()->node_size() - 1).name()) {
99  throw cms::Exception("BadInitialization") << "provided list of decoder graphs have different output nodes";
100  }
101  }
102  }
103 
104  // check that the appropriate number of links have been specified
105  if (linkToGraphMap_.size() <= maxNumberOfLinks_) {
106  throw cms::Exception("BadInitialization")
107  << "Autoencoder graph number must be specified for all link allocation possibilities. Only "
108  << linkToGraphMap_.size() << " values specified while " << maxNumberOfLinks_ << "links are possible";
109  }
110 
111  // check that all graph indices specified exist in the model file lists
112  for (const auto& graphNumber : linkToGraphMap_) {
113  if (graphNumber >= modelFilePaths_.size()) {
114  throw cms::Exception("BadInitialization")
115  << "Autoencoder graph number " << graphNumber << " is larger than the size of the provided list of graphs "
116  << modelFilePaths_.size();
117  }
118  }
119 }
120 
122  const std::vector<l1t::HGCalTriggerCell>& trigCellVecInput,
123  std::vector<l1t::HGCalTriggerCell>& trigCellVecOutput,
124  std::vector<l1t::HGCalConcentratorData>& ae_encodedLayer_Output) {
125  std::array<double, nTriggerCells_> mipPt;
126  std::array<double, nTriggerCells_> uncompressedCharge;
127  std::array<double, nTriggerCells_> compressedCharge;
128  std::array<double, maxAEInputSize_> ae_inputArray;
129  std::array<double, nTriggerCells_> ae_outputArray;
130 
131  //reset inputs to 0 to account for zero suppressed trigger cells
132  mipPt.fill(0);
133  uncompressedCharge.fill(0);
134  compressedCharge.fill(0);
135  ae_inputArray.fill(0);
136  ae_outputArray.fill(0);
137 
138  double modSum = 0;
139 
140  int bitsPerOutput = outputBitsPerLink_.at(nLinks);
141 
142  // largest expected input and output values, used for bit truncation
143  // values of -1 for the number of bits used to keep full precision, in which case the MaxIntSize variables are not used
144  double inputMaxIntSize = 1;
145  if (bitsPerInput_ > 0)
146  inputMaxIntSize = 1 << bitsPerInput_;
147  double outputMaxIntSize = 1;
148  if (bitsPerOutput > 0)
149  outputMaxIntSize = 1 << bitsPerOutput;
150  double outputMaxIntSizeGlobal = 1;
151  if (maxBitsPerOutput_ > 0)
152  outputMaxIntSizeGlobal = 1 << maxBitsPerOutput_;
153 
154  for (const auto& trigCell : trigCellVecInput) {
155  if (triggerTools_.isScintillator(trigCell.detId()))
156  return; //currently, only silicon modules are setup to work (mapping of scinillators would be different, and needs to be investigated)
157 
158  HGCalTriggerDetId id(trigCell.detId());
159  uint cellu = id.triggerCellU();
160  uint cellv = id.triggerCellV();
161  int inputIndex = cellUVremap_[cellu][cellv];
162  if (inputIndex < 0) {
163  throw cms::Exception("BadInitialization")
164  << "Invalid index provided for trigger cell u=" << cellu << " v=" << cellv << " in cellUVRemap[" << cellu
165  << "][" << cellv << "]";
166  }
167 
168  mipPt[inputIndex] = trigCell.mipPt();
169  uncompressedCharge[inputIndex] = trigCell.uncompressedCharge();
170  compressedCharge[inputIndex] = trigCell.compressedCharge();
171 
172  modSum += trigCell.mipPt();
173  }
174 
175  double normalization = modSum;
176  if (modSum > 0) {
177  //Use a bit shift normalization like will be implemented in ECON, rather than floating point sum
178  //Normalizes to the MSB value of the module sum
180  int msb = int(log2(modSum));
181  normalization = pow(2, msb);
182  }
183 
184  //normalize inputs to module sum
185  for (uint i = 0; i < nInputs_; i++) {
186  int remapIndex = cellRemap_[i];
187  if (remapIndex < 0)
188  continue;
189  ae_inputArray[i] = mipPt[remapIndex] / normalization;
190  //round to precision of input, if bitsPerInput_ is -1 keep full precision
191  if (bitsPerInput_ > 0) {
192  ae_inputArray[i] = std::round(ae_inputArray[i] * inputMaxIntSize) / inputMaxIntSize;
193  }
194  }
195  }
196 
197  tensorflow::Tensor encoder_input(tensorflow::DT_FLOAT,
199 
200  float* d = encoder_input.flat<float>().data();
201  for (uint i = 0; i < nInputs_; i++, d++) {
202  *d = ae_inputArray[i];
203  }
204 
205  int graphIndex = linkToGraphMap_.at(nLinks);
206 
207  std::vector<tensorflow::Tensor> encoder_outputs;
208  tensorflow::run(session_encoder_.at(graphIndex).get(),
209  {{inputTensorName_encoder_, encoder_input}},
211  &encoder_outputs);
212 
213  if (encoder_outputs.empty()) {
214  throw cms::Exception("BadInitialization") << "Autoencoder graph returning empty output vector";
215  }
216 
217  d = encoder_outputs[0].flat<float>().data();
218  for (int i = 0; i < encoder_outputs[0].NumElements(); i++, d++) {
219  ae_encodedLayer_[i] = *d;
220  //truncate the encoded layer bits
221  if (bitsPerOutput > 0 && maxBitsPerOutput_ > 0) {
222  ae_encodedLayer_[i] = std::round(ae_encodedLayer_[i] * outputMaxIntSize) / outputMaxIntSize;
223  }
224  }
225 
226  tensorflow::Tensor decoder_input(tensorflow::DT_FLOAT, {decoderShape_[0], decoderShape_[1]});
227  d = decoder_input.flat<float>().data();
228  for (int i = 0; i < nEncodedLayerNodes_; i++, d++) {
229  *d = ae_encodedLayer_[i];
230  }
231 
232  std::vector<tensorflow::Tensor> decoder_outputs;
233  tensorflow::run(session_decoder_.at(graphIndex).get(),
234  {{inputTensorName_decoder_, decoder_input}},
236  &decoder_outputs);
237 
238  double outputSum = 0.;
239 
240  d = decoder_outputs[0].flat<float>().data();
241  for (uint i = 0; i < nInputs_; i++, d++) {
242  int remapIndex = cellRemapNoDuplicates_[i];
243  if (remapIndex < 0)
244  continue;
245  outputSum += *d * normalization;
246  ae_outputArray[remapIndex] = *d;
247  }
248 
249  double renormalizationFactor = 1.;
250  if (preserveModuleSum_) {
251  renormalizationFactor = modSum / outputSum;
252  }
253 
254  // Add data back into trigger cells
255  if (modSum > 0) {
256  //get detID for everything but cell, take first entry detID and subtract off cellU and cellV contribution
257  HGCalTriggerDetId id(trigCellVecInput.at(0).detId());
258  int subdet = id.subdet();
259  int zp = id.zside();
260  int type = id.type();
261  int layer = id.layer();
262  int waferU = id.waferU();
263  int waferV = id.waferV();
264 
265  //use first TC to find mipPt conversions to Et and ADC
266  float mipPtToEt_conv = trigCellVecInput[0].et() / trigCellVecInput[0].mipPt();
267  float mipToADC_conv = trigCellVecInput[0].hwPt() / (trigCellVecInput[0].mipPt() * cosh(trigCellVecInput[0].eta()));
268 
269  for (int i = 0; i < nTriggerCells_; i++) {
270  if (ae_outputArray[i] > 0) {
271  int cellU = ae_outputCellU_[i];
272  int cellV = ae_outputCellV_[i];
273 
274  HGCalTriggerDetId id(subdet, zp, type, layer, waferU, waferV, cellU, cellV);
275 
277  continue;
278 
280 
281  double mipPt = ae_outputArray[i] * normalization * renormalizationFactor;
282  double adc = mipPt * cosh(point.eta()) * mipToADC_conv;
283  double et = mipPt * mipPtToEt_conv;
284 
285  if (mipPt < zeroSuppresionThreshold_)
286  continue;
287 
288  l1t::HGCalTriggerCell triggerCell(reco::LeafCandidate::LorentzVector(), adc, 0, 0, 0, id);
289  //Keep the pre-autoencoder charge for this cell
290  triggerCell.setUncompressedCharge(uncompressedCharge[i]);
291  triggerCell.setCompressedCharge(compressedCharge[i]);
292  triggerCell.setMipPt(mipPt);
293 
294  math::PtEtaPhiMLorentzVector p4(et, point.eta(), point.phi(), 0.);
295 
296  triggerCell.setP4(p4);
297  triggerCell.setPosition(point);
298 
299  trigCellVecOutput.push_back(triggerCell);
300  }
301  }
302 
303  if (saveEncodedValues_) {
304  id = HGCalTriggerDetId(subdet, zp, type, layer, waferU, waferV, 0, 0);
305  for (int i = 0; i < nEncodedLayerNodes_; i++) {
306  l1t::HGCalConcentratorData encodedLayerData(ae_encodedLayer_[i] * outputMaxIntSizeGlobal, i, id);
307  ae_encodedLayer_Output.push_back(encodedLayerData);
308  }
309  }
310  }
311 }
bool isScintillator(const DetId &id) const
int32_t waferU(const int32_t index)
GraphDef * loadGraphDef(const std::string &pbFile)
Definition: TensorFlow.cc:129
static constexpr int ae_outputCellU_[nTriggerCells_]
const HGCalTriggerGeometryBase * getTriggerGeometry() const
std::vector< std::unique_ptr< tensorflow::Session > > session_decoder_
std::unique_ptr< tensorflow::GraphDef > graphDef_encoder_
std::vector< std::unique_ptr< tensorflow::Session > > session_encoder_
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
virtual bool validTriggerCell(const unsigned trigger_cell_id) const =0
HGCalConcentratorAutoEncoderImpl(const edm::ParameterSet &conf)
GlobalPoint getTCPosition(const DetId &id) const
static constexpr int cellUVremap_[cellUVSize_][cellUVSize_]
void setMipPt(double value)
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:281
void setCompressedCharge(uint32_t value)
d
Definition: ztail.py:151
Session * createSession()
Definition: TensorFlow.cc:146
std::unique_ptr< tensorflow::GraphDef > graphDef_decoder_
void setLogging(const std::string &level="3")
Definition: TensorFlow.cc:90
void setUncompressedCharge(uint32_t value)
static constexpr unsigned int maxNumberOfLinks_
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
int32_t waferV(const int32_t index)
void select(unsigned nLinks, const std::vector< l1t::HGCalTriggerCell > &trigCellVecInput, std::vector< l1t::HGCalTriggerCell > &trigCellVecOutput, std::vector< l1t::HGCalConcentratorData > &ae_EncodedOutput)
static constexpr int ae_outputCellV_[nTriggerCells_]
void setPosition(const GlobalPoint &position)
void setP4(const LorentzVector &p4) final
set 4-momentum
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
std::vector< edm::ParameterSet > modelFilePaths_
*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
uint16_t *__restrict__ uint16_t const *__restrict__ adc
std::array< double, nEncodedLayerNodes_ > ae_encodedLayer_