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&);
44  ~HcalDigiToRawuHTR() override;
45 
46  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
47 
48  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
49 
50 private:
53 
59 
60  bool premix_;
61 };
62 
64  _verbosity(iConfig.getUntrackedParameter<int>("Verbosity", 0)),
65  electronicsMapLabel_(iConfig.getParameter<std::string>("ElectronicsMap")),
66  tok_QIE10DigiCollection_(consumes<HcalDataFrameContainer<QIE10DataFrame> >(iConfig.getParameter<edm::InputTag>("QIE10"))),
67  tok_QIE11DigiCollection_(consumes<HcalDataFrameContainer<QIE11DataFrame> >(iConfig.getParameter<edm::InputTag>("QIE11"))),
68  tok_HBHEDigiCollection_(consumes<HBHEDigiCollection >(iConfig.getParameter<edm::InputTag>("HBHEqie8"))),
69  tok_HFDigiCollection_(consumes<HFDigiCollection>(iConfig.getParameter<edm::InputTag>("HFqie8"))),
70  tok_TPDigiCollection_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("TP"))),
71  premix_(iConfig.getParameter<bool>("premix"))
72 {
73  produces<FEDRawDataCollection>("");
74 }
75 
77 
79 
80  using namespace edm;
81 
84  const HcalElectronicsMap* readoutMap = item.product();
85 
86  //collection to be inserted into event
87  std::unique_ptr<FEDRawDataCollection> fed_buffers(new FEDRawDataCollection());
88 
89  //
90  // Extracting All the Collections containing useful Info
91  edm::Handle<QIE10DigiCollection> qie10DigiCollection;
92  edm::Handle<QIE11DigiCollection> qie11DigiCollection;
93  edm::Handle<HBHEDigiCollection> hbheDigiCollection;
94  edm::Handle<HFDigiCollection> hfDigiCollection;
96  iEvent.getByToken(tok_QIE10DigiCollection_,qie10DigiCollection);
97  iEvent.getByToken(tok_QIE11DigiCollection_,qie11DigiCollection);
98  iEvent.getByToken(tok_HBHEDigiCollection_,hbheDigiCollection);
99  iEvent.getByToken(tok_HFDigiCollection_,hfDigiCollection);
100  iEvent.getByToken(tok_TPDigiCollection_,tpDigiCollection);
101 
102  // first argument is the fedid (minFEDID+crateId)
103  map<int,unique_ptr<HCalFED> > fedMap;
104 
105  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
106  // QIE10 precision data
107  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
108  UHTRpacker uhtrs;
109  // loop over each digi and allocate memory for each
110  if( qie10DigiCollection.isValid() ){
111  const QIE10DigiCollection& qie10dc=*(qie10DigiCollection);
112  for (unsigned int j=0; j < qie10dc.size(); j++){
113  QIE10DataFrame qiedf = static_cast<QIE10DataFrame>(qie10dc[j]);
114  DetId detid = qiedf.detid();
115  HcalElectronicsId eid(readoutMap->lookup(detid));
116  int crateId = eid.crateId();
117  int slotId = eid.slot();
118  int uhtrIndex = ((slotId&0xF)<<8) | (crateId&0xFF);
119  int presamples = qiedf.presamples();
120 
121  /* Defining a custom index that will encode only
122  the information about the crate and slot of a
123  given channel: crate: bits 0-7
124  slot: bits 8-12 */
125 
126  if( ! uhtrs.exist( uhtrIndex ) ){
127  uhtrs.newUHTR( uhtrIndex , presamples );
128  }
129  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,_verbosity);
130  }
131  }
132  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
133  // QIE11 precision data
134  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
135  //UHTRpacker uhtrs;
136  // loop over each digi and allocate memory for each
137  if( qie11DigiCollection.isValid() ){
138  const QIE11DigiCollection& qie11dc=*(qie11DigiCollection);
139  for (unsigned int j=0; j < qie11dc.size(); j++){
140  QIE11DataFrame qiedf = static_cast<QIE11DataFrame>(qie11dc[j]);
141  DetId detid = qiedf.detid();
142  HcalElectronicsId eid(readoutMap->lookup(detid));
143  int crateId = eid.crateId();
144  int slotId = eid.slot();
145  int uhtrIndex = ((slotId&0xF)<<8) | (crateId&0xFF);
146  int presamples = qiedf.presamples();
147 
148  if( ! uhtrs.exist(uhtrIndex) ){
149  uhtrs.newUHTR( uhtrIndex , presamples );
150  }
151  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,_verbosity);
152  }
153  }
154  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
155  // HF (QIE8) precision data
156  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
157  // loop over each digi and allocate memory for each
158  if(hfDigiCollection.isValid()){
159  const HFDigiCollection& qie8hfdc=*(hfDigiCollection);
160  for(HFDigiCollection::const_iterator qiedf=qie8hfdc.begin();qiedf!=qie8hfdc.end();qiedf++){
161  DetId detid = qiedf->id();
162 
163  HcalElectronicsId eid(readoutMap->lookup(detid));
164  int crateId = eid.crateId();
165  int slotId = eid.slot();
166  int uhtrIndex = (crateId&0xFF) | ((slotId&0xF)<<8) ;
167  int presamples = qiedf->presamples();
168 
169  if( ! uhtrs.exist(uhtrIndex) ){
170  uhtrs.newUHTR( uhtrIndex , presamples );
171  }
172  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,premix_,_verbosity);
173  }
174  }
175  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
176  // HBHE (QIE8) precision data
177  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
178  // loop over each digi and allocate memory for each
179  if(hbheDigiCollection.isValid()){
180  const HBHEDigiCollection& qie8hbhedc=*(hbheDigiCollection);
181  for(HBHEDigiCollection::const_iterator qiedf=qie8hbhedc.begin();qiedf!=qie8hbhedc.end();qiedf++){
182  DetId detid = qiedf->id();
183 
184  HcalElectronicsId eid(readoutMap->lookup(detid));
185  int crateId = eid.crateId();
186  int slotId = eid.slot();
187  int uhtrIndex = (crateId&0xFF) | ((slotId&0xF)<<8) ;
188  int presamples = qiedf->presamples();
189 
190  if( ! uhtrs.exist(uhtrIndex) ){
191  uhtrs.newUHTR( uhtrIndex , presamples );
192  }
193  uhtrs.addChannel(uhtrIndex,qiedf,readoutMap,premix_,_verbosity);
194  }
195  }
196  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
197  // TP data
198  // - - - - - - - - - - - - - - - - - - - - - - - - - - -
199  // loop over each digi and allocate memory for each
200  if(tpDigiCollection.isValid()){
202  for(HcalTrigPrimDigiCollection::const_iterator qiedf=qietpdc.begin();qiedf!=qietpdc.end();qiedf++){
203  DetId detid = qiedf->id();
204  HcalElectronicsId eid(readoutMap->lookupTrigger(detid));
205 
206  int crateId = eid.crateId();
207  int slotId = eid.slot();
208  int uhtrIndex = (crateId&0xFF) | ((slotId&0xF)<<8);
209  int ilink = eid.fiberIndex();
210  int itower = eid.fiberChanId();
211  int channelid = (itower&0xF) | ((ilink&0xF)<<4);
212  int presamples = qiedf->presamples();
213 
214  if( ! uhtrs.exist(uhtrIndex) ){
215  uhtrs.newUHTR( uhtrIndex , presamples );
216  }
217  uhtrs.addChannel(uhtrIndex,qiedf,channelid,_verbosity);
218  }
219  }
220  // -----------------------------------------------------
221  // -----------------------------------------------------
222  // loop over each uHTR and format data
223  // -----------------------------------------------------
224  // -----------------------------------------------------
225  // loop over each uHTR and format data
226  int idxuhtr =-1;
227  for( UHTRpacker::UHTRMap::iterator uhtr = uhtrs.uhtrs.begin() ; uhtr != uhtrs.uhtrs.end() ; ++uhtr){
228 
229  idxuhtr ++;
230 
231  uint64_t crateId = (uhtr->first)&0xFF;
232  uint64_t slotId = (uhtr->first&0xF00)>>8;
233 
234  uhtrs.finalizeHeadTail(&(uhtr->second),_verbosity);
235  int fedId = FEDNumbering::MINHCALuTCAFEDID + crateId;
236  if( fedMap.find(fedId) == fedMap.end() ){
237  /* QUESTION: where should the orbit number come from? */
238  fedMap[fedId] = std::unique_ptr<HCalFED>(new HCalFED(fedId,iEvent.id().event(),iEvent.orbitNumber(),iEvent.bunchCrossing()));
239  }
240  fedMap[fedId]->addUHTR(uhtr->second,crateId,slotId);
241  }// end loop over uhtr containers
242 
243  /* ------------------------------------------------------
244  ------------------------------------------------------
245  putting together the FEDRawDataCollection
246  ------------------------------------------------------
247  ------------------------------------------------------ */
248  for( map<int,unique_ptr<HCalFED> >::iterator fed = fedMap.begin() ; fed != fedMap.end() ; ++fed ){
249 
250  int fedId = fed->first;
251 
252  auto & rawData = fed_buffers->FEDData(fedId);
253  fed->second->formatFEDdata(rawData);
254 
255  FEDHeader hcalFEDHeader(rawData.data());
256  hcalFEDHeader.set(rawData.data(), 1, iEvent.id().event(), iEvent.bunchCrossing(), fedId);
257  FEDTrailer hcalFEDTrailer(rawData.data()+(rawData.size()-8));
258  hcalFEDTrailer.set(rawData.data()+(rawData.size()-8), rawData.size()/8, evf::compute_crc(rawData.data(),rawData.size()), 0, 0);
259 
260  }// end loop over FEDs with data
261 
262  iEvent.put(std::move(fed_buffers));
263 
264 }
265 
266 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
268  //The following says we do not know what parameters are allowed so do no validation
269  // Please change this to state exactly what you do use, even if it is no parameters
271  desc.addUntracked<int>("Verbosity", 0);
272  desc.add<std::string>("ElectronicsMap", "");
273  desc.add<edm::InputTag>("QIE10", edm::InputTag("simHcalDigis", "HFQIE10DigiCollection"));
274  desc.add<edm::InputTag>("QIE11", edm::InputTag("simHcalDigis", "HBHEQIE11DigiCollection"));
275  desc.add<edm::InputTag>("HBHEqie8", edm::InputTag("simHcalDigis"));
276  desc.add<edm::InputTag>("HFqie8", edm::InputTag("simHcalDigis"));
277  desc.add<edm::InputTag>("TP", edm::InputTag("simHcalTriggerPrimitiveDigis"));
278  desc.add<bool>("premix", false);
279  descriptions.add("hcalDigiToRawuHTR",desc);
280  descriptions.addDefault(desc);
281 }
282 
283 //define this as a plug-in
DetId detid() const
Get the detector id.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
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:136
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void finalizeHeadTail(uhtrData *uhtr, bool verbosity)
Definition: PackerHelp.h:476
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int presamples() const
for backward compatibility
uhtrData * newUHTR(int uhtrIndex, int ps=0, int orn=0, int bcn=0, uint64_t evt=0)
Definition: PackerHelp.h:435
std::vector< HFDataFrame >::const_iterator const_iterator
int bunchCrossing() const
Definition: EventBase.h:66
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
static void set(unsigned char *trailer, uint32_t lenght, uint16_t crc, uint8_t evt_stat, uint8_t tts, bool moreTrailers=false)
Set all fields in the trailer.
Definition: FEDTrailer.cc:60
edm::EDGetTokenT< HcalDataFrameContainer< QIE10DataFrame > > tok_QIE10DigiCollection_
std::string electronicsMapLabel_
unsigned short compute_crc(unsigned char *buffer, unsigned int bufSize)
Definition: CRC16.h:67
edm::EDGetTokenT< HFDigiCollection > tok_HFDigiCollection_
~HcalDigiToRawuHTR() override
DetId detid() const
Get the detector id.
int orbitNumber() const
Definition: EventBase.h:67
HcalDigiToRawuHTR(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
const_iterator end() const
UHTRMap uhtrs
Definition: PackerHelp.h:314
Definition: DetId.h:18
unsigned long long uint64_t
Definition: Time.h:15
const T & get() const
Definition: EventSetup.h:58
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void addChannel(int uhtrIndex, edm::SortedCollection< HFDataFrame >::const_iterator &qiedf, const HcalElectronicsMap *readoutMap, bool premix, int verbosity=0)
Definition: PackerHelp.h:499
int presamples() const
for backward compatibility
bool exist(int uhtrIndex)
Definition: PackerHelp.h:345
edm::EventID id() const
Definition: EventBase.h:60
HLT enums.
const DetId lookupTrigger(HcalElectronicsId fId) const
brief lookup the trigger logical detid associated with the given electronics id
static void set(unsigned char *header, uint8_t triggerType, uint32_t lvl1ID, uint16_t bxID, uint16_t sourceID, uint8_t version=0, bool moreHeaders=false)
Set all fields in the header.
Definition: FEDHeader.cc:47
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_