CMS 3D CMS Logo

Classes | Namespaces | Typedefs | Functions | Variables
PtAssignmentNNRegression.cc File Reference
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/MuonDetId/interface/CSCDetId.h"
#include "L1Trigger/L1TMuonOverlapPhase2/interface/PtAssignmentNNRegression.h"
#include <boost/archive/text_oarchive.hpp>
#include <boost/archive/text_iarchive.hpp>
#include <boost/algorithm/string/replace.hpp>
#include <sstream>
#include <fstream>

Go to the source code of this file.

Classes

struct  OmtfHit
 

Namespaces

 lutNN
 

Typedefs

typedef LutNetworkFixedPointRegression2Outputs< input_I, input_F, networkInputSize, layer1_lut_I, layer1_lut_F, layer1_neurons, layer1_output_I, layer2_input_I, layer2_lut_I, layer2_lut_F, layer2_neurons, layer3_input_I, layer3_0_inputCnt, layer3_0_lut_I, layer3_0_lut_F, output0_I, output0_F, layer3_1_inputCnt, layer3_1_lut_I, layer3_1_lut_F, output1_I, output1_F > lutNN::LutNetworkFP
 

Functions

bool omtfHitToEventInput (OmtfHit &hit, std::vector< float > &inputs, unsigned int omtfRefLayer, bool print)
 

Variables

static constexpr int lutNN::input_F = 4
 
static constexpr int lutNN::input_I = 10
 
static constexpr int lutNN::layer1_lut_F = 13
 
static constexpr int lutNN::layer1_lut_I = 3
 
static constexpr int lutNN::layer1_neurons = 16
 
static constexpr int lutNN::layer1_output_I = 4
 
static constexpr int lutNN::layer2_input_I = 8
 
static constexpr int lutNN::layer2_lut_F = 11
 
static constexpr int lutNN::layer2_lut_I = 5
 
static constexpr int lutNN::layer2_neurons = 9
 
static constexpr int lutNN::layer3_0_inputCnt = 8
 
static constexpr int lutNN::layer3_0_lut_F = 11
 
static constexpr int lutNN::layer3_0_lut_I = 5
 
static constexpr int lutNN::layer3_1_inputCnt = 1
 
static constexpr int lutNN::layer3_1_lut_F = 11
 
static constexpr int lutNN::layer3_1_lut_I = 4
 
static constexpr int lutNN::layer3_input_I = 5
 
static constexpr std::size_t lutNN::networkInputSize = 18
 
static constexpr int lutNN::output0_F = 2
 
static constexpr int lutNN::output0_I = 8
 
static constexpr int lutNN::output1_F = 0
 
static constexpr int lutNN::output1_I = 8
 

Function Documentation

◆ omtfHitToEventInput()

bool omtfHitToEventInput ( OmtfHit hit,
std::vector< float > &  inputs,
unsigned int  omtfRefLayer,
bool  print 
)

TODO <<<<<<<<<<<<<<<<<<<<<<<<<<<<<

Definition at line 110 of file PtAssignmentNNRegression.cc.

References funct::abs(), TauDecayModes::dec, nano_mu_digi_cff::float, ALPAKA_ACCELERATOR_NAMESPACE::caPixelDoublets::if(), PixelMapPlotter::inputs, createfilelist::int, hltrates_dqm_sourceclient-live_cfg::offset, print(), and hit::z.

Referenced by PtAssignmentNNRegression::getPts().

110  {
111  float offset = (omtfRefLayer << 7) + 64;
112 
113  if (hit.valid) {
114  if ((hit.layer == 1 || hit.layer == 3 || hit.layer == 5) && hit.quality < 4)
115  return false;
116 
117  int rangeFactor = 2; //rangeFactor scales the hit.phiDist such that the event->inputs is smaller then 63
118  if (hit.layer == 1) {
119  rangeFactor = 8;
120  }
121  /*else if(hit.layer == 8 || hit.layer == 17) {
122  rangeFactor = 4;
123  }*/
124  else if (hit.layer == 3) {
125  rangeFactor = 4;
126  } else if (hit.layer == 9) {
127  rangeFactor = 1;
128  }
129  /*else {
130  rangeFactor = 2;
131  }
132  */
133 
134  rangeFactor *= 2; //TODO !!!!!!!!!!!!!!!!!!!
135 
136  if (std::abs(hit.phiDist) >= (63 * rangeFactor)) {
137  edm::LogImportant("OMTFReconstruction") //<<" muonPt "<<omtfEvent.muonPt<<" omtfPt "<<omtfEvent.omtfPt
138  << " RefLayer " << omtfRefLayer << " layer " << int(hit.layer) << " hit.phiDist " << hit.phiDist << " valid "
139  << ((short)hit.valid) << " !!!!!!!!!!!!!!!!!!!!!!!!" << endl;
140  hit.phiDist = copysign(63 * rangeFactor, hit.phiDist);
141  }
142 
143  inputs.at(hit.layer) = (float)hit.phiDist / (float)rangeFactor + offset;
144 
145  if (inputs.at(hit.layer) >=
146  1022) //the last address i.e. 1023 is reserved for the no-hit value, so interpolation between the 1022 and 1023 has no sense
147  inputs.at(hit.layer) = 1022;
148 
149  if (print || inputs.at(hit.layer) < 0) {
150  edm::LogImportant("OMTFReconstruction") //<<"rawData "<<hex<<setw(16)<<hit.rawData
151  << " layer " << dec << int(hit.layer);
152  edm::LogImportant("OMTFReconstruction")
153  << " phiDist " << hit.phiDist << " inputVal " << inputs.at(hit.layer) << " hit.z " << int(hit.z) << " valid "
154  << ((short)hit.valid) << " quality " << (short)hit.quality << " omtfRefLayer " << omtfRefLayer;
155  if (inputs.at(hit.layer) < 0)
156  edm::LogImportant("OMTFReconstruction") << " event->inputs.at(hit.layer) < 0 !!!!!!!!!!!!!!!!!" << endl;
157  edm::LogImportant("OMTFReconstruction") << endl;
158  }
159 
160  if (inputs[hit.layer] >= 1024) { //TODO should be the size of the LUT of the first layer
161  edm::LogImportant("OMTFReconstruction") << " event->inputs[hit.layer] >= 1024 !!!!!!!!!!!!!!!!!" << endl;
162  }
163  return true;
164  }
165 
166  return false;
167 }
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Log< level::Error, true > LogImportant
if(threadIdxLocalY==0 &&threadIdxLocalX==0)