CMS 3D CMS Logo

HcalTriggerPrimitiveAlgo.h
Go to the documentation of this file.
1 #ifndef HcalSimAlgos_HcalTriggerPrimitiveAlgo_h
2 #define HcalSimAlgos_HcalTriggerPrimitiveAlgo_h
3 
6 
11 
13 
15 
18 
19 #include <map>
20 #include <vector>
21 
22 class CaloGeometry;
23 class IntegerCaloSamples;
24 
26 public:
27  HcalTriggerPrimitiveAlgo(bool pf, const std::vector<double>& w, int latency,
28  uint32_t FG_threshold, uint32_t FG_HF_threshold, uint32_t ZS_threshold,
31  uint32_t minSignalThreshold=0, uint32_t PMT_NoiseThreshold=0);
33 
34  template<typename... Digis>
35  void run(const HcalTPGCoder* incoder,
36  const HcalTPGCompressor* outcoder,
37  const HcalDbService* conditions,
39  const HcalTrigTowerGeometry* trigTowerGeometry,
40  float rctlsb, const HcalFeatureBit* LongvrsShortCut,
41  const Digis&... digis);
42 
43  template<typename T, typename... Args>
44  void addDigis(const T& collection, const Args&... digis) {
45  addDigis(collection);
46  addDigis(digis...);
47  };
48 
49  template<typename T>
50  void addDigis(const T& collection) {
51  for (const auto& digi: collection) {
52  addSignal(digi);
53  }
54  };
55 
56  template<typename D>
58  for (auto i = collection.begin(); i != collection.end(); ++i) {
59  D digi(*i);
60  addSignal(digi);
61  }
62  };
63 
65  void runFEFormatError(const FEDRawDataCollection* rawraw,
66  const HcalElectronicsMap* emap,
68  void setPeakFinderAlgorithm(int algo);
69  void setNCTScaleShift(int);
70  void setRCTScaleShift(int);
71 
72  void setUpgradeFlags(bool hb, bool he, bool hf);
73  void overrideParameters(const edm::ParameterSet& ps);
74 
75  private:
76 
78  void addSignal(const HBHEDataFrame & frame);
79  void addSignal(const HFDataFrame & frame);
80  void addSignal(const QIE10DataFrame& frame);
81  void addSignal(const QIE11DataFrame& frame);
82  void addSignal(const IntegerCaloSamples & samples);
83  void addFG(const HcalTrigTowerDetId& id, std::vector<bool>& msb);
84  void addUpgradeFG(const HcalTrigTowerDetId& id, int depth, const std::vector<std::bitset<2>>& bits);
85 
86  bool validUpgradeFG(const HcalTrigTowerDetId& id, int depth) const;
87  bool validChannel(const QIE10DataFrame& digi, int ts) const;
88  bool needLegacyFG(const HcalTrigTowerDetId& id) const;
89 
91  void analyze(IntegerCaloSamples & samples, HcalTriggerPrimitiveDigi & result);
92  // 2017: QIE11
93  void analyze2017(IntegerCaloSamples& samples, HcalTriggerPrimitiveDigi& result, const HcalFinegrainBit& fg_algo);
94  // Version 0: RCT
95  void analyzeHF(IntegerCaloSamples & samples, HcalTriggerPrimitiveDigi & result, const int hf_lumi_shift);
96  // Version 1: 1x1
97  void analyzeHF2016(
98  const IntegerCaloSamples& SAMPLES,
100  const int HF_LUMI_SHIFT,
101  const HcalFeatureBit* HCALFEM
102  );
103  // With dual anode readout
104  void analyzeHF2017(
105  const IntegerCaloSamples& SAMPLES,
106  HcalTriggerPrimitiveDigi& result,
107  const int HF_LUMI_SHIFT,
108  const HcalFeatureBit* HCALFEM
109  );
110 
111  // Member initialized by constructor
115  double theThreshold;
116  bool peakfind_;
117  std::vector<double> weights_;
118  int latency_;
119  uint32_t FG_threshold_;
121  uint32_t ZS_threshold_;
131 
132  // Algo1: isPeak = TS[i-1] < TS[i] && TS[i] >= TS[i+1]
133  // Algo2: isPeak = TSS[i-1] < TSS[i] && TSS[i] >= TSS[i+1],
134  // TSS[i] = TS[i] + TS[i+1]
135  // Default: Algo2
137 
138  // Member not initialzed
139  //std::vector<HcalTrigTowerDetId> towerIds(const HcalDetId & id) const;
140 
142 
143  typedef std::map<HcalTrigTowerDetId, IntegerCaloSamples> SumMap;
144  SumMap theSumMap;
145 
146  struct HFDetails {
151  };
152  typedef std::map<HcalTrigTowerDetId, std::map<uint32_t, HFDetails>> HFDetailMap;
153  HFDetailMap theHFDetailMap;
154 
158  std::vector<bool> validity;
159  std::vector<bool> fgbit;
160  };
161  typedef std::map<HcalTrigTowerDetId, std::map<uint32_t, std::array<HFUpgradeDetails, 4>>> HFUpgradeDetailMap;
162  HFUpgradeDetailMap theHFUpgradeDetailMap;
163 
164  typedef std::vector<IntegerCaloSamples> SumFGContainer;
165  typedef std::map< HcalTrigTowerDetId, SumFGContainer > TowerMapFGSum;
166  TowerMapFGSum theTowerMapFGSum;
167 
168  // ==============================
169  // = HF Veto
170  // ==============================
171  // Sum = Long + Short;" // intermediate calculation.
172  // if ((Short < MinSignalThresholdET OR Long < MinSignalThresholdET)
173  // AND Sum > PMTNoiseThresholdET) VetoedSum = 0;
174  // else VetoedSum = Sum;
175  // ==============================
176  // Map from FG id to veto booleans
178  typedef std::map<uint32_t, std::vector<bool> > TowerMapVeto;
179  TowerMapVeto HF_Veto;
180 
181  typedef std::map<HcalTrigTowerDetId, std::vector<bool> > FGbitMap;
182  FGbitMap fgMap_;
183 
184  typedef std::vector<HcalFinegrainBit::Tower> FGUpgradeContainer;
185  typedef std::map<HcalTrigTowerDetId, FGUpgradeContainer> FGUpgradeMap;
186  FGUpgradeMap fgUpgradeMap_;
187 
188  bool upgrade_hb_ = false;
189  bool upgrade_he_ = false;
190  bool upgrade_hf_ = false;
191 
193 
194  bool override_adc_hf_ = false;
196  bool override_tdc_hf_ = false;
197  unsigned long long override_tdc_hf_value_;
198 
199  // HE constants
200  static const int HBHE_OVERLAP_TOWER = 16;
201  static const int LAST_FINEGRAIN_DEPTH = 6;
202  static const int LAST_FINEGRAIN_TOWER = 28;
203 
204  // Fine-grain in HF ignores tower 29, and starts with 30
205  static const int FIRST_FINEGRAIN_TOWER = 30;
206 
210  // Consider CaloTPGTranscoderULUT.h for values
211  static const int QIE10_MAX_LINEARIZATION_ET = 0x7FF;
212  static const int QIE11_MAX_LINEARIZATION_ET = 0x7FF;
213 };
214 
215 template<typename... Digis>
217  const HcalTPGCompressor* outcoder,
218  const HcalDbService* conditions,
220  const HcalTrigTowerGeometry* trigTowerGeometry,
221  float rctlsb, const HcalFeatureBit* LongvrsShortCut,
222  const Digis&... digis) {
223  theTrigTowerGeometry = trigTowerGeometry;
224 
225  incoder_ = dynamic_cast<const HcaluLUTTPGCoder*>(incoder);
226  outcoder_ = outcoder;
227  conditions_ = conditions;
228 
229  theSumMap.clear();
230  theTowerMapFGSum.clear();
231  HF_Veto.clear();
232  fgMap_.clear();
233  fgUpgradeMap_.clear();
234  theHFDetailMap.clear();
235  theHFUpgradeDetailMap.clear();
236 
237  // Add all digi collections
238  addDigis(digis...);
239 
240  // Prepare the fine-grain calculation algorithm for HB/HE
241  int version = 0;
244  if (override_parameters_.exists("FGVersionHBHE"))
245  version = override_parameters_.getParameter<uint32_t>("FGVersionHBHE");
246  HcalFinegrainBit fg_algo(version);
247 
248  // VME produces additional bits on the front used by lumi but not the
249  // trigger, this shift corrects those out by right shifting over them.
250  for (auto& item: theSumMap) {
251  result.push_back(HcalTriggerPrimitiveDigi(item.first));
252  HcalTrigTowerDetId detId(item.second.id());
253  if(detId.ietaAbs() >= theTrigTowerGeometry->firstHFTower(detId.version())) {
254  if (detId.version() == 0) {
255  analyzeHF(item.second, result.back(), RCTScaleShift);
256  } else if (detId.version() == 1) {
257  if (upgrade_hf_)
258  analyzeHF2017(item.second, result.back(), NCTScaleShift, LongvrsShortCut);
259  else
260  analyzeHF2016(item.second, result.back(), NCTScaleShift, LongvrsShortCut);
261  } else {
262  // Things are going to go poorly
263  }
264  }
265  else {
266  // Determine which energy reconstruction path to take based on the
267  // fine-grain availability:
268  // * QIE8 TP add entries into fgMap_
269  // * QIE11 TP add entries into fgUpgradeMap_
270  // (not for tower 16 unless HB is upgraded, too)
271  if (fgMap_.find(item.first) != fgMap_.end()) {
272  analyze(item.second, result.back());
273  } else if (fgUpgradeMap_.find(item.first) != fgUpgradeMap_.end()) {
274  analyze2017(item.second, result.back(), fg_algo);
275  }
276  }
277  }
278 
279  // Free up some memory
280  theSumMap.clear();
281  theTowerMapFGSum.clear();
282  HF_Veto.clear();
283  fgMap_.clear();
284  fgUpgradeMap_.clear();
285  theHFDetailMap.clear();
286  theHFUpgradeDetailMap.clear();
287 
288  return;
289 }
290 
291 #endif
T getParameter(std::string const &) const
void analyze(IntegerCaloSamples &samples, HcalTriggerPrimitiveDigi &result)
adds the actual RecHits
void runFEFormatError(const FEDRawDataCollection *rawraw, const HcalElectronicsMap *emap, HcalTrigPrimDigiCollection &result)
std::map< uint32_t, std::vector< bool > > TowerMapVeto
HcalTriggerPrimitiveAlgo(bool pf, const std::vector< double > &w, int latency, uint32_t FG_threshold, uint32_t FG_HF_threshold, uint32_t ZS_threshold, int numberOfSamples, int numberOfPresamples, int numberOfSamplesHF, int numberOfPresamplesHF, uint32_t minSignalThreshold=0, uint32_t PMT_NoiseThreshold=0)
HFUpgradeDetailMap theHFUpgradeDetailMap
void analyzeHF(IntegerCaloSamples &samples, HcalTriggerPrimitiveDigi &result, const int hf_lumi_shift)
const double w
Definition: UKUtility.cc:23
static const int QIE10_MAX_LINEARIZATION_ET
void analyzeHF2016(const IntegerCaloSamples &SAMPLES, HcalTriggerPrimitiveDigi &result, const int HF_LUMI_SHIFT, const HcalFeatureBit *HCALFEM)
static const int QIE11_LUT_BITMASK
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision bits
bool exists(std::string const &parameterName) const
checks if a parameter exists
void push_back(T const &t)
const HcalTPGCompressor * outcoder_
std::map< HcalTrigTowerDetId, FGUpgradeContainer > FGUpgradeMap
void addFG(const HcalTrigTowerDetId &id, std::vector< bool > &msb)
const_iterator begin() const
static const int QIE8_LUT_BITMASK
bool needLegacyFG(const HcalTrigTowerDetId &id) const
numberOfSamples
threshold for setting fine grain bit
bool validChannel(const QIE10DataFrame &digi, int ts) const
std::vector< IntegerCaloSamples > SumFGContainer
const HcalTrigTowerGeometry * theTrigTowerGeometry
const HcalDbService * conditions_
void run(const HcalTPGCoder *incoder, const HcalTPGCompressor *outcoder, const HcalDbService *conditions, HcalTrigPrimDigiCollection &result, const HcalTrigTowerGeometry *trigTowerGeometry, float rctlsb, const HcalFeatureBit *LongvrsShortCut, const Digis &...digis)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
int getFGVersionHBHE() const
get FineGrain Algorithm Version for HBHE
void runZS(HcalTrigPrimDigiCollection &tp)
static const int QIE10_LUT_BITMASK
std::map< HcalTrigTowerDetId, std::vector< bool > > FGbitMap
std::vector< HcalFinegrainBit::Tower > FGUpgradeContainer
void addDigis(const T &collection, const Args &...digis)
void addDigis(const T &collection)
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:151
std::map< HcalTrigTowerDetId, SumFGContainer > TowerMapFGSum
unsigned long long override_tdc_hf_value_
bool validUpgradeFG(const HcalTrigTowerDetId &id, int depth) const
static const int QIE11_MAX_LINEARIZATION_ET
std::map< HcalTrigTowerDetId, IntegerCaloSamples > SumMap
const HcalTPParameters * getHcalTPParameters() const
ZS_threshold
threshold for setting fine grain bit
void addUpgradeFG(const HcalTrigTowerDetId &id, int depth, const std::vector< std::bitset< 2 >> &bits)
const_iterator end() const
std::map< HcalTrigTowerDetId, std::map< uint32_t, std::array< HFUpgradeDetails, 4 > > > HFUpgradeDetailMap
const HcaluLUTTPGCoder * incoder_
void analyzeHF2017(const IntegerCaloSamples &SAMPLES, HcalTriggerPrimitiveDigi &result, const int HF_LUMI_SHIFT, const HcalFeatureBit *HCALFEM)
long double T
latency
hardware algo
void addSignal(const HBHEDataFrame &frame)
adds the signal to the map
std::map< HcalTrigTowerDetId, std::map< uint32_t, HFDetails > > HFDetailMap
void setUpgradeFlags(bool hb, bool he, bool hf)
void addDigis(const HcalDataFrameContainer< D > &collection)
void overrideParameters(const edm::ParameterSet &ps)
const_reference back() const
void analyze2017(IntegerCaloSamples &samples, HcalTriggerPrimitiveDigi &result, const HcalFinegrainBit &fg_algo)
int firstHFTower(int version) const