CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCDigiToRaw.cc
Go to the documentation of this file.
1 
13 #include "boost/dynamic_bitset.hpp"
14 #include "boost/foreach.hpp"
23 #include <algorithm>
24 
25 using namespace edm;
26 using namespace std;
27 
28 namespace cscd2r
29 {
31 CSCDetId chamberID(const CSCDetId & cscDetId)
32 {
33  CSCDetId chamberId = cscDetId.chamberId();
34  if (chamberId.ring() ==4)
35  {
36  chamberId = CSCDetId(chamberId.endcap(), chamberId.station(), 1, chamberId.chamber(), 0);
37  }
38  return chamberId;
39 }
40 
41 template<typename LCTCollection>
42 bool accept(const CSCDetId & cscId, const LCTCollection & lcts,
43  int bxMin, int bxMax)
44 {
45  if (bxMin == -999) return true;
46  int nominalBX = 6;
47  CSCDetId chamberId = chamberID(cscId);
48  typename LCTCollection::Range lctRange = lcts.get(chamberId);
49  bool result = false;
50  for (typename LCTCollection::const_iterator lctItr = lctRange.first;
51  lctItr != lctRange.second; ++lctItr)
52  {
53  int bx = lctItr->getBX() - nominalBX;
54  if (bx >= bxMin && bx <= bxMax)
55  {
56  result = true;
57  break;
58  }
59  }
60  return result;
61 }
62 
63 // need to specialize for pretriggers, since they don't have a getBX()
64 template<>
65 bool accept(const CSCDetId & cscId, const CSCCLCTPreTriggerCollection & lcts,
66  int bxMin, int bxMax)
67 {
68  if (bxMin == -999) return true;
69  int nominalBX = 6;
70  CSCDetId chamberId = chamberID(cscId);
71  CSCCLCTPreTriggerCollection::Range lctRange = lcts.get(chamberId);
72  bool result = false;
73  for (CSCCLCTPreTriggerCollection::const_iterator lctItr = lctRange.first;
74  lctItr != lctRange.second; ++lctItr)
75  {
76  int bx = *lctItr - nominalBX;
77  if (bx >= bxMin && bx <= bxMax)
78  {
79  result = true;
80  break;
81  }
82  }
83  return result;
84 }
85 
86 }
87 
88 
90  : alctWindowMin_(pset.getParameter<int>("alctWindowMin")),
91  alctWindowMax_(pset.getParameter<int>("alctWindowMax")),
92  clctWindowMin_(pset.getParameter<int>("clctWindowMin")),
93  clctWindowMax_(pset.getParameter<int>("clctWindowMax")),
94  preTriggerWindowMin_(pset.getParameter<int>("preTriggerWindowMin")),
95  preTriggerWindowMax_(pset.getParameter<int>("preTriggerWindowMax")),
96  theFormatVersion(2005),
97  usePreTriggers(true)
98 {}
99 
101 {
102  theChamberDataMap.clear();
104 }
105 
106 
108 {
109  CSCDetId chamberId = cscd2r::chamberID(cscDetId);
110  // find the entry into the map
111  map<CSCDetId, CSCEventData>::iterator chamberMapItr = theChamberDataMap.find(chamberId);
112  if (chamberMapItr == theChamberDataMap.end())
113  {
114  // make an entry, telling it the correct chamberType
115  int chamberType = chamberId.iChamberType();
116  chamberMapItr = theChamberDataMap.insert(pair<CSCDetId, CSCEventData>(chamberId, CSCEventData(chamberType, theFormatVersion))).first;
117  }
118  CSCEventData & cscData = chamberMapItr->second;
119  cscData.dmbHeader()->setCrateAddress(theElectronicsMap->crate(cscDetId), theElectronicsMap->dmb(cscDetId));
120 
121  if (theFormatVersion == 2013)
122  {
123  // Set DMB version field to distinguish between ME11s and other chambers (ME11 - 2, others - 1)
124  bool me11 = ((chamberId.station()==1) && (chamberId.ring()==4)) || ((chamberId.station()==1) && (chamberId.ring()==1));
125  if (me11)
126  {
127  cscData.dmbHeader()->setdmbVersion(2);
128  }
129  else
130  {
131  cscData.dmbHeader()->setdmbVersion(1);
132  }
133 
134  }
135  return cscData;
136 }
137 
138 
139 
141  const CSCCLCTPreTriggerCollection & preTriggers)
142 { //iterate over chambers with strip digis in them
143  for (CSCStripDigiCollection::DigiRangeIterator j=stripDigis.begin(); j!=stripDigis.end(); ++j)
144  {
145  CSCDetId cscDetId=(*j).first;
146  // only digitize if there are pre-triggers
147 
148  /* !!! Testing. Uncomment for production */
149  if (!usePreTriggers ||
151  {
152  bool me1a = (cscDetId.station()==1) && (cscDetId.ring()==4);
153  bool zplus = (cscDetId.endcap() == 1);
154  bool me1b = (cscDetId.station()==1) && (cscDetId.ring()==1);
155 
156  CSCEventData & cscData = findEventData(cscDetId);
157 
158  std::vector<CSCStripDigi>::const_iterator digiItr = (*j).second.first;
159  std::vector<CSCStripDigi>::const_iterator last = (*j).second.second;
160  for ( ; digiItr != last; ++digiItr)
161  {
162  CSCStripDigi digi = *digiItr;
163  int strip = digi.getStrip();
164  if (theFormatVersion == 2013)
165  {
166  if ( me1a && zplus )
167  {
168  digi.setStrip(49-strip); // 1-48 -> 48-1
169  }
170  if ( me1b && !zplus)
171  {
172  digi.setStrip(65-strip); // 1-64 -> 64-1
173  }
174  if ( me1a )
175  {
176  strip = digi.getStrip(); // reset back 1-16 to 65-80 digi
177  digi.setStrip(strip+64);
178  }
179 
180  }
181  else
182  {
183  if ( me1a && zplus )
184  {
185  digi.setStrip(17-strip); // 1-16 -> 16-1
186  }
187  if ( me1b && !zplus)
188  {
189  digi.setStrip(65-strip); // 1-64 -> 64-1
190  }
191  if ( me1a )
192  {
193  strip = digi.getStrip(); // reset back 1-16 to 65-80 digi
194  digi.setStrip(strip+64);
195  }
196  }
197  cscData.add(digi, cscDetId.layer() );
198  }
199  }
200  }
201 }
202 
203 
205  const CSCALCTDigiCollection & alctDigis)
206 {
207  add(alctDigis);
208  for (CSCWireDigiCollection::DigiRangeIterator j=wireDigis.begin(); j!=wireDigis.end(); ++j)
209  {
210  CSCDetId cscDetId=(*j).first;
211  if (cscd2r::accept(cscDetId, alctDigis, alctWindowMin_, alctWindowMax_))
212  {
213  CSCEventData & cscData = findEventData(cscDetId);
214  std::vector<CSCWireDigi>::const_iterator digiItr = (*j).second.first;
215  std::vector<CSCWireDigi>::const_iterator last = (*j).second.second;
216  for ( ; digiItr != last; ++digiItr)
217  {
218  cscData.add(*digiItr, cscDetId.layer() );
219  }
220  }
221  }
222 
223 }
224 
225 void CSCDigiToRaw::add(const CSCComparatorDigiCollection & comparatorDigis,
226  const CSCCLCTDigiCollection & clctDigis)
227 {
228  add(clctDigis);
229  for (CSCComparatorDigiCollection::DigiRangeIterator j=comparatorDigis.begin(); j!=comparatorDigis.end(); ++j)
230  {
231  CSCDetId cscDetId=(*j).first;
232  CSCEventData & cscData = findEventData(cscDetId);
233  if (cscd2r::accept(cscDetId, clctDigis, clctWindowMin_, clctWindowMax_))
234  {
235  bool me1a = (cscDetId.station()==1) && (cscDetId.ring()==4);
236 
237 
238  BOOST_FOREACH(CSCComparatorDigi digi, (*j).second)
239  {
240  if (theFormatVersion == 2013)
241  {
242  // Move ME1/A comparators from CFEB=0 to CFEB=4 if this has not
243  // been done already.
244  if (me1a && digi.getStrip() <= 48)
245  {
246  CSCComparatorDigi digi_corr(64+digi.getStrip(),
247  digi.getComparator(),
248  digi.getTimeBinWord());
249  cscData.add(digi_corr, cscDetId); // This version does ME11 strips swapping
250  // cscData.add(digi_corr, cscDetId.layer()); // This one doesn't
251  }
252  else
253  {
254  cscData.add(digi, cscDetId); // This version does ME11 strips swapping
255  // cscData.add(digi, cscDetId.layer()); // This one doesn't
256  }
257  }
258  else
259  {
260  // Move ME1/A comparators from CFEB=0 to CFEB=4 if this has not
261  // been done already.
262  if (me1a && digi.getStrip() <= 16)
263  {
264  CSCComparatorDigi digi_corr(64+digi.getStrip(),
265  digi.getComparator(),
266  digi.getTimeBinWord());
267  cscData.add(digi_corr, cscDetId.layer());
268  }
269  else
270  {
271  cscData.add(digi, cscDetId.layer());
272  }
273  }
274  }
275  }
276  }
277 }
278 
280 {
281  for (CSCALCTDigiCollection::DigiRangeIterator j=alctDigis.begin(); j!=alctDigis.end(); ++j)
282  {
283  CSCDetId cscDetId=(*j).first;
284  CSCEventData & cscData = findEventData(cscDetId);
285 
286  cscData.add(std::vector<CSCALCTDigi>((*j).second.first, (*j).second.second));
287  }
288 }
289 
291 {
292  for (CSCCLCTDigiCollection::DigiRangeIterator j=clctDigis.begin(); j!=clctDigis.end(); ++j)
293  {
294  CSCDetId cscDetId=(*j).first;
295  CSCEventData & cscData = findEventData(cscDetId);
296 
297  cscData.add(std::vector<CSCCLCTDigi>((*j).second.first, (*j).second.second));
298  }
299 }
300 
302 {
303  for (CSCCorrelatedLCTDigiCollection::DigiRangeIterator j=corrLCTDigis.begin(); j!=corrLCTDigis.end(); ++j)
304  {
305  CSCDetId cscDetId=(*j).first;
306  CSCEventData & cscData = findEventData(cscDetId);
307 
308  cscData.add(std::vector<CSCCorrelatedLCTDigi>((*j).second.first, (*j).second.second));
309  }
310 
311 }
312 
313 
315  const CSCWireDigiCollection& wireDigis,
316  const CSCComparatorDigiCollection& comparatorDigis,
317  const CSCALCTDigiCollection& alctDigis,
318  const CSCCLCTDigiCollection& clctDigis,
319  const CSCCLCTPreTriggerCollection & preTriggers,
320  const CSCCorrelatedLCTDigiCollection& correlatedLCTDigis,
321  FEDRawDataCollection& fed_buffers,
322  const CSCChamberMap* mapping,
323  Event & e, uint16_t format_version, bool use_pre_triggers)
324 {
325 
326  theFormatVersion = format_version;
327  usePreTriggers = use_pre_triggers;
328  //bits of code from ORCA/Muon/METBFormatter - thanks, Rick:)!
329 
330  //get fed object from fed_buffers
331  // make a map from the index of a chamber to the event data from it
332  beginEvent(mapping);
333  add(stripDigis, preTriggers);
334  add(wireDigis, alctDigis);
335  add(comparatorDigis, clctDigis);
336  add(correlatedLCTDigis);
337 
338  int l1a=e.id().event(); //need to add increments or get it from lct digis
339  int bx = l1a;//same as above
340  //int startingFED = FEDNumbering::MINCSCFEDID;
341 
342  if (theFormatVersion == 2005)
343  {
344  std::map<int, CSCDCCEventData> dccMap;
345  for (int idcc=FEDNumbering::MINCSCFEDID;
346  idcc<=FEDNumbering::MAXCSCFEDID; ++idcc)
347  {
348  //idcc goes from startingFed to startingFED+7
349  // @@ if ReadoutMapping changes, this'll have to change
350  // DCCs 1,2,4,5 have 5 DDUs. Otherwise, 4
351  //int nDDUs = (idcc < 2) || (idcc ==4) || (idcc ==5)
352  // ? 5 : 4;
353  //@@ WARNING some DCCs only have 4 DDUs, but I'm giving them all 5, for now
354  int nDDUs = 5;
355  dccMap.insert(std::pair<int, CSCDCCEventData>(idcc, CSCDCCEventData(idcc, nDDUs, bx, l1a) ) );
356  }
357 
358  for (int idcc=FEDNumbering::MINCSCFEDID;
359  idcc<=FEDNumbering::MAXCSCFEDID; ++idcc)
360  {
361 
362  // for every chamber with data, add to a DDU in this DCC Event
363  for (map<CSCDetId, CSCEventData>::iterator chamberItr = theChamberDataMap.begin();
364  chamberItr != theChamberDataMap.end(); ++chamberItr)
365  {
366  int indexDCC = mapping->slink(chamberItr->first);
367  if (indexDCC == idcc)
368  {
369  //FIXME (What does this mean? Is something wrong?)
370  std::map<int, CSCDCCEventData>::iterator dccMapItr = dccMap.find(indexDCC);
371  if (dccMapItr == dccMap.end())
372  {
373  throw cms::Exception("CSCDigiToRaw") << "Bad DCC number:" << indexDCC;
374  }
375  // get id's based on ChamberId from mapping
376 
377  int dduId = mapping->ddu(chamberItr->first);
378  int dduSlot = mapping->dduSlot(chamberItr->first);
379  int dduInput = mapping->dduInput(chamberItr->first);
380  int dmbId = mapping->dmb(chamberItr->first);
381  dccMapItr->second.addChamber(chamberItr->second, dduId, dduSlot, dduInput, dmbId, theFormatVersion);
382  }
383  }
384  }
385 
386  // FIXME: FEDRawData size set to 2*64 to add FED header and trailer
387  for (std::map<int, CSCDCCEventData>::iterator dccMapItr = dccMap.begin();
388  dccMapItr != dccMap.end(); ++dccMapItr)
389  {
390  boost::dynamic_bitset<> dccBits = dccMapItr->second.pack();
391  FEDRawData & fedRawData = fed_buffers.FEDData(dccMapItr->first);
392  fedRawData.resize(dccBits.size());
393  //fill data with dccEvent
394  bitset_utilities::bitsetToChar(dccBits, fedRawData.data());
395  FEDTrailer cscFEDTrailer(fedRawData.data()+(fedRawData.size()-8));
396  cscFEDTrailer.set(fedRawData.data()+(fedRawData.size()-8),
397  fedRawData.size()/8,
398  evf::compute_crc(fedRawData.data(),fedRawData.size()), 0, 0);
399  }
400 
401  }
402  else if (theFormatVersion == 2013)
403  {
404 
405  std::map<int, CSCDDUEventData> dduMap;
406  unsigned int ddu_fmt_version = 0x7;
407  const unsigned postLS1_map [] = { 841, 842, 843, 844, 845, 846, 847, 848, 849,
408  831, 832, 833, 834, 835, 836, 837, 838, 839,
409  861, 862, 863, 864, 865, 866, 867, 868, 869,
410  851, 852, 853, 854, 855, 856, 857, 858, 859 };
411 
412 
414  for (unsigned int i = 0; i < 36; i++)
415  {
416  unsigned int iddu = postLS1_map[i];
417  // make a new one
418  CSCDDUHeader newDDUHeader(bx,l1a, iddu, ddu_fmt_version);
419 
420  dduMap.insert(std::pair<int, CSCDDUEventData>(iddu, CSCDDUEventData(newDDUHeader) ) );
421  }
422 
423 
425  for (unsigned int i = 0; i < 36; i++)
426  {
427  unsigned int iddu = postLS1_map[i];
428  // for every chamber with data, add to a DDU in this DCC Event
429  for (map<CSCDetId, CSCEventData>::iterator chamberItr = theChamberDataMap.begin();
430  chamberItr != theChamberDataMap.end(); ++chamberItr)
431  {
432  unsigned int indexDDU = mapping->slink(chamberItr->first);
433 
435  if ( (indexDDU >= FEDNumbering::MINCSCFEDID) && (indexDDU <= FEDNumbering::MAXCSCFEDID) ) // Still preLS1 DCC FEDs mapping
436  {
437 
438  int dduID = mapping->ddu(chamberItr->first); // try to switch to DDU ID mapping
439  if ((dduID >= FEDNumbering::MINCSCDDUFEDID) && (dduID <= FEDNumbering::MAXCSCDDUFEDID) ) // DDU ID is in expectedi post-LS1 FED ID range
440  {
441  indexDDU = dduID;
442  }
443  else // Messy
444  {
445  // Lets try to change pre-LS1 1-36 ID range to post-LS1 MINCSCDDUFEDID - MAXCSCDDUFEDID range
446  dduID &= 0xFF;
447  if ((dduID <= 36) && (dduID > 0)) indexDDU = postLS1_map[dduID -1]; // indexDDU = FEDNumbering::MINCSCDDUFEDID + dduID-1;
448  }
449 
450  }
451 
452  if (indexDDU == iddu)
453  {
454  std::map<int, CSCDDUEventData>::iterator dduMapItr = dduMap.find(indexDDU);
455  if (dduMapItr == dduMap.end())
456  {
457  throw cms::Exception("CSCDigiToRaw") << "Bad DDU number:" << indexDDU;
458  }
459  // get id's based on ChamberId from mapping
460 
461  int dduInput = mapping->dduInput(chamberItr->first);
462  int dmbId = mapping->dmb(chamberItr->first);
463  // int crateId = mapping->crate(chamberItr->first);
464  dduMapItr->second.add( chamberItr->second, dmbId, dduInput, theFormatVersion );
465  }
466  }
467  }
468 
469  // FIXME: FEDRawData size set to 2*64 to add FED header and trailer
470  for (std::map<int, CSCDDUEventData>::iterator dduMapItr = dduMap.begin();
471  dduMapItr != dduMap.end(); ++dduMapItr)
472  {
473  boost::dynamic_bitset<> dduBits = dduMapItr->second.pack();
474  FEDRawData & fedRawData = fed_buffers.FEDData(dduMapItr->first);
475  fedRawData.resize(dduBits.size()/8);
476  //fill data with dduEvent
477  bitset_utilities::bitsetToChar(dduBits, fedRawData.data());
478  FEDTrailer cscFEDTrailer(fedRawData.data()+(fedRawData.size()-8));
479  cscFEDTrailer.set(fedRawData.data()+(fedRawData.size()-8),
480  fedRawData.size()/8,
481  evf::compute_crc(fedRawData.data(),fedRawData.size()), 0, 0);
482 
483  }
484  }
485 }
486 
487 
01/20/05 A.Tumanov
int chamber() const
Definition: CSCDetId.h:81
EventNumber_t event() const
Definition: EventID.h:44
void add(const CSCStripDigiCollection &stripDigis, const CSCCLCTPreTriggerCollection &preTriggers)
int i
Definition: DBlmapReader.cc:9
uint16_t theFormatVersion
Definition: CSCDigiToRaw.h:67
int getStrip() const
Get the strip number.
bool usePreTriggers
Definition: CSCDigiToRaw.h:68
int ddu(const CSCDetId &) const
ddu id for given DetId
int clctWindowMin_
Definition: CSCDigiToRaw.h:63
int dmb(const CSCDetId &) const
dmb id for given DetId
void setdmbVersion(unsigned int version)
Definition: CSCDMBHeader.h:43
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:26
int getComparator() const
Get Comparator readings.
int layer() const
Definition: CSCDetId.h:74
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
int getStrip() const
Definition: CSCStripDigi.h:51
int endcap() const
Definition: CSCDetId.h:106
int alctWindowMax_
Definition: CSCDigiToRaw.h:62
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
void setStrip(int istrip)
Definition: CSCStripDigi.h:75
const CSCChamberMap * theElectronicsMap
Definition: CSCDigiToRaw.h:60
void resize(size_t newsize)
Definition: FEDRawData.cc:32
int slink(const CSCDetId &) const
slink id for given DetId
unsigned short compute_crc(unsigned char *buffer, unsigned int bufSize)
Definition: CRC16.h:67
tuple result
Definition: query.py:137
int alctWindowMin_
Definition: CSCDigiToRaw.h:61
int dduInput(const CSCDetId &) const
ddu input for given DetId
int j
Definition: DBlmapReader.cc:9
CSCDetId chamberId() const
Definition: CSCDetId.h:66
CSCDetId chamberID(const CSCDetId &cscDetId)
takes layer ID, converts to chamber ID, switching ME1A to ME11
Definition: CSCDigiToRaw.cc:31
CSCDigiToRaw(const edm::ParameterSet &pset)
Constructor.
Definition: CSCDigiToRaw.cc:89
int getTimeBinWord() const
Return the word with each bit corresponding to a time bin.
int clctWindowMax_
Definition: CSCDigiToRaw.h:64
void beginEvent(const CSCChamberMap *electronicsMap)
void setCrateAddress(int crate, int dmbId)
Definition: CSCDMBHeader.h:41
int preTriggerWindowMin_
Definition: CSCDigiToRaw.h:65
unsigned short iChamberType()
Definition: CSCDetId.h:120
int preTriggerWindowMax_
Definition: CSCDigiToRaw.h:66
int ring() const
Definition: CSCDetId.h:88
PixelRecoRange< float > Range
std::map< CSCDetId, CSCEventData > theChamberDataMap
Definition: CSCDigiToRaw.h:59
void bitsetToChar(const boost::dynamic_bitset<> &bs, unsigned char *result)
this method takes bitset obj and returns char * array
int crate(const CSCDetId &) const
Interface required use in digi-to-raw.
std::vector< CSCCLCTPreTrigger >::const_iterator const_iterator
int dduSlot(const CSCDetId &) const
ddu slot for given DetId
string bx
Definition: rpc-layouts.py:13
CSCEventData & findEventData(const CSCDetId &cscDetId)
pick out the correct data object for this chamber
void add(const CSCStripDigi &, int layer)
routines to add digis to the data
edm::EventID id() const
Definition: EventBase.h:56
EcalElectronicsMapping const * electronicsMap(0)
bool accept(const CSCDetId &cscId, const LCTCollection &lcts, int bxMin, int bxMax)
Definition: CSCDigiToRaw.cc:42
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
int station() const
Definition: CSCDetId.h:99
std::pair< const_iterator, const_iterator > Range
const CSCDMBHeader * dmbHeader() const
the DAQ motherboard header. A good place for event and chamber info
Definition: CSCEventData.h:93
void createFedBuffers(const CSCStripDigiCollection &stripDigis, const CSCWireDigiCollection &wireDigis, const CSCComparatorDigiCollection &comparatorDigis, const CSCALCTDigiCollection &alctDigis, const CSCCLCTDigiCollection &clctDigis, const CSCCLCTPreTriggerCollection &preTriggers, const CSCCorrelatedLCTDigiCollection &correlatedLCTDigis, FEDRawDataCollection &fed_buffers, const CSCChamberMap *theMapping, edm::Event &e, uint16_t theFormatVersion=2005, bool usePreTriggers=true)
Take a vector of digis and fill the FEDRawDataCollection.