CMS 3D CMS Logo

HcalDigiToRawuHTR.cc
Go to the documentation of this file.
6 
9 
12 
14 
22 
25 
26 #include "PackerHelp.h"
27 
28 #include <iostream>
29 #include <fstream>
30 #include <sstream>
31 #include <string>
32 
33 /* QUESTION: what about dual FED readout? */
34 /* QUESTION: what do I do if the number of 16-bit words
35  are not divisible by 4? -- these need to
36  fit into the 64-bit words of the FEDRawDataFormat */
37 
38 
39 using namespace std;
40 
42 public:
43  explicit HcalDigiToRawuHTR(const edm::ParameterSet&);
45 
46  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
47 
48 private:
49  virtual void produce(edm::Event&, const edm::EventSetup&) override;
50  void getData(const edm::Event&, const edm::EventSetup&);
51 
54 
65 
66  edm::InputTag qie10Tag_, qie11Tag_, hbheqie8Tag_, hfqie8Tag_, trigTag_;
67 };
68 
70  _verbosity(iConfig.getUntrackedParameter<int>("Verbosity", 0)),
71  electronicsMapLabel_(iConfig.getParameter<std::string>("ElectronicsMap")),
72  qie10Tag_(iConfig.getParameter<edm::InputTag>("QIE10")),
73  qie11Tag_(iConfig.getParameter<edm::InputTag>("QIE11")),
74  hbheqie8Tag_(iConfig.getParameter<edm::InputTag>("HBHEqie8")),
75  hfqie8Tag_(iConfig.getParameter<edm::InputTag>("HFqie8")),
76  trigTag_(iConfig.getParameter<edm::InputTag>("TP"))
77 {
78  produces<FEDRawDataCollection>("");
79  tok_QIE10DigiCollection_ = consumes<HcalDataFrameContainer<QIE10DataFrame> >(qie10Tag_);
80  tok_QIE11DigiCollection_ = consumes<HcalDataFrameContainer<QIE11DataFrame> >(qie11Tag_);
81  tok_HBHEDigiCollection_ = consumes<HBHEDigiCollection >(hbheqie8Tag_);
82  tok_HFDigiCollection_ = consumes<HFDigiCollection>(hfqie8Tag_);
83  tok_TPDigiCollection_ = consumes<HcalTrigPrimDigiCollection>(trigTag_);
84 }
85 
87 
89 
90  using namespace edm;
91 
94  const HcalElectronicsMap* readoutMap = item.product();
95 
96  //collection to be inserted into event
97  std::unique_ptr<FEDRawDataCollection> fed_buffers(new FEDRawDataCollection());
98 
99  //
100  // Extracting All the Collections containing useful Info
106 
107  // first argument is the fedid (minFEDID+crateId)
108  map<int,unique_ptr<HCalFED> > fedMap;
109 
110  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
111  // QIE10 precision data
112  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
113  UHTRpacker uhtrs;
114  // loop over each digi and allocate memory for each
115  if( qie10DigiCollection.isValid() ){
116  const QIE10DigiCollection& qie10dc=*(qie10DigiCollection);
117  for (unsigned int j=0; j < qie10dc.size(); j++){
118  QIE10DataFrame qiedf = static_cast<QIE10DataFrame>(qie10dc[j]);
119  DetId detid = qiedf.detid();
120  HcalElectronicsId eid(readoutMap->lookup(detid));
121  int crateId = eid.crateId();
122  int slotId = eid.slot();
123  int uhtrIndex = ((slotId&0xF)<<8) | (crateId&0xFF);
124  int presamples = qiedf.presamples();
125 
126  /* Defining a custom index that will encode only
127  the information about the crate and slot of a
128  given channel: crate: bits 0-7
129  slot: bits 8-12 */
130 
131  if( ! uhtrs.exist( uhtrIndex ) ){
132  uhtrs.newUHTR( uhtrIndex , presamples );
133  }
134  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,_verbosity);
135  }
136  }
137  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
138  // QIE11 precision data
139  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
140  //UHTRpacker uhtrs;
141  // loop over each digi and allocate memory for each
142  if( qie11DigiCollection.isValid() ){
143  const QIE11DigiCollection& qie11dc=*(qie11DigiCollection);
144  for (unsigned int j=0; j < qie11dc.size(); j++){
145  QIE11DataFrame qiedf = static_cast<QIE11DataFrame>(qie11dc[j]);
146  DetId detid = qiedf.detid();
147  HcalElectronicsId eid(readoutMap->lookup(detid));
148  int crateId = eid.crateId();
149  int slotId = eid.slot();
150  int uhtrIndex = ((slotId&0xF)<<8) | (crateId&0xFF);
151  int presamples = qiedf.presamples();
152 
153  if( ! uhtrs.exist(uhtrIndex) ){
154  uhtrs.newUHTR( uhtrIndex , presamples );
155  }
156  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,_verbosity);
157  }
158  }
159  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
160  // HF (QIE8) precision data
161  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
162  // loop over each digi and allocate memory for each
163  if(hfDigiCollection.isValid()){
164  const HFDigiCollection& qie8hfdc=*(hfDigiCollection);
165  for(HFDigiCollection::const_iterator qiedf=qie8hfdc.begin();qiedf!=qie8hfdc.end();qiedf++){
166  DetId detid = qiedf->id();
167 
168  HcalElectronicsId eid(readoutMap->lookup(detid));
169  int crateId = eid.crateId();
170  int slotId = eid.slot();
171  int uhtrIndex = (crateId&0xFF) | ((slotId&0xF)<<8) ;
172  int presamples = qiedf->presamples();
173 
174  if( ! uhtrs.exist(uhtrIndex) ){
175  uhtrs.newUHTR( uhtrIndex , presamples );
176  }
177  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,_verbosity);
178  }
179  }
180  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
181  // HBHE (QIE8) precision data
182  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
183  // loop over each digi and allocate memory for each
184  if(hbheDigiCollection.isValid()){
185  const HBHEDigiCollection& qie8hbhedc=*(hbheDigiCollection);
186  for(HBHEDigiCollection::const_iterator qiedf=qie8hbhedc.begin();qiedf!=qie8hbhedc.end();qiedf++){
187  DetId detid = qiedf->id();
188 
189  HcalElectronicsId eid(readoutMap->lookup(detid));
190  int crateId = eid.crateId();
191  int slotId = eid.slot();
192  int uhtrIndex = (crateId&0xFF) | ((slotId&0xF)<<8) ;
193  int presamples = qiedf->presamples();
194 
195  if( ! uhtrs.exist(uhtrIndex) ){
196  uhtrs.newUHTR( uhtrIndex , presamples );
197  }
198  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,_verbosity);
199  }
200  }
201  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
202  // TP data
203  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
204  // loop over each digi and allocate memory for each
205  if(tpDigiCollection.isValid()){
207  for(HcalTrigPrimDigiCollection::const_iterator qiedf=qietpdc.begin();qiedf!=qietpdc.end();qiedf++){
208  DetId detid = qiedf->id();
209  HcalElectronicsId eid(readoutMap->lookupTrigger(detid));
210 
211  int crateId = eid.crateId();
212  int slotId = eid.slot();
213  int uhtrIndex = (crateId&0xFF) | ((slotId&0xF)<<8);
214  int ilink = eid.fiberIndex();
215  int itower = eid.fiberChanId();
216  int channelid = (itower&0xF) | ((ilink&0xF)<<4);
217  int presamples = qiedf->presamples();
218 
219  if( ! uhtrs.exist(uhtrIndex) ){
220  uhtrs.newUHTR( uhtrIndex , presamples );
221  }
222  uhtrs.addChannel(uhtrIndex,qiedf,channelid,_verbosity);
223  }
224  }
225  // -----------------------------------------------------
226  // -----------------------------------------------------
227  // loop over each uHTR and format data
228  // -----------------------------------------------------
229  // -----------------------------------------------------
230  // loop over each uHTR and format data
231  int idxuhtr =-1;
232  for( UHTRpacker::UHTRMap::iterator uhtr = uhtrs.uhtrs.begin() ; uhtr != uhtrs.uhtrs.end() ; ++uhtr){
233 
234  idxuhtr ++;
235 
236  uint64_t crateId = (uhtr->first)&0xFF;
237  uint64_t slotId = (uhtr->first&0xF00)>>8;
238 
239  uhtrs.finalizeHeadTail(&(uhtr->second),_verbosity);
240  int fedId = FEDNumbering::MINHCALuTCAFEDID + crateId;
241  if( fedMap.find(fedId) == fedMap.end() ){
242  /* QUESTION: where should the orbit number come from? */
243  fedMap[fedId] = std::unique_ptr<HCalFED>(new HCalFED(fedId,iEvent.id().event(),iEvent.orbitNumber(),iEvent.bunchCrossing()));
244  }
245  fedMap[fedId]->addUHTR(uhtr->second,crateId,slotId);
246  }// end loop over uhtr containers
247 
248  /* ------------------------------------------------------
249  ------------------------------------------------------
250  putting together the FEDRawDataCollection
251  ------------------------------------------------------
252  ------------------------------------------------------ */
253  for( map<int,unique_ptr<HCalFED> >::iterator fed = fedMap.begin() ; fed != fedMap.end() ; ++fed ){
254 
255  int fedId = fed->first;
256 
257  auto & rawData = fed_buffers->FEDData(fedId);
258  fed->second->formatFEDdata(rawData);
259 
260  FEDHeader hcalFEDHeader(rawData.data());
261  hcalFEDHeader.set(rawData.data(), 1, iEvent.id().event(), iEvent.bunchCrossing(), fedId);
262  FEDTrailer hcalFEDTrailer(rawData.data()+(rawData.size()-8));
263  hcalFEDTrailer.set(rawData.data()+(rawData.size()-8), rawData.size()/8, evf::compute_crc(rawData.data(),rawData.size()), 0, 0);
264 
265  }// end loop over FEDs with data
266 
267  iEvent.put(std::move(fed_buffers));
268 
269 }
270 
271 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
273  //The following says we do not know what parameters are allowed so do no validation
274  // Please change this to state exactly what you do use, even if it is no parameters
276  desc.addUntracked<int>("Verbosity", 0);
277  desc.add<std::string>("ElectronicsMap", "");
278  desc.add<edm::InputTag>("QIE10", edm::InputTag("simHcalDigis", "HFQIE10DigiCollection"));
279  desc.add<edm::InputTag>("QIE11", edm::InputTag("simHcalDigis", "HBHEQIE11DigiCollection"));
280  desc.add<edm::InputTag>("HBHEqie8", edm::InputTag("simHcalDigis"));
281  desc.add<edm::InputTag>("HFqie8", edm::InputTag("simHcalDigis"));
282  desc.add<edm::InputTag>("TP", edm::InputTag("simHcalTriggerPrimitiveDigis"));
283  descriptions.add("hcalDigiToRawuHTR",desc);
284  descriptions.addDefault(desc);
285 }
286 
287 //define this as a plug-in
DetId detid() const
Get the detector id.
EventNumber_t event() const
Definition: EventID.h:41
edm::EDGetTokenT< HBHEDigiCollection > tok_HBHEDigiCollection_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void finalizeHeadTail(uhtrData *uhtr, bool verbosity)
Definition: PackerHelp.h:439
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int presamples() const
for backward compatibility
edm::InputTag hfqie8Tag_
uhtrData * newUHTR(int uhtrIndex, int ps=0, int orn=0, int bcn=0, uint64_t evt=0)
Definition: PackerHelp.h:398
virtual void produce(edm::Event &, const edm::EventSetup &) override
std::vector< HFDataFrame >::const_iterator const_iterator
int bunchCrossing() const
Definition: EventBase.h:64
edm::InputTag qie11Tag_
def getData(doc, options, dataset, site)
edm::InputTag qie10Tag_
static void set(unsigned char *trailer, int evt_lgth, int crc, int evt_stat, int tts, bool T=false)
Set all fields in the trailer.
Definition: FEDTrailer.cc:42
edm::Handle< HcalTrigPrimDigiCollection > tpDigiCollection
void addChannel(int uhtrIndex, edm::SortedCollection< HFDataFrame >::const_iterator &qiedf, const HcalElectronicsMap *readoutMap, int verbosity=0)
Definition: PackerHelp.h:462
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< HcalDataFrameContainer< QIE10DataFrame > > tok_QIE10DigiCollection_
std::string electronicsMapLabel_
static void set(unsigned char *header, int evt_ty, int lvl1_ID, int bx_ID, int source_ID, int version=0, bool H=false)
Set all fields in the header.
Definition: FEDHeader.cc:40
unsigned short compute_crc(unsigned char *buffer, unsigned int bufSize)
Definition: CRC16.h:67
edm::EDGetTokenT< HFDigiCollection > tok_HFDigiCollection_
DetId detid() const
Get the detector id.
edm::Handle< HFDigiCollection > hfDigiCollection
int orbitNumber() const
Definition: EventBase.h:65
HcalDigiToRawuHTR(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::Handle< QIE10DigiCollection > qie10DigiCollection
edm::InputTag trigTag_
const_iterator end() const
UHTRMap uhtrs
Definition: PackerHelp.h:299
Definition: DetId.h:18
unsigned long long uint64_t
Definition: Time.h:15
edm::InputTag hbheqie8Tag_
const T & get() const
Definition: EventSetup.h:56
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int presamples() const
for backward compatibility
bool exist(int uhtrIndex)
Definition: PackerHelp.h:330
edm::EventID id() const
Definition: EventBase.h:58
HLT enums.
edm::Handle< QIE11DigiCollection > qie11DigiCollection
edm::Handle< HBHEDigiCollection > hbheDigiCollection
const DetId lookupTrigger(HcalElectronicsId fId) const
brief lookup the trigger logical detid associated with the given electronics id
T const * product() const
Definition: ESHandle.h:86
Readout chain identification for Hcal.
edm::EDGetTokenT< HcalTrigPrimDigiCollection > tok_TPDigiCollection_
const DetId lookup(HcalElectronicsId fId) const
lookup the logical detid associated with the given electronics id
def move(src, dest)
Definition: eostools.py:510
const_iterator begin() const
edm::EDGetTokenT< HcalDataFrameContainer< QIE11DataFrame > > tok_QIE11DigiCollection_