CMS 3D CMS Logo

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