CMS 3D CMS Logo

L1TZDCProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1Trigger/L1TZDC
4 // Class: L1TZDC
5 //
14 //
15 // Original Author: James Brooke
16 // Created: Thu, 05 Dec 2013 17:39:27 GMT
17 //
18 // Copied for ZDC by: Chris McGinn
19 // Copy Made: Wed, 03 Aug 2023
20 // Contact: christopher.mc.ginn@cern.ch or
21 // cfmcginn on github for bugs/issues
22 //
23 #include <memory>
24 #include <iostream>
25 // user include files
26 
39 
43 
47 
53 
55 
56 //
57 // class declaration
58 //
59 
60 using namespace l1t;
61 
63 public:
64  explicit L1TZDCProducer(const edm::ParameterSet& ps);
65  ~L1TZDCProducer() override = default;
66 
67  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
68 
69 private:
70  void produce(edm::Event&, const edm::EventSetup&) override;
71 
72  void beginRun(edm::Run const&, edm::EventSetup const&) override;
73  // void endRun(edm::Run const&, edm::EventSetup const&) override;
74 
75  int zdcLUTIndexHelper(int iDetPos, int iBXPos);
76 
77  // ----------member data ---------------------------
78 
79  // input tokens
80  //Add the ZDC token, candidateToken for caloParamsHelper
83 
84  //input ints
85  int bxFirst_;
86  int bxLast_;
88 
89  // put tokens
91 
92  //Following the L1TStage2Layer2Producer
93  std::unique_ptr<CaloParamsHelper> params_;
95 };
96 
98  // register what you produce
99  etToken_ = produces<EtSumBxCollection>();
100 
101  // register what you consume and keep token for later access:
102  zdcToken_ = consumes<QIE10DigiCollection>(ps.getParameter<edm::InputTag>("zdcDigis"));
103  candidateToken_ = esConsumes<CaloParams, L1TCaloParamsRcd, edm::Transition::BeginRun>();
104 
105  bxFirst_ = ps.getParameter<int>("bxFirst");
106  bxLast_ = ps.getParameter<int>("bxLast");
107  sampleToCenterBX_ = ps.getParameter<int>("sampleToCenterBX");
108 }
109 
110 // ------------ method called to produce the data ------------
112  using namespace edm;
113  using namespace l1t;
114 
115  LogDebug("L1TZDCProducer") << "L1TZDCProducer::produce function called..." << std::endl;
116 
117  // reduced collection to be emplaced in output
118  EtSumBxCollection etsumsReduced(0, bxFirst_, bxLast_);
119 
120  //inputs
121  Handle<QIE10DigiCollection> zdcDigiCollection;
122  iEvent.getByToken(zdcToken_, zdcDigiCollection);
123 
124  //Produce ZDC EtSums IF input zdcDigiCollection is valid
125  if (zdcDigiCollection.isValid()) {
126  //In lieu of bxFirst, bxLast, use the number of the timeslice samples
127  const QIE10DataFrame& frametest = (*zdcDigiCollection)[0];
128  int nSamples = frametest.samples();
129 
130  //Check if the sample requested as center bx is valid
131  int sampleToCenterBXChecked = sampleToCenterBX_;
132 
133  if (sampleToCenterBXChecked < 0) {
134  edm::LogWarning("L1TZDCProducer") << "sampleToCenterBX LT 0; Set bx 0 to sample 0 (minimum allowed)" << std::endl;
135  sampleToCenterBXChecked = 0;
136  } else if (sampleToCenterBXChecked >= nSamples) {
137  edm::LogWarning("L1TZDCProducer")
138  << "sampleToCenterBX GE nsamples; Set bx 0 to sample nsamples-1 (maximum allowed)" << std::endl;
139  sampleToCenterBXChecked = nSamples - 1;
140  }
141 
142  //rawadc[detector index][time slices]
143  unsigned short rawadc[18][10];
144 
145  // the loop below loops over all the elements of the QIE10DigiCollection. Each entry corresponds to one channel
146  for (QIE10DigiCollection::const_iterator it = zdcDigiCollection->begin(); it != zdcDigiCollection->end(); it++) {
147  const QIE10DataFrame& frame(*it);
148  HcalZDCDetId cell = frame.id();
149  int zside = cell.zside();
150  int section = cell.section();
151  int channel = cell.channel();
152 
153  if (zside != -1 && zside != 1)
154  continue;
155  if (section != 1 && section != 2)
156  continue;
157  if (section == 1 && (channel < 1 || channel > 5))
158  continue;
159  if (section == 2 && (channel < 1 || channel > 4))
160  continue;
161 
162  int ihitid = (zside == 1 ? 9 : 0) + (section == 2 ? 5 : 0) + (channel - 1);
163  //the loop below iterates over the time slices
164  for (int iTS = 0; iTS < nSamples; iTS++) {
165  unsigned short adc = (unsigned short)frame[iTS].adc();
166  rawadc[ihitid][iTS] = adc;
167  } // end of loop over iTS
168  } //end of loop over channels
169 
170  for (int ibx = 0; ibx < nSamples; ibx++) {
171  double cEMP = 0, cEMM = 0, cHDP = 0, cHDM = 0;
172  double sumcEMP = 0, sumcEMM = 0, sumcHDP = 0, sumcHDM = 0;
173  //idet=0-4 correpond to the EM channels
174  for (int idet = 0; idet < 5; idet++) {
175  unsigned short EMP = rawadc[idet + 9][ibx];
176  unsigned short EMM = rawadc[idet][ibx];
177 
178  int cEMP_LUTIndex = zdcLUTIndexHelper(idet + 9, (int)EMP);
179  int cEMM_LUTIndex = zdcLUTIndexHelper(idet, (int)EMM);
180 
181  cEMP = ((double)params_->zdcLUT()->data(cEMP_LUTIndex)) / ((double)scaleFactor_);
182  cEMM = ((double)params_->zdcLUT()->data(cEMM_LUTIndex)) / ((double)scaleFactor_);
183 
184  sumcEMP = sumcEMP + cEMP;
185  sumcEMM = sumcEMM + cEMM;
186  }
187  //idet=5-8 correspond to HAD channels
188  for (int idet = 5; idet < 9; idet++) {
189  unsigned short HDP = rawadc[idet + 9][ibx];
190  unsigned short HDM = rawadc[idet][ibx];
191 
192  int cHDP_LUTIndex = zdcLUTIndexHelper(idet + 9, (int)HDP);
193  int cHDM_LUTIndex = zdcLUTIndexHelper(idet, (int)HDM);
194 
195  cHDP = ((double)params_->zdcLUT()->data(cHDP_LUTIndex)) / ((double)scaleFactor_);
196  cHDM = ((double)params_->zdcLUT()->data(cHDM_LUTIndex)) / ((double)scaleFactor_);
197 
198  sumcHDP = sumcHDP + cHDP;
199  sumcHDM = sumcHDM + cHDM;
200  }
201  double sumM = sumcEMM + sumcHDM;
202  double sumP = sumcEMP + sumcHDP;
203 
204  if (ibx == 4) {
205  edm::LogInfo("L1TZDCProducer") << ", sumM= " << sumM << std::endl;
206  edm::LogInfo("L1TZDCProducer") << ", sumP= " << sumP << std::endl;
207  }
208  l1t::EtSum tempEtM = l1t::EtSum();
209  tempEtM.setHwPt(sumM);
210  tempEtM.setHwEta(-1.);
211  tempEtM.setHwPhi(0.);
213 
214  l1t::EtSum tempEtP = l1t::EtSum();
215  tempEtP.setHwPt(sumP);
216  tempEtP.setHwEta(1.);
217  tempEtP.setHwPhi(0.);
219 
220  if (ibx >= sampleToCenterBXChecked + bxFirst_ && ibx <= sampleToCenterBXChecked + bxLast_) {
221  etsumsReduced.push_back(ibx - sampleToCenterBXChecked, CaloTools::etSumP4Demux(tempEtP));
222  etsumsReduced.push_back(ibx - sampleToCenterBXChecked, CaloTools::etSumP4Demux(tempEtM));
223  }
224  } // end of loop over bunch crossings
225  } // end if(zdcDigiCollection.isValid())
226  else {
227  // If the collection is not valid issue a warning before putting an empty collection
228  edm::LogWarning("L1TZDCProducer") << "zdcDigis not valid; return empty ZDC Et Sum BXCollection" << std::endl;
229  }
230 
231  // Emplace even if !zdcDigiCollection.isValid()
232  // Output in this case will be an empty collection
233  iEvent.emplace(etToken_, std::move(etsumsReduced));
234 }
235 
236 // ------------ method called when starting to processes a run ------------
237 void L1TZDCProducer::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
238  edm::ESHandle<CaloParams> candidateHandle = iSetup.getHandle(candidateToken_);
239 
240  params_ = std::make_unique<CaloParamsHelper>(*candidateHandle.product());
241  scaleFactor_ = params_->zdcLUT()->data(0); //First position is the integer scaling factor
242  edm::LogInfo("L1TZDCProducer") << "SCALE FACTOR FOR LUT: " << scaleFactor_ << std::endl;
243 }
244 
245 // LUT HELPER METHOD
246 int L1TZDCProducer::zdcLUTIndexHelper(int iDetPos, int iBxPos) { return 1 + iDetPos * 256 + iBxPos; }
247 
248 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
251  desc.add<edm::InputTag>("zdcDigis", edm::InputTag("hcalDigis", "ZDC"));
252  desc.add<int>("bxFirst", -2);
253  desc.add<int>("bxLast", 2);
254  desc.add<int>("sampleToCenterBX", 4);
255  descriptions.add("l1tZDCProducer", desc);
256 }
257 
258 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
Section section() const
get the section
Definition: HcalZDCDetId.cc:44
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &, const edm::EventSetup &) override
delete x;
Definition: CaloConfig.h:22
int zside(DetId const &)
int zdcLUTIndexHelper(int iDetPos, int iBXPos)
void setType(EtSumType type)
Definition: EtSum.cc:11
edm::ESGetToken< CaloParams, L1TCaloParamsRcd > candidateToken_
edm::EDPutTokenT< EtSumBxCollection > etToken_
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
L1TZDCProducer(const edm::ParameterSet &ps)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
const_iterator end() const
Log< level::Info, false > LogInfo
std::unique_ptr< CaloParamsHelper > params_
void setHwPhi(int phi)
Definition: L1Candidate.h:30
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const_iterator begin() const
The iterator returned can not safely be used across threads.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< QIE10DigiCollection > zdcToken_
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
void setHwPt(int pt)
Definition: L1Candidate.h:28
void setHwEta(int eta)
Definition: L1Candidate.h:29
void beginRun(edm::Run const &, edm::EventSetup const &) override
Log< level::Warning, false > LogWarning
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalZDCDetId.h:39
static l1t::EtSum etSumP4Demux(l1t::EtSum &)
Definition: CaloTools.cc:364
def move(src, dest)
Definition: eostools.py:511
constexpr int samples() const
total number of samples in the digi
int channel() const
get the channel
Definition: HcalZDCDetId.cc:63
Definition: Run.h:45
uint16_t *__restrict__ uint16_t const *__restrict__ adc
void push_back(int bx, T object)
#define LogDebug(id)