CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CSCDigiToRaw.cc
Go to the documentation of this file.
1 
15 #include "boost/dynamic_bitset.hpp"
23 
24 #include <algorithm>
25 
26 using namespace edm;
27 using namespace std;
28 
30  : alctWindowMin_(pset.getParameter<int>("alctWindowMin")),
31  alctWindowMax_(pset.getParameter<int>("alctWindowMax")),
32  clctWindowMin_(pset.getParameter<int>("clctWindowMin")),
33  clctWindowMax_(pset.getParameter<int>("clctWindowMax")),
34  preTriggerWindowMin_(pset.getParameter<int>("preTriggerWindowMin")),
35  preTriggerWindowMax_(pset.getParameter<int>("preTriggerWindowMax")),
36  // pre-LS1 - '2005'. post-LS1 - '2013'
37  formatVersion_(pset.getParameter<unsigned>("formatVersion")),
38  // don't check for consistency with trig primitives
39  // overrides usePreTriggers
40  packEverything_(pset.getParameter<bool>("packEverything")),
41  usePreTriggers_(pset.getParameter<bool>("usePreTriggers")),
42  packByCFEB_(pset.getParameter<bool>("packByCFEB")) {}
43 
45  CSCDetId chamberId = CSCDigiToRawAccept::chamberID(cscDetId);
46  // find the entry into the map
47  map<CSCDetId, CSCEventData>::iterator chamberMapItr = info.theChamberDataMap.find(chamberId);
48  if (chamberMapItr == info.theChamberDataMap.end()) {
49  // make an entry, telling it the correct chamberType
50  int chamberType = chamberId.iChamberType();
51  chamberMapItr = info.theChamberDataMap
52  .insert(pair<CSCDetId, CSCEventData>(chamberId, CSCEventData(chamberType, info.formatVersion_)))
53  .first;
54  }
55  CSCEventData& cscData = chamberMapItr->second;
56  cscData.dmbHeader()->setCrateAddress(info.theElectronicsMap->crate(cscDetId), info.theElectronicsMap->dmb(cscDetId));
57 
58  if (info.formatVersion_ == 2013) {
59  // Set DMB version field to distinguish between ME11s and other chambers (ME11 - 2, others - 1)
60  bool me11 = ((chamberId.station() == 1) && (chamberId.ring() == 4)) ||
61  ((chamberId.station() == 1) && (chamberId.ring() == 1));
62  if (me11) {
63  cscData.dmbHeader()->setdmbVersion(2);
64  } else {
65  cscData.dmbHeader()->setdmbVersion(1);
66  }
67  }
68  return cscData;
69 }
70 
72  const CSCCLCTPreTriggerCollection* preTriggers,
73  const CSCCLCTPreTriggerDigiCollection* preTriggerDigis,
74  FindEventDataInfo& fedInfo) const {
75  //iterate over chambers with strip digis in them
76  for (CSCStripDigiCollection::DigiRangeIterator j = stripDigis.begin(); j != stripDigis.end(); ++j) {
77  CSCDetId cscDetId = (*j).first;
78  // only digitize if there are pre-triggers
79 
80  // determine where the pretriggers are
81  std::vector<bool> preTriggerInCFEB;
82  preTriggerInCFEB.resize(CSCConstants::MAX_CFEBS_RUN2);
83 
84  // pretrigger flag must be set and the pretrigger collection must be nonzero!
85  const bool usePreTriggers = usePreTriggers_ and preTriggers != nullptr;
86  if (!usePreTriggers || packEverything_ ||
87  (usePreTriggers && CSCDigiToRawAccept::accept(cscDetId,
88  *preTriggerDigis,
92  preTriggerInCFEB))) {
93  bool me1a = (cscDetId.station() == 1) && (cscDetId.ring() == 4);
94  bool zplus = (cscDetId.endcap() == 1);
95  bool me1b = (cscDetId.station() == 1) && (cscDetId.ring() == 1);
96 
97  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
98 
99  std::vector<CSCStripDigi>::const_iterator digiItr = (*j).second.first;
100  std::vector<CSCStripDigi>::const_iterator last = (*j).second.second;
101  for (; digiItr != last; ++digiItr) {
102  CSCStripDigi digi = *digiItr;
103  int strip = digi.getStrip();
104  int cfeb = digi.getCFEB();
105  // CSC strip digis in ME1/a have CFEB number 0, 1, or 2
106  // But a pretrigger in ME1/a has CFEB number 4, 5, or 6 (+4)
107  if (me1a)
109 
110  // At this point, if we are packing by CFEBs and there is no
111  // pretrigger in this CFEB, ignore this strip digi
112  if (packByCFEB_ and not preTriggerInCFEB[cfeb])
113  continue;
114 
115  // From LS1 on ME1a strips are unganged
116  if (fedInfo.formatVersion_ == 2013) {
117  if (me1a && zplus) {
118  digi.setStrip(CSCConstants::NUM_STRIPS_ME1A_UNGANGED + 1 - strip); // 1-48 -> 48-1
119  }
120  if (me1b && !zplus) {
121  digi.setStrip(CSCConstants::NUM_STRIPS_ME1B + 1 - strip); // 1-64 -> 64-1
122  }
123  if (me1a) {
124  strip = digi.getStrip(); // reset back 1-16 to 65-80 digi
126  }
127  }
128  // During Run-1 ME1a strips are triple-ganged
129  else {
130  if (me1a && zplus) {
131  digi.setStrip(CSCConstants::NUM_STRIPS_ME1A_GANGED + 1 - strip); // 1-16 -> 16-1
132  }
133  if (me1b && !zplus) {
134  digi.setStrip(CSCConstants::NUM_STRIPS_ME1B + 1 - strip); // 1-64 -> 64-1
135  }
136  if (me1a) {
137  strip = digi.getStrip(); // reset back 1-16 to 65-80 digi
139  }
140  }
141  cscData.add(digi, cscDetId.layer());
142  }
143  }
144  }
145 }
146 
148  const CSCALCTDigiCollection& alctDigis,
149  FindEventDataInfo& fedInfo) const {
150  add(alctDigis, fedInfo);
151  for (CSCWireDigiCollection::DigiRangeIterator j = wireDigis.begin(); j != wireDigis.end(); ++j) {
152  CSCDetId cscDetId = (*j).first;
155  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
156  std::vector<CSCWireDigi>::const_iterator digiItr = (*j).second.first;
157  std::vector<CSCWireDigi>::const_iterator last = (*j).second.second;
158  for (; digiItr != last; ++digiItr) {
159  cscData.add(*digiItr, cscDetId.layer());
160  }
161  }
162  }
163 }
164 
165 void CSCDigiToRaw::add(const CSCComparatorDigiCollection& comparatorDigis,
166  const CSCCLCTDigiCollection& clctDigis,
167  FindEventDataInfo& fedInfo) const {
168  add(clctDigis, fedInfo);
169  for (auto const& j : comparatorDigis) {
170  CSCDetId cscDetId = j.first;
171  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
174  bool me1a = (cscDetId.station() == 1) && (cscDetId.ring() == 4);
175 
176  /*
177  Add the comparator digi to the fedInfo.
178  Move ME1/A comparators from CFEB=0 to CFEB=4 if this has not
179  been done already.
180  Consider the case for triple-ganged and unganged ME1A strips
181  */
182  for (auto digi = j.second.first; digi != j.second.second; ++digi) {
183  if (fedInfo.formatVersion_ == 2013) {
184  // unganged case
185  if (me1a && digi->getStrip() <= CSCConstants::NUM_STRIPS_ME1A_UNGANGED) {
186  CSCComparatorDigi digi_corr(
187  CSCConstants::NUM_STRIPS_ME1B + digi->getStrip(), digi->getComparator(), digi->getTimeBinWord());
188  cscData.add(digi_corr, cscDetId);
189  } else {
190  cscData.add(*digi, cscDetId);
191  }
192  }
193  // triple-ganged case
194  else {
195  if (me1a && digi->getStrip() <= CSCConstants::NUM_STRIPS_ME1A_GANGED) {
196  CSCComparatorDigi digi_corr(
197  CSCConstants::NUM_STRIPS_ME1B + digi->getStrip(), digi->getComparator(), digi->getTimeBinWord());
198  cscData.add(digi_corr, cscDetId.layer());
199  } else {
200  cscData.add(*digi, cscDetId.layer());
201  }
202  }
203  }
204  }
205  }
206 }
207 
208 void CSCDigiToRaw::add(const CSCALCTDigiCollection& alctDigis, FindEventDataInfo& fedInfo) const {
209  for (CSCALCTDigiCollection::DigiRangeIterator j = alctDigis.begin(); j != alctDigis.end(); ++j) {
210  CSCDetId cscDetId = (*j).first;
211  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
212 
213  cscData.add(std::vector<CSCALCTDigi>((*j).second.first, (*j).second.second));
214  }
215 }
216 
217 void CSCDigiToRaw::add(const CSCCLCTDigiCollection& clctDigis, FindEventDataInfo& fedInfo) const {
218  for (CSCCLCTDigiCollection::DigiRangeIterator j = clctDigis.begin(); j != clctDigis.end(); ++j) {
219  CSCDetId cscDetId = (*j).first;
220  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
221 
222  bool me11a = cscDetId.station() == 1 && cscDetId.ring() == 4;
223  //CLCTs are packed by chamber not by A/B parts in ME11
224  //me11a appears only in simulation with SLHC algorithm settings
225  //without the shift, it's impossible to distinguish A and B parts
226  if (me11a && fedInfo.formatVersion_ == 2013) {
227  std::vector<CSCCLCTDigi> shiftedDigis((*j).second.first, (*j).second.second);
228  for (std::vector<CSCCLCTDigi>::iterator iC = shiftedDigis.begin(); iC != shiftedDigis.end(); ++iC) {
229  if (iC->getCFEB() < CSCConstants::NUM_CFEBS_ME1A_UNGANGED) { //sanity check, mostly
230  (*iC) = CSCCLCTDigi(iC->isValid(),
231  iC->getQuality(),
232  iC->getPattern(),
233  iC->getStripType(),
234  iC->getBend(),
235  iC->getStrip(),
236  iC->getCFEB() + CSCConstants::NUM_CFEBS_ME1B,
237  iC->getBX(),
238  iC->getTrknmb(),
239  iC->getFullBX());
240  }
241  }
242  cscData.add(shiftedDigis);
243  } else {
244  cscData.add(std::vector<CSCCLCTDigi>((*j).second.first, (*j).second.second));
245  }
246  }
247 }
248 
249 void CSCDigiToRaw::add(const CSCCorrelatedLCTDigiCollection& corrLCTDigis, FindEventDataInfo& fedInfo) const {
250  for (CSCCorrelatedLCTDigiCollection::DigiRangeIterator j = corrLCTDigis.begin(); j != corrLCTDigis.end(); ++j) {
251  CSCDetId cscDetId = (*j).first;
252  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
253 
254  bool me11a = cscDetId.station() == 1 && cscDetId.ring() == 4;
255  //LCTs are packed by chamber not by A/B parts in ME11
256  //me11a appears only in simulation with SLHC algorithm settings
257  //without the shift, it's impossible to distinguish A and B parts
258  if (me11a && fedInfo.formatVersion_ == 2013) {
259  std::vector<CSCCorrelatedLCTDigi> shiftedDigis((*j).second.first, (*j).second.second);
260  for (std::vector<CSCCorrelatedLCTDigi>::iterator iC = shiftedDigis.begin(); iC != shiftedDigis.end(); ++iC) {
261  if (iC->getStrip() < CSCConstants::NUM_HALF_STRIPS_ME1A_UNGANGED) { //sanity check, mostly
262  (*iC) = CSCCorrelatedLCTDigi(iC->getTrknmb(),
263  iC->isValid(),
264  iC->getQuality(),
265  iC->getKeyWG(),
266  iC->getStrip() + CSCConstants::NUM_HALF_STRIPS_ME1B,
267  iC->getPattern(),
268  iC->getBend(),
269  iC->getBX(),
270  iC->getMPCLink(),
271  iC->getBX0(),
272  iC->getSyncErr(),
273  iC->getCSCID());
274  }
275  }
276  cscData.add(shiftedDigis);
277  } else {
278  cscData.add(std::vector<CSCCorrelatedLCTDigi>((*j).second.first, (*j).second.second));
279  }
280  }
281 }
282 
283 void CSCDigiToRaw::add(const CSCShowerDigiCollection& cscShowerDigis, FindEventDataInfo& fedInfo) const {
284  for (const auto& shower : cscShowerDigis) {
285  const CSCDetId& cscDetId = shower.first;
286  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
287 
288  cscData.add(std::vector<CSCShowerDigi>(shower.second.first, shower.second.second));
289  }
290 }
291 
292 void CSCDigiToRaw::add(const GEMPadDigiClusterCollection& gemPadClusters, FindEventDataInfo& fedInfo) const {
293  for (const auto& jclus : gemPadClusters) {
294  const GEMDetId& gemDetId = jclus.first;
295 
296  const int zendcap = gemDetId.region() == 1 ? 1 : 2;
297  CSCDetId cscDetId(zendcap, gemDetId.station(), 1, gemDetId.chamber(), 0);
298  CSCEventData& cscData = findEventData(cscDetId, fedInfo);
299 
300  cscData.add(std::vector<GEMPadDigiCluster>(jclus.second.first, jclus.second.second), gemDetId);
301  }
302 }
303 
305  const CSCWireDigiCollection& wireDigis,
306  const CSCComparatorDigiCollection& comparatorDigis,
307  const CSCALCTDigiCollection& alctDigis,
308  const CSCCLCTDigiCollection& clctDigis,
309  const CSCCLCTPreTriggerCollection* preTriggers,
310  const CSCCLCTPreTriggerDigiCollection* preTriggerDigis,
311  const CSCCorrelatedLCTDigiCollection& correlatedLCTDigis,
312  const CSCShowerDigiCollection* cscShowerDigis,
313  const GEMPadDigiClusterCollection* gemPadDigiClusters,
314  FEDRawDataCollection& fed_buffers,
315  const CSCChamberMap* mapping,
316  const EventID& eid) const {
317  //get fed object from fed_buffers
318  // make a map from the index of a chamber to the event data from it
319  FindEventDataInfo fedInfo{mapping, formatVersion_};
320  add(stripDigis, preTriggers, preTriggerDigis, fedInfo);
321  add(wireDigis, alctDigis, fedInfo);
322  add(comparatorDigis, clctDigis, fedInfo);
323  add(correlatedLCTDigis, fedInfo);
324 
325  // Starting Run-3, the CSC DAQ will pack/unpack CSC showers
326  if (cscShowerDigis) {
327  add(*cscShowerDigis, fedInfo);
328  }
329 
330  // Starting Run-3, the CSC DAQ will pack/unpack GEM clusters
331  if (gemPadDigiClusters) {
332  add(*gemPadDigiClusters, fedInfo);
333  }
334  int l1a = eid.event(); //need to add increments or get it from lct digis
335  int bx = l1a; //same as above
336 
337  if (fedInfo.formatVersion_ == 2005)
338  {
339  std::map<int, CSCDCCEventData> dccMap;
340  //idcc goes from startingFed to startingFED+7
341  for (int idcc = FEDNumbering::MINCSCFEDID; idcc <= FEDNumbering::MAXCSCFEDID; ++idcc) {
342  dccMap.insert(std::pair<int, CSCDCCEventData>(idcc, CSCDCCEventData(idcc, CSCConstants::NUM_DDUS, bx, l1a)));
343  }
344 
345  for (int idcc = FEDNumbering::MINCSCFEDID; idcc <= FEDNumbering::MAXCSCFEDID; ++idcc) {
346  // for every chamber with data, add to a DDU in this DCC Event
347  for (map<CSCDetId, CSCEventData>::iterator chamberItr = fedInfo.theChamberDataMap.begin();
348  chamberItr != fedInfo.theChamberDataMap.end();
349  ++chamberItr) {
350  int indexDCC = mapping->slink(chamberItr->first);
351  if (indexDCC == idcc) {
352  //FIXME (What does this mean? Is something wrong?)
353  std::map<int, CSCDCCEventData>::iterator dccMapItr = dccMap.find(indexDCC);
354  if (dccMapItr == dccMap.end()) {
355  throw cms::Exception("CSCDigiToRaw") << "Bad DCC number:" << indexDCC;
356  }
357  // get id's based on ChamberId from mapping
358 
359  int dduId = mapping->ddu(chamberItr->first);
360  int dduSlot = mapping->dduSlot(chamberItr->first);
361  int dduInput = mapping->dduInput(chamberItr->first);
362  int dmbId = mapping->dmb(chamberItr->first);
363  dccMapItr->second.addChamber(chamberItr->second, dduId, dduSlot, dduInput, dmbId, formatVersion_);
364  }
365  }
366  }
367 
368  // FIXME: FEDRawData size set to 2*64 to add FED header and trailer
369  for (std::map<int, CSCDCCEventData>::iterator dccMapItr = dccMap.begin(); dccMapItr != dccMap.end(); ++dccMapItr) {
370  boost::dynamic_bitset<> dccBits = dccMapItr->second.pack();
371  FEDRawData& fedRawData = fed_buffers.FEDData(dccMapItr->first);
372  fedRawData.resize(dccBits.size());
373  //fill data with dccEvent
374  bitset_utilities::bitsetToChar(dccBits, fedRawData.data());
375  FEDTrailer cscFEDTrailer(fedRawData.data() + (fedRawData.size() - 8));
376  cscFEDTrailer.set(fedRawData.data() + (fedRawData.size() - 8),
377  fedRawData.size() / 8,
378  evf::compute_crc(fedRawData.data(), fedRawData.size()),
379  0,
380  0);
381  }
382  }
384  else if (formatVersion_ == 2013) {
385  std::map<int, CSCDDUEventData> dduMap;
386  unsigned int ddu_fmt_version = 0x7;
387  const unsigned postLS1_map[] = {841, 842, 843, 844, 845, 846, 847, 848, 849, 831, 832, 833,
388  834, 835, 836, 837, 838, 839, 861, 862, 863, 864, 865, 866,
389  867, 868, 869, 851, 852, 853, 854, 855, 856, 857, 858, 859};
390 
392  for (unsigned int i = 0; i < 36; i++) {
393  unsigned int iddu = postLS1_map[i];
394  // make a new one
395  CSCDDUHeader newDDUHeader(bx, l1a, iddu, ddu_fmt_version);
396 
397  dduMap.insert(std::pair<int, CSCDDUEventData>(iddu, CSCDDUEventData(newDDUHeader)));
398  }
399 
401  for (unsigned int i = 0; i < 36; i++) {
402  unsigned int iddu = postLS1_map[i];
403  // for every chamber with data, add to a DDU in this DCC Event
404  for (map<CSCDetId, CSCEventData>::iterator chamberItr = fedInfo.theChamberDataMap.begin();
405  chamberItr != fedInfo.theChamberDataMap.end();
406  ++chamberItr) {
407  unsigned int indexDDU = mapping->slink(chamberItr->first);
408 
410  // Still preLS1 DCC FEDs mapping
411  if ((indexDDU >= FEDNumbering::MINCSCFEDID) && (indexDDU <= FEDNumbering::MAXCSCFEDID)) {
412  int dduID = mapping->ddu(chamberItr->first); // try to switch to DDU ID mapping
413  // DDU ID is in expectedi post-LS1 FED ID range
414  if ((dduID >= FEDNumbering::MINCSCDDUFEDID) && (dduID <= FEDNumbering::MAXCSCDDUFEDID)) {
415  indexDDU = dduID;
416  } else // Messy
417  {
418  // Lets try to change pre-LS1 1-36 ID range to post-LS1 MINCSCDDUFEDID - MAXCSCDDUFEDID range
419  dduID &= 0xFF;
420  // indexDDU = FEDNumbering::MINCSCDDUFEDID + dduID-1;
421  if ((dduID <= 36) && (dduID > 0))
422  indexDDU = postLS1_map[dduID - 1];
423  }
424  }
425 
426  if (indexDDU == iddu) {
427  std::map<int, CSCDDUEventData>::iterator dduMapItr = dduMap.find(indexDDU);
428  if (dduMapItr == dduMap.end()) {
429  throw cms::Exception("CSCDigiToRaw") << "Bad DDU number:" << indexDDU;
430  }
431  // get id's based on ChamberId from mapping
432 
433  int dduInput = mapping->dduInput(chamberItr->first);
434  int dmbId = mapping->dmb(chamberItr->first);
435  // int crateId = mapping->crate(chamberItr->first);
436  dduMapItr->second.add(chamberItr->second, dmbId, dduInput, formatVersion_);
437  }
438  }
439  }
440 
441  // FIXME: FEDRawData size set to 2*64 to add FED header and trailer
442  for (std::map<int, CSCDDUEventData>::iterator dduMapItr = dduMap.begin(); dduMapItr != dduMap.end(); ++dduMapItr) {
443  boost::dynamic_bitset<> dduBits = dduMapItr->second.pack();
444  FEDRawData& fedRawData = fed_buffers.FEDData(dduMapItr->first);
445  fedRawData.resize(dduBits.size() / 8);
446  //fill data with dduEvent
447  bitset_utilities::bitsetToChar(dduBits, fedRawData.data());
448  FEDTrailer cscFEDTrailer(fedRawData.data() + (fedRawData.size() - 8));
449  cscFEDTrailer.set(fedRawData.data() + (fedRawData.size() - 8),
450  fedRawData.size() / 8,
451  evf::compute_crc(fedRawData.data(), fedRawData.size()),
452  0,
453  0);
454  }
455  }
456 }
01/20/05 A.Tumanov
EventNumber_t event() const
Definition: EventID.h:40
static const TGPicture * info(bool iBackgroundIsBlack)
int getCFEB() const
Get the CFEB number. Counts from 0.
Definition: CSCStripDigi.cc:25
bool accept(const CSCDetId &cscId, const LCTCollection &lcts, int bxMin, int bxMax, int nominalBX)
uint16_t formatVersion_
Definition: CSCDigiToRaw.h:81
bool usePreTriggers_
Definition: CSCDigiToRaw.h:83
const int clctWindowMin_
Definition: CSCDigiToRaw.h:76
void createFedBuffers(const CSCStripDigiCollection &stripDigis, const CSCWireDigiCollection &wireDigis, const CSCComparatorDigiCollection &comparatorDigis, const CSCALCTDigiCollection &alctDigis, const CSCCLCTDigiCollection &clctDigis, const CSCCLCTPreTriggerCollection *preTriggers, const CSCCLCTPreTriggerDigiCollection *preTriggerDigis, const CSCCorrelatedLCTDigiCollection &correlatedLCTDigis, const CSCShowerDigiCollection *showerDigis, const GEMPadDigiClusterCollection *padDigiClusters, FEDRawDataCollection &fed_buffers, const CSCChamberMap *theMapping, const edm::EventID &eid) const
Take a vector of digis and fill the FEDRawDataCollection.
bool packEverything_
Definition: CSCDigiToRaw.h:82
int ddu(const CSCDetId &) const
ddu id for given DetId
const int clctWindowMax_
Definition: CSCDigiToRaw.h:77
int dmb(const CSCDetId &) const
dmb id for given DetId
void setdmbVersion(unsigned int version)
Definition: CSCDMBHeader.h:32
int layer() const
Definition: CSCDetId.h:56
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:45
int getStrip() const
Definition: CSCStripDigi.h:41
int endcap() const
Definition: CSCDetId.h:85
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
void setStrip(int istrip)
Definition: CSCStripDigi.h:65
CSCEventData & findEventData(const CSCDetId &cscDetId, FindEventDataInfo &) const
pick out the correct data object for this chamber
Definition: CSCDigiToRaw.cc:44
void resize(size_t newsize)
Definition: FEDRawData.cc:28
int slink(const CSCDetId &) const
slink id for given DetId
const int preTriggerWindowMin_
Definition: CSCDigiToRaw.h:78
constexpr int region() const
Definition: GEMDetId.h:171
unsigned short compute_crc(unsigned char *buffer, unsigned int bufSize)
Definition: CRC16.h:46
const int preTriggerWindowMax_
Definition: CSCDigiToRaw.h:79
int dduInput(const CSCDetId &) const
ddu input for given DetId
CSCDetId chamberID(const CSCDetId &cscDetId)
CSCDigiToRaw(const edm::ParameterSet &pset)
Constructor.
Definition: CSCDigiToRaw.cc:29
unsigned short iChamberType() const
Definition: CSCDetId.h:96
const CSCChamberMap * theElectronicsMap
Definition: CSCDigiToRaw.h:52
void setCrateAddress(int crate, int dmbId)
Definition: CSCDMBHeader.h:30
const int alctWindowMin_
Definition: CSCDigiToRaw.h:74
int ring() const
Definition: CSCDetId.h:68
const int alctWindowMax_
Definition: CSCDigiToRaw.h:75
constexpr int chamber() const
Definition: GEMDetId.h:183
void bitsetToChar(const boost::dynamic_bitset<> &bs, unsigned char *result)
this method takes bitset obj and returns char * array
void add(const CSCStripDigiCollection &stripDigis, const CSCCLCTPreTriggerCollection *preTriggers, const CSCCLCTPreTriggerDigiCollection *preTriggerDigis, FindEventDataInfo &) const
Definition: CSCDigiToRaw.cc:71
int crate(const CSCDetId &) const
Interface required use in digi-to-raw.
constexpr int station() const
Definition: GEMDetId.h:179
int dduSlot(const CSCDetId &) const
ddu slot for given DetId
void add(const CSCStripDigi &, int layer)
routines to add digis to the data
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:24
int station() const
Definition: CSCDetId.h:79
tuple last
Definition: dqmdumpme.py:56
const CSCDMBHeader * dmbHeader() const
the DAQ motherboard header. A good place for event and chamber info
Definition: CSCEventData.h:73
bool packByCFEB_
Definition: CSCDigiToRaw.h:84
A container for a generic type of digis indexed by some index, implemented with a map&lt;IndexType...