CMS 3D CMS Logo

HcalTTPDigiProducer.cc
Go to the documentation of this file.
15 
16 #include <cstdio>
17 
19 public:
20  explicit HcalTTPDigiProducer(const edm::ParameterSet& ps);
21  ~HcalTTPDigiProducer() override = default;
22 
23  void produce(edm::Event& e, const edm::EventSetup& c) override;
24 
25 private:
26  bool isMasked(HcalDetId id);
27  bool decision(int nP, int nM, int bit);
28 
31  std::vector<unsigned int> maskedChannels_;
33  int calc_[4];
34  int nHits_[4], nHFp_[4], nHFm_[4];
35  char pReq_[4], mReq_[4], pmLogic_[4];
37  int fwAlgo_;
39  unsigned int threshold_;
40 
41  int SoI_;
42 
43  static const int inputs_[];
44 };
45 
46 // DO NOT MODIFY: Mapping between iphi (array index) and TTP input (value) for HF
47 const int HcalTTPDigiProducer::inputs_[] = {30, 66, 4, 44, 4, 44, 0, 68, 0, 68, 16, 48, 16, 48, 6, 46, 6, 46,
48  2, 70, 2, 70, 18, 50, 18, 50, 12, 40, 12, 40, 8, 52, 8, 52, 20, 36,
49  20, 36, 14, 42, 14, 42, 10, 54, 10, 54, 22, 38, 22, 38, 24, 56, 24, 56,
50  32, 60, 32, 60, 28, 64, 28, 64, 26, 58, 26, 58, 34, 62, 34, 62, 30, 66};
51 
53  tok_hf_ = consumes<HFDigiCollection>(ps.getParameter<edm::InputTag>("HFDigiCollection"));
54  tok_tpgCoder_ = esConsumes<HcalTPGCoder, HcalTPGRecord>();
55  maskedChannels_ = ps.getParameter<std::vector<unsigned int> >("maskedChannels");
56  bit_[0] = ps.getParameter<std::string>("defTT8");
57  bit_[1] = ps.getParameter<std::string>("defTT9");
58  bit_[2] = ps.getParameter<std::string>("defTT10");
59  bit_[3] = ps.getParameter<std::string>("defTTLocal");
60 
61  for (int i = 0; i < 4; i++) {
62  nHits_[i] = -1;
63  nHFp_[i] = -1;
64  nHFm_[i] = -1;
65  pReq_[i] = ' ';
66  mReq_[i] = ' ';
67  pmLogic_[i] = ' ';
68  calc_[i] = sscanf(bit_[i].c_str(),
69  "hits>=%d:hfp%c=%d%chfm%c=%d",
70  &(nHits_[i]),
71  &(pReq_[i]),
72  &(nHFp_[i]),
73  &(pmLogic_[i]),
74  &(mReq_[i]),
75  &(nHFm_[i]));
76  if (calc_[i] == 1) {
77  if (nHits_[i] < 0)
78  throw cms::Exception("HcalTTPDigiProducer") << "Unable to read logic for technical trigger";
79  } else if (calc_[i] == 6) {
80  if (nHits_[i] < 0 || nHFp_[i] < 0 || nHFm_[i] < 0)
81  throw cms::Exception("HcalTTPDigiProducer") << "Unable to read logic for technical trigger";
82  if ((pReq_[i] != '>' && pReq_[i] != '<') || (mReq_[i] != '>' && mReq_[i] != '<') ||
83  (pmLogic_[i] != ':' && pmLogic_[i] != '|'))
84  throw cms::Exception("HcalTTPDigiProducer") << "Technical Trigger logic must obey the following format:\n"
85  "\"hits>=[A1]:hfp[B1]=[A2][C]hfm[B2]=[A3]\",\n"
86  "or \"hits>=[A1]\",\n"
87  "with A# >= 0, B# = (</>) and C = (:/|)";
88  } else {
89  throw cms::Exception("HcalTTPDigiProducer") << "Unable to read logic for technical trigger";
90  }
91  }
92 
93  id_ = ps.getUntrackedParameter<int>("id", -1);
94  samples_ = ps.getParameter<int>("samples");
95  presamples_ = ps.getParameter<int>("presamples");
96  iEtaMin_ = ps.getParameter<int>("iEtaMin");
97  iEtaMax_ = ps.getParameter<int>("iEtaMax");
98  threshold_ = ps.getParameter<unsigned int>("threshold");
99  fwAlgo_ = ps.getParameter<int>("fwAlgorithm");
100 
101  SoI_ = ps.getParameter<int>("HFSoI");
102 
103  if (samples_ > 8) {
104  samples_ = 8;
105  edm::LogWarning("HcalTTPDigiProducer") << "Samples forced to maximum value of 8";
106  }
107  if (presamples_ - SoI_ > 0) { // Too many presamples
108  presamples_ = SoI_;
109  edm::LogWarning("HcalTTPDigiProducer") << "Presamples reset to HF SoI value";
110  }
111 
112  produces<HcalTTPDigiCollection>();
113 }
114 
116  for (unsigned int i = 0; i < maskedChannels_.size(); i++)
117  if (id.rawId() == maskedChannels_.at(i))
118  return true;
119  return false;
120 }
121 
122 bool HcalTTPDigiProducer::decision(int nP, int nM, int bit) {
123  bool pOK = false;
124  bool mOK = false;
125  if ((nP + nM) < nHits_[bit])
126  return false;
127  if (calc_[bit] == 1)
128  return ((nP + nM) >= nHits_[bit]);
129 
130  if (pReq_[bit] == '>')
131  pOK = (nP >= nHFp_[bit]);
132  else if (pReq_[bit] == '<')
133  pOK = (nP <= nHFp_[bit]);
134 
135  if (mReq_[bit] == '>')
136  mOK = (nM >= nHFm_[bit]);
137  else if (mReq_[bit] == '<')
138  mOK = (nM <= nHFm_[bit]);
139 
140  if (pmLogic_[bit] == ':')
141  return (pOK && mOK);
142  else if (pmLogic_[bit] == '|')
143  return (pOK || mOK);
144 
145  // Should not ever get here...need to create a warning message
146  edm::LogWarning("HcalTTPDigiProducer") << "Trigger logic exhausted. Returning false";
147  return false;
148 }
149 
151  // Step A: Get Inputs
152  edm::Handle<HFDigiCollection> hfDigiCollection;
153  e.getByToken(tok_hf_, hfDigiCollection);
154  edm::ESHandle<HcalTPGCoder> inputCoder = eventSetup.getHandle(tok_tpgCoder_);
155 
156  // Step B: Create empty output
157  std::unique_ptr<HcalTTPDigiCollection> ttpResult(new HcalTTPDigiCollection());
158 
159  // Step C: Compute TTP inputs
160  uint16_t trigInputs[40];
161  int nP[8];
162  int nM[8];
163  for (int i = 0; i < 8; i++) {
164  nP[i] = 0;
165  nM[i] = 0;
166  for (int j = 0; j < 5; j++)
167  trigInputs[j * 8 + i] = 0;
168  }
169  for (HFDigiCollection::const_iterator theDigi = hfDigiCollection->begin(); theDigi != hfDigiCollection->end();
170  theDigi++) {
171  HcalDetId id = HcalDetId(theDigi->id());
172  if (isMasked(id))
173  continue;
174  if (id.ietaAbs() < iEtaMin_ || id.ietaAbs() > iEtaMax_)
175  continue;
176 
177  IntegerCaloSamples samples(id, theDigi->size());
178  inputCoder->adc2Linear(*theDigi, samples);
179 
180  for (int relSample = -presamples_; relSample < (samples_ - presamples_); relSample++) {
181  if (samples[SoI_ + relSample] >= threshold_) {
182  int linSample = presamples_ + relSample;
183  int offset = (-1 + id.zside()) / 2;
184  int shift = inputs_[id.iphi() + offset];
185  int group = 0;
186  while (shift >= 16) {
187  shift -= 16;
188  group++;
189  }
190  if (!(trigInputs[(linSample * 8) + group] & (1 << shift)))
191  (id.ieta() > 0) ? (nP[linSample]++) : (nM[linSample]++);
192  trigInputs[(linSample * 8) + group] |= (1 << shift);
193  }
194  }
195  }
196 
197  // Step D: Compute trigger decision and fill TTP digi
198  uint8_t trigOutput[8];
199  uint32_t algoDepBits[8];
200  HcalTTPDigi ttpDigi(id_, samples_, presamples_, 0, fwAlgo_, 0);
201  for (int linSample = 0; linSample < 8; linSample++) {
202  trigOutput[linSample] = 0;
203  algoDepBits[linSample] = 0;
204  if (linSample < samples_) {
205  for (int j = 0; j < 4; j++)
206  trigOutput[linSample] |= (decision(nP[linSample], nM[linSample], j) << j);
207  int nT = nP[linSample] + nM[linSample];
208 
209  // Algorithm Dependent bits for FW flavor = 1
210  // NOTE: this disagrees with the fw var. names that implies (LSB) T,M,P (MSB)
211  if (fwAlgo_ == 1)
212  algoDepBits[linSample] = (nT & 0x7F) | ((nP[linSample] & 0x3F) << 7) | ((nM[linSample] & 0x3F) << 13);
213  ttpDigi.setSample(
214  (linSample - presamples_), &trigInputs[linSample * 8], algoDepBits[linSample], trigOutput[linSample]);
215  }
216  }
217  ttpResult->push_back(ttpDigi);
218 
219  // Step E: Put outputs into event
220  e.put(std::move(ttpResult));
221 }
222 
225 
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ESGetToken< HcalTPGCoder, HcalTPGRecord > tok_tpgCoder_
std::vector< T >::const_iterator const_iterator
void produce(edm::Event &e, const edm::EventSetup &c) override
std::vector< unsigned int > maskedChannels_
bool isMasked(HcalDetId id)
T getUntrackedParameter(std::string const &, T const &) const
static const int inputs_[]
~HcalTTPDigiProducer() override=default
void setSample(int relativeSample, const uint16_t *triggerInputs, const uint32_t algodep, const uint8_t outputTrigger)
Definition: HcalTTPDigi.cc:28
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const_iterator begin() const
bool decision(int nP, int nM, int bit)
const_iterator end() const
virtual void adc2Linear(const HBHEDataFrame &df, IntegerCaloSamples &ics) const =0
edm::SortedCollection< HcalTTPDigi > HcalTTPDigiCollection
edm::EDGetTokenT< HFDigiCollection > tok_hf_
static unsigned int const shift
Log< level::Warning, false > LogWarning
HcalTTPDigiProducer(const edm::ParameterSet &ps)
def move(src, dest)
Definition: eostools.py:511