CMS 3D CMS Logo

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