CMS 3D CMS Logo

ScCALORawToDigi.cc
Go to the documentation of this file.
2 
4  using namespace edm;
5  using namespace l1ScoutingRun3;
6  srcInputTag = iConfig.getParameter<InputTag>("srcInputTag");
7  enableAllSums_ = iConfig.getUntrackedParameter<bool>("enableAllSums", false);
8  debug_ = iConfig.getUntrackedParameter<bool>("debug", false);
9 
10  orbitBufferJets_ = std::vector<std::vector<Jet>>(3565);
11  orbitBufferEGammas_ = std::vector<std::vector<EGamma>>(3565);
12  orbitBufferTaus_ = std::vector<std::vector<Tau>>(3565);
13  orbitBufferEtSums_ = std::vector<std::vector<BxSums>>(3565);
14 
15  nJetsOrbit_ = 0;
16  nEGammasOrbit_ = 0;
17  nTausOrbit_ = 0;
18  nEtSumsOrbit_ = 0;
19 
20  produces<JetOrbitCollection>().setBranchAlias("JetOrbitCollection");
21  produces<TauOrbitCollection>().setBranchAlias("TauOrbitCollection");
22  produces<EGammaOrbitCollection>().setBranchAlias("EGammaOrbitCollection");
23  produces<BxSumsOrbitCollection>().setBranchAlias("BxSumsOrbitCollection");
24 
25  rawToken = consumes<SDSRawDataCollection>(srcInputTag);
26 }
27 
29 
31  using namespace edm;
32  using namespace l1ScoutingRun3;
33 
34  Handle<SDSRawDataCollection> ScoutingRawDataCollection;
35  iEvent.getByToken(rawToken, ScoutingRawDataCollection);
36 
37  const FEDRawData& sourceRawData = ScoutingRawDataCollection->FEDData(SDSNumbering::CaloSDSID);
38  size_t orbitSize = sourceRawData.size();
39 
40  std::unique_ptr<JetOrbitCollection> unpackedJets(new JetOrbitCollection);
41  std::unique_ptr<TauOrbitCollection> unpackedTaus(new TauOrbitCollection);
42  std::unique_ptr<EGammaOrbitCollection> unpackedEGammas(new EGammaOrbitCollection);
43  std::unique_ptr<BxSumsOrbitCollection> unpackedEtSums(new BxSumsOrbitCollection);
44 
45  if ((sourceRawData.size() == 0) && debug_) {
46  std::cout << "No raw data for CALO source\n";
47  }
48 
49  // unpack current orbit and store data into the orbitBufferr
50  unpackOrbit(sourceRawData.data(), orbitSize);
51 
52  // fill orbit collection and clear the Bx buffer vector
53  unpackedJets->fillAndClear(orbitBufferJets_, nJetsOrbit_);
54  unpackedEGammas->fillAndClear(orbitBufferEGammas_, nEGammasOrbit_);
55  unpackedTaus->fillAndClear(orbitBufferTaus_, nTausOrbit_);
56  unpackedEtSums->fillAndClear(orbitBufferEtSums_, nEtSumsOrbit_);
57 
58  // store collections in the event
59  iEvent.put(std::move(unpackedJets));
60  iEvent.put(std::move(unpackedTaus));
61  iEvent.put(std::move(unpackedEGammas));
62  iEvent.put(std::move(unpackedEtSums));
63 }
64 
65 void ScCaloRawToDigi::unpackOrbit(const unsigned char* buf, size_t len) {
66  using namespace l1ScoutingRun3;
67 
68  // reset counters
69  nJetsOrbit_ = 0;
70  nEGammasOrbit_ = 0;
71  nTausOrbit_ = 0;
72  nEtSumsOrbit_ = 0;
73 
74  size_t pos = 0;
75 
76  while (pos < len) {
77  assert(pos + sizeof(demux::block) <= len);
78 
79  demux::block* bl = (demux::block*)(buf + pos);
80  pos += sizeof(demux::block);
81 
82  assert(pos <= len);
83  uint32_t orbit = bl->orbit & 0x7FFFFFFF;
84  uint32_t bx = bl->bx;
85 
86  if (debug_) {
87  std::cout << "CALO Orbit " << orbit << ", BX -> " << bx << std::endl;
88  }
89 
90  // unpack jets from first link
91  if (debug_)
92  std::cout << "--- Jets link 1 ---\n";
93  unpackLinkJets(bl->jet1, bx);
94 
95  // unpack jets from second link
96  if (debug_)
97  std::cout << "--- Jets link 2 ---\n";
98  unpackLinkJets(bl->jet2, bx);
99 
100  // unpack eg from first link
101  if (debug_)
102  std::cout << "--- E/g link 1 ---\n";
104 
105  // unpack eg from second link link
106  if (debug_)
107  std::cout << "--- E/g link 2 ---\n";
109 
110  // unpack taus from first link
111  if (debug_)
112  std::cout << "--- Taus link 1 ---\n";
113  unpackLinkTaus(bl->tau1, bx);
114 
115  // unpack taus from second link
116  if (debug_)
117  std::cout << "--- Taus link 2 ---\n";
118  unpackLinkTaus(bl->tau2, bx);
119 
120  // unpack et sums
121  if (debug_)
122  std::cout << "--- Sums ---\n";
123  unpackEtSums(bl->sum, bx);
124 
125  } // end of bx objects
126 }
127 
128 void ScCaloRawToDigi::unpackLinkJets(uint32_t* dataBlock, int bx) {
129  using namespace l1ScoutingRun3;
130 
131  int32_t ET(0), Eta(0), Phi(0), Qual(0);
132  for (uint32_t i = 0; i < 6; i++) {
133  ET = ((dataBlock[i] >> demux::shiftsJet::ET) & demux::masksJet::ET);
134 
135  if (ET != 0) {
136  Eta = ((dataBlock[i] >> demux::shiftsJet::eta) & demux::masksJet::eta);
137  Phi = ((dataBlock[i] >> demux::shiftsJet::phi) & demux::masksJet::phi);
138  Qual = ((dataBlock[i] >> demux::shiftsJet::qual) & demux::masksJet::qual);
139 
140  if (Eta > 127)
141  Eta = Eta - 256;
142 
143  Jet jet(ET, Eta, Phi, Qual);
144  orbitBufferJets_[bx].push_back(jet);
145  nJetsOrbit_++;
146 
147  if (debug_) {
148  std::cout << "Jet " << i << std::endl;
149  std::cout << " Raw: 0x" << std::hex << dataBlock[i] << std::dec << std::endl;
150  printJet(jet);
151  }
152  }
153  } // end link jets unpacking loop
154 }
155 
156 void ScCaloRawToDigi::unpackLinkEGammas(uint32_t* dataBlock, int bx) {
157  using namespace l1ScoutingRun3;
158 
159  int32_t ET(0), Eta(0), Phi(0), Iso(0);
160  for (uint32_t i = 0; i < 6; i++) {
161  ET = ((dataBlock[i] >> demux::shiftsEGamma::ET) & demux::masksEGamma::ET);
162  if (ET != 0) {
163  Eta = ((dataBlock[i] >> demux::shiftsEGamma::eta) & demux::masksEGamma::eta);
165  Iso = ((dataBlock[i] >> demux::shiftsEGamma::iso) & demux::masksEGamma::iso);
166 
167  if (Eta > 127)
168  Eta = Eta - 256;
169 
170  EGamma eGamma(ET, Eta, Phi, Iso);
171  orbitBufferEGammas_[bx].push_back(eGamma);
172  nEGammasOrbit_++;
173 
174  if (debug_) {
175  std::cout << "E/g " << i << std::endl;
176  std::cout << " Raw: 0x" << std::hex << dataBlock[i] << std::dec << std::endl;
177  printEGamma(eGamma);
178  }
179  }
180  } // end link e/gammas unpacking loop
181 }
182 
183 void ScCaloRawToDigi::unpackLinkTaus(uint32_t* dataBlock, int bx) {
184  using namespace l1ScoutingRun3;
185 
186  int32_t ET(0), Eta(0), Phi(0), Iso(0);
187  for (uint32_t i = 0; i < 6; i++) {
188  ET = ((dataBlock[i] >> demux::shiftsTau::ET) & demux::masksTau::ET);
189  if (ET != 0) {
190  Eta = ((dataBlock[i] >> demux::shiftsTau::eta) & demux::masksTau::eta);
191  Phi = ((dataBlock[i] >> demux::shiftsTau::phi) & demux::masksTau::phi);
192  Iso = ((dataBlock[i] >> demux::shiftsTau::iso) & demux::masksTau::iso);
193 
194  if (Eta > 127)
195  Eta = Eta - 256;
196 
197  Tau tau(ET, Eta, Phi, Iso);
198  orbitBufferTaus_[bx].push_back(tau);
199  nTausOrbit_++;
200 
201  if (debug_) {
202  std::cout << "Tau " << i << std::endl;
203  std::cout << " Raw: 0x" << std::hex << dataBlock[i] << std::dec << std::endl;
204  printTau(tau);
205  }
206  }
207  } // end link taus unpacking loop
208 }
209 
210 void ScCaloRawToDigi::unpackEtSums(uint32_t* dataBlock, int bx) {
211  using namespace l1ScoutingRun3;
212 
213  BxSums bxSums;
214 
215  int32_t ETEt(0), ETEttem(0), ETMinBiasHFP0(0); // ET
216  int32_t HTEt(0), HTtowerCount(0), HTMinBiasHFM0(0); // HT
217  int32_t ETmissEt(0), ETmissPhi(0), ETmissASYMET(0), ETmissMinBiasHFP1(0); //ETMiss
218  int32_t HTmissEt(0), HTmissPhi(0), HTmissASYMHT(0), HTmissMinBiasHFM1(0); // HTMiss
219  int32_t ETHFmissEt(0), ETHFmissPhi(0), ETHFmissASYMETHF(0), ETHFmissCENT(0); // ETHFMiss
220  int32_t HTHFmissEt(0), HTHFmissPhi(0), HTHFmissASYMHTHF(0), HTHFmissCENT(0); // HTHFMiss
221 
222  // ET block
223  ETEt = ((dataBlock[0] >> demux::shiftsESums::ETEt) & demux::masksESums::ETEt);
224  ETEttem = ((dataBlock[0] >> demux::shiftsESums::ETEttem) & demux::masksESums::ETEttem);
225 
226  bxSums.setHwTotalEt(ETEt);
227  bxSums.setHwTotalEtEm(ETEttem);
228 
229  // HT block
230  HTEt = ((dataBlock[1] >> demux::shiftsESums::HTEt) & demux::masksESums::HTEt);
231 
232  bxSums.setHwTotalHt(HTEt);
233 
234  // ETMiss block
235  ETmissEt = ((dataBlock[2] >> demux::shiftsESums::ETmissEt) & demux::masksESums::ETmissEt);
236  ETmissPhi = ((dataBlock[2] >> demux::shiftsESums::ETmissPhi) & demux::masksESums::ETmissPhi);
237 
238  if (ETmissEt > 0) {
239  bxSums.setHwMissEt(ETmissEt);
240  bxSums.setHwMissEtPhi(ETmissPhi);
241  }
242 
243  // HTMiss block
244  HTmissEt = ((dataBlock[3] >> demux::shiftsESums::HTmissEt) & demux::masksESums::HTmissEt);
245  HTmissPhi = ((dataBlock[3] >> demux::shiftsESums::HTmissPhi) & demux::masksESums::HTmissPhi);
246 
247  if (HTmissEt > 0) {
248  bxSums.setHwMissHt(HTmissEt);
249  bxSums.setHwMissHtPhi(HTmissPhi);
250  }
251 
252  // ETHFMiss block
253  ETHFmissEt = ((dataBlock[4] >> demux::shiftsESums::ETHFmissEt) & demux::masksESums::ETHFmissEt);
254  ETHFmissPhi = ((dataBlock[4] >> demux::shiftsESums::ETHFmissPhi) & demux::masksESums::ETHFmissPhi);
255 
256  if (ETHFmissEt > 0) {
257  bxSums.setHwMissEtHF(ETHFmissEt);
258  bxSums.setHwMissEtHFPhi(ETHFmissPhi);
259  }
260 
261  // HTHFMiss block
262  HTHFmissEt = ((dataBlock[5] >> demux::shiftsESums::ETHFmissEt) & demux::masksESums::ETHFmissEt);
263  HTHFmissPhi = ((dataBlock[5] >> demux::shiftsESums::ETHFmissPhi) & demux::masksESums::ETHFmissPhi);
264 
265  if (HTHFmissEt > 0) {
266  bxSums.setHwMissHtHF(HTHFmissEt);
267  bxSums.setHwMissHtHFPhi(HTHFmissPhi);
268  }
269 
270  // Insert additional sums
271  if (enableAllSums_) {
272  // ET block
273  ETMinBiasHFP0 = ((dataBlock[0] >> demux::shiftsESums::ETMinBiasHF) & demux::masksESums::ETMinBiasHF);
274  bxSums.setMinBiasHFP0(ETMinBiasHFP0);
275 
276  // HT block
277  HTtowerCount = ((dataBlock[1] >> demux::shiftsESums::HTtowerCount) & demux::masksESums::HTtowerCount);
278  HTMinBiasHFM0 = ((dataBlock[1] >> demux::shiftsESums::HTMinBiasHF) & demux::masksESums::HTMinBiasHF);
279 
280  bxSums.setTowerCount(HTtowerCount);
281  bxSums.setMinBiasHFM0(HTMinBiasHFM0);
282 
283  // ET Miss block
284  ETmissASYMET = ((dataBlock[2] >> demux::shiftsESums::ETmissASYMET) & demux::masksESums::ETmissASYMET);
285  ETmissMinBiasHFP1 = ((dataBlock[2] >> demux::shiftsESums::ETmissMinBiasHF) & demux::masksESums::ETmissMinBiasHF);
286  bxSums.setHwAsymEt(ETmissASYMET);
287  bxSums.setMinBiasHFP1(ETmissMinBiasHFP1);
288 
289  // HT Miss block
290  HTmissASYMHT = ((dataBlock[3] >> demux::shiftsESums::HTmissASYMHT) & demux::masksESums::HTmissASYMHT);
291  HTmissMinBiasHFM1 = ((dataBlock[3] >> demux::shiftsESums::HTmissMinBiasHF) & demux::masksESums::HTmissMinBiasHF);
292 
293  bxSums.setHwAsymHt(HTmissASYMHT);
294  bxSums.setMinBiasHFM1(HTmissMinBiasHFM1);
295 
296  // ETHFMiss
297  ETHFmissASYMETHF = ((dataBlock[4] >> demux::shiftsESums::ETHFmissASYMETHF) & demux::masksESums::ETHFmissASYMETHF);
298  ETHFmissCENT = ((dataBlock[4] >> demux::shiftsESums::ETHFmissCENT) & demux::masksESums::ETHFmissCENT);
299 
300  bxSums.setHwAsymEtHF(ETHFmissASYMETHF);
301 
302  // HTHFMiss
303  HTHFmissASYMHTHF = ((dataBlock[5] >> demux::shiftsESums::ETHFmissASYMETHF) & demux::masksESums::ETHFmissASYMETHF);
304  HTHFmissCENT = ((dataBlock[5] >> demux::shiftsESums::ETHFmissCENT) & demux::masksESums::ETHFmissCENT);
305 
306  bxSums.setHwAsymHtHF(HTHFmissASYMHTHF);
307  bxSums.setCentrality((HTHFmissCENT << 4) + ETHFmissCENT);
308  }
309 
310  orbitBufferEtSums_[bx].push_back(bxSums);
311  nEtSumsOrbit_ += 1;
312 
313  if (debug_) {
314  std::cout << "Raw frames:\n";
315  for (int frame = 0; frame < 6; frame++) {
316  std::cout << " frame " << frame << ": 0x" << std::hex << dataBlock[frame] << std::dec << std::endl;
317  }
318  printBxSums(bxSums);
319  }
320 }
321 
324  desc.setUnknown();
325  descriptions.addDefault(desc);
326 }
327 
void unpackOrbit(const unsigned char *buf, size_t len)
void setHwTotalEtEm(int hwTotalEtEm)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void setHwMissEtHFPhi(int hwMissEtHFPhi)
void produce(edm::Event &, const edm::EventSetup &) override
void setMinBiasHFM0(int minBiasHFM0)
void setMinBiasHFP0(int minBiasHFP0)
void unpackLinkJets(uint32_t *dataBlock, int bx)
void setHwAsymHtHF(int hwAsymHtHF)
~ScCaloRawToDigi() override
const FEDRawData & FEDData(int sourceId) const
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:48
void printJet(const Jet &jet, std::ostream &outs=std::cout)
void printEGamma(const EGamma &eGamma, std::ostream &outs=std::cout)
assert(be >=bs)
void printBxSums(const BxSums &sums, std::ostream &outs=std::cout)
std::vector< std::vector< l1ScoutingRun3::Tau > > orbitBufferTaus_
T getUntrackedParameter(std::string const &, T const &) const
ScCaloRawToDigi(const edm::ParameterSet &)
void setHwMissEtHF(int hwMissEtHF)
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
void setHwAsymHt(int hwAsymHt)
Definition: Jet.py:1
void setHwMissHtHF(int hwMissHtHF)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void unpackEtSums(uint32_t *dataBlock, int bx)
void setHwTotalHt(int hwTotalHt)
std::vector< std::vector< l1ScoutingRun3::BxSums > > orbitBufferEtSums_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static constexpr int CaloSDSID
Definition: SDSNumbering.h:18
void setHwAsymEtHF(int hwAsymEtHF)
void setCentrality(int centrality)
std::vector< std::vector< l1ScoutingRun3::EGamma > > orbitBufferEGammas_
void setTowerCount(int towerCount)
void setHwMissHt(int hwMissHt)
Definition: Tau.py:1
void setHwAsymEt(int hwAsymEt)
void setMinBiasHFM1(int minBiasHFM1)
edm::EDGetToken rawToken
void setMinBiasHFP1(int minBiasHFP1)
void printTau(const Tau &tau, std::ostream &outs=std::cout)
edm::InputTag srcInputTag
std::vector< std::vector< l1ScoutingRun3::Jet > > orbitBufferJets_
void unpackLinkTaus(uint32_t *dataBlock, int bx)
void setHwMissEtPhi(int hwMissEtPhi)
HLT enums.
void unpackLinkEGammas(uint32_t *dataBlock, int bx)
void setHwTotalEt(int hwTotalEt)
void setHwMissHtPhi(int hwMissHtPhi)
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:24
#define ET
void setHwMissEt(int hwMissEt)
void setHwMissHtHFPhi(int hwMissHtHFPhi)
def move(src, dest)
Definition: eostools.py:511