CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalRawToDigi.cc
Go to the documentation of this file.
1 using namespace std;
10 #include <iostream>
11 
13  dataTag_(conf.getParameter<edm::InputTag>("InputLabel")),
14  unpacker_(conf.getUntrackedParameter<int>("HcalFirstFED",int(FEDNumbering::MINHCALFEDID)),conf.getParameter<int>("firstSample"),conf.getParameter<int>("lastSample")),
15  filter_(conf.getParameter<bool>("FilterDataQuality"),conf.getParameter<bool>("FilterDataQuality"),
16  false,
17  0, 0,
18  -1),
19  fedUnpackList_(conf.getUntrackedParameter<std::vector<int> >("FEDs", std::vector<int>())),
20  firstFED_(conf.getUntrackedParameter<int>("HcalFirstFED",FEDNumbering::MINHCALFEDID)),
21  unpackCalib_(conf.getUntrackedParameter<bool>("UnpackCalib",false)),
22  unpackZDC_(conf.getUntrackedParameter<bool>("UnpackZDC",false)),
23  unpackTTP_(conf.getUntrackedParameter<bool>("UnpackTTP",false)),
24  silent_(conf.getUntrackedParameter<bool>("silent",true)),
25  complainEmptyData_(conf.getUntrackedParameter<bool>("ComplainEmptyData",false)),
26  expectedOrbitMessageTime_(conf.getUntrackedParameter<int>("ExpectedOrbitMessageTime",-1))
27 {
28  if (fedUnpackList_.empty()) {
30  fedUnpackList_.push_back(i);
31  }
32 
34  std::ostringstream ss;
35  for (unsigned int i=0; i<fedUnpackList_.size(); i++)
36  ss << fedUnpackList_[i] << " ";
37  edm::LogInfo("HCAL") << "HcalRawToDigi will unpack FEDs ( " << ss.str() << ")";
38 
39  // products produced...
40  produces<HBHEDigiCollection>();
41  produces<HFDigiCollection>();
42  produces<HODigiCollection>();
43  produces<HcalTrigPrimDigiCollection>();
44  produces<HOTrigPrimDigiCollection>();
45  produces<HcalUnpackerReport>();
46  if (unpackCalib_)
47  produces<HcalCalibDigiCollection>();
48  if (unpackZDC_)
49  produces<ZDCDigiCollection>();
50  if (unpackTTP_)
51  produces<HcalTTPDigiCollection>();
52 
53  memset(&stats_,0,sizeof(stats_));
54 
55 }
56 
57 // Virtual destructor needed.
59 
60 // Functions that gets called by framework every event
62 {
63  // Step A: Get Inputs
65  e.getByLabel(dataTag_,rawraw);
66  // get the mapping
68  es.get<HcalDbRecord>().get( pSetup );
69  const HcalElectronicsMap* readoutMap=pSetup->getHcalMapping();
70 
71  // Step B: Create empty output : three vectors for three classes...
72  std::vector<HBHEDataFrame> hbhe;
73  std::vector<HODataFrame> ho;
74  std::vector<HFDataFrame> hf;
75  std::vector<HcalTriggerPrimitiveDigi> htp;
76  std::vector<HcalCalibDataFrame> hc;
77  std::vector<ZDCDataFrame> zdc;
78  std::vector<HcalTTPDigi> ttp;
79  std::vector<HOTriggerPrimitiveDigi> hotp;
80  std::auto_ptr<HcalUnpackerReport> report(new HcalUnpackerReport);
81 
82  // Heuristics: use ave+(max-ave)/8
83  if (stats_.max_hbhe>0) hbhe.reserve(stats_.ave_hbhe+(stats_.max_hbhe-stats_.ave_hbhe)/8);
84  if (stats_.max_ho>0) ho.reserve(stats_.ave_ho+(stats_.max_ho-stats_.ave_ho)/8);
85  if (stats_.max_hf>0) hf.reserve(stats_.ave_hf+(stats_.max_hf-stats_.ave_hf)/8);
87  if (stats_.max_tp>0) htp.reserve(stats_.ave_tp+(stats_.max_tp-stats_.ave_tp)/8);
88  if (stats_.max_tpho>0) hotp.reserve(stats_.ave_tpho+(stats_.max_tpho-stats_.ave_tpho)/8);
89 
90  if (unpackZDC_) zdc.reserve(24);
91 
92 
94  colls.hbheCont=&hbhe;
95  colls.hoCont=&ho;
96  colls.hfCont=&hf;
97  colls.tpCont=&htp;
98  colls.tphoCont=&hotp;
99  colls.calibCont=&hc;
100  colls.zdcCont=&zdc;
101  if (unpackTTP_) colls.ttp=&ttp;
102 
103  // Step C: unpack all requested FEDs
104  for (std::vector<int>::const_iterator i=fedUnpackList_.begin(); i!=fedUnpackList_.end(); i++) {
105  const FEDRawData& fed = rawraw->FEDData(*i);
106  if (fed.size()==0) {
107  if (complainEmptyData_) {
108  if (!silent_) edm::LogWarning("EmptyData") << "No data for FED " << *i;
109  report->addError(*i);
110  }
111  } else if (fed.size()<8*3) {
112  if (!silent_) edm::LogWarning("EmptyData") << "Tiny data " << fed.size() << " for FED " << *i;
113  report->addError(*i);
114  } else {
115  try {
116  unpacker_.unpack(fed,*readoutMap,colls, *report,silent_);
117  report->addUnpacked(*i);
118  } catch (cms::Exception& e) {
119  if (!silent_) edm::LogWarning("Unpacking error") << e.what();
120  report->addError(*i);
121  } catch (...) {
122  if (!silent_) edm::LogWarning("Unpacking exception");
123  report->addError(*i);
124  }
125  }
126  }
127 
128 
129  // gather statistics
130  stats_.max_hbhe=std::max(stats_.max_hbhe,(int)hbhe.size());
131  stats_.ave_hbhe=(stats_.ave_hbhe*stats_.n+hbhe.size())/(stats_.n+1);
132  stats_.max_ho=std::max(stats_.max_ho,(int)ho.size());
133  stats_.ave_ho=(stats_.ave_ho*stats_.n+ho.size())/(stats_.n+1);
134  stats_.max_hf=std::max(stats_.max_hf,(int)hf.size());
135  stats_.ave_hf=(stats_.ave_hf*stats_.n+hf.size())/(stats_.n+1);
136  stats_.max_tp=std::max(stats_.max_tp,(int)htp.size());
137  stats_.ave_tp=(stats_.ave_tp*stats_.n+htp.size())/(stats_.n+1);
138  stats_.max_tpho=std::max(stats_.max_tpho,(int)hotp.size());
139  stats_.ave_tpho=(stats_.ave_tpho*stats_.n+hotp.size())/(stats_.n+1);
140  stats_.max_calib=std::max(stats_.max_calib,(int)hc.size());
141  stats_.ave_calib=(stats_.ave_calib*stats_.n+hc.size())/(stats_.n+1);
142 
143 
144  stats_.n++;
145 
146  // Step B: encapsulate vectors in actual collections
147  std::auto_ptr<HBHEDigiCollection> hbhe_prod(new HBHEDigiCollection());
148  std::auto_ptr<HFDigiCollection> hf_prod(new HFDigiCollection());
149  std::auto_ptr<HODigiCollection> ho_prod(new HODigiCollection());
150  std::auto_ptr<HcalTrigPrimDigiCollection> htp_prod(new HcalTrigPrimDigiCollection());
151  std::auto_ptr<HOTrigPrimDigiCollection> hotp_prod(new HOTrigPrimDigiCollection());
152 
153  hbhe_prod->swap_contents(hbhe);
154  hf_prod->swap_contents(hf);
155  ho_prod->swap_contents(ho);
156  htp_prod->swap_contents(htp);
157  hotp_prod->swap_contents(hotp);
158 
159  // Step C2: filter FEDs, if required
160  if (filter_.active()) {
161  HBHEDigiCollection filtered_hbhe=filter_.filter(*hbhe_prod,*report);
162  HODigiCollection filtered_ho=filter_.filter(*ho_prod,*report);
163  HFDigiCollection filtered_hf=filter_.filter(*hf_prod,*report);
164 
165  hbhe_prod->swap(filtered_hbhe);
166  ho_prod->swap(filtered_ho);
167  hf_prod->swap(filtered_hf);
168  }
169 
170 
171  // Step D: Put outputs into event
172  // just until the sorting is proven
173  hbhe_prod->sort();
174  ho_prod->sort();
175  hf_prod->sort();
176  htp_prod->sort();
177  hotp_prod->sort();
178 
179  e.put(hbhe_prod);
180  e.put(ho_prod);
181  e.put(hf_prod);
182  e.put(htp_prod);
183  e.put(hotp_prod);
184 
186  if (unpackCalib_) {
187  std::auto_ptr<HcalCalibDigiCollection> hc_prod(new HcalCalibDigiCollection());
188  hc_prod->swap_contents(hc);
189 
190  if (filter_.active()) {
191  HcalCalibDigiCollection filtered_calib=filter_.filter(*hc_prod,*report);
192  hc_prod->swap(filtered_calib);
193  }
194 
195  hc_prod->sort();
196  e.put(hc_prod);
197  }
198 
200  if (unpackZDC_) {
201  std::auto_ptr<ZDCDigiCollection> prod(new ZDCDigiCollection());
202  prod->swap_contents(zdc);
203 
204  if (filter_.active()) {
205  ZDCDigiCollection filtered_zdc=filter_.filter(*prod,*report);
206  prod->swap(filtered_zdc);
207  }
208 
209  prod->sort();
210  e.put(prod);
211  }
212 
213  if (unpackTTP_) {
214  std::auto_ptr<HcalTTPDigiCollection> prod(new HcalTTPDigiCollection());
215  prod->swap_contents(ttp);
216 
217  prod->sort();
218  e.put(prod);
219  }
220  e.put(report);
221 
222 
223 }
224 
225 
virtual char const * what() const
Definition: Exception.cc:141
int i
Definition: DBlmapReader.cc:9
HcalDataFrameFilter filter_
Definition: HcalRawToDigi.h:36
std::vector< HcalTTPDigi > * ttp
Definition: HcalUnpacker.h:31
edm::SortedCollection< HcalTriggerPrimitiveDigi > HcalTrigPrimDigiCollection
std::vector< HFDataFrame > * hfCont
Definition: HcalUnpacker.h:26
void swap(SortedCollection &other)
virtual void produce(edm::Event &e, const edm::EventSetup &c)
std::vector< HBHEDataFrame > * hbheCont
Definition: HcalUnpacker.h:24
std::vector< HOTriggerPrimitiveDigi > * tphoCont
Definition: HcalUnpacker.h:30
edm::SortedCollection< HOTriggerPrimitiveDigi > HOTrigPrimDigiCollection
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:49
edm::SortedCollection< ZDCDataFrame > ZDCDigiCollection
tuple report
Definition: zeeHLT_cff.py:9
edm::SortedCollection< HODataFrame > HODigiCollection
bool complainEmptyData_
Definition: HcalRawToDigi.h:40
const T & max(const T &a, const T &b)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
virtual ~HcalRawToDigi()
std::vector< HcalTriggerPrimitiveDigi > * tpCont
Definition: HcalUnpacker.h:29
struct HcalRawToDigi::Statistics stats_
int expectedOrbitMessageTime_
Definition: HcalRawToDigi.h:41
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
HcalRawToDigi(const edm::ParameterSet &ps)
tuple conf
Definition: dbtoconf.py:185
void setExpectedOrbitMessageTime(int time)
Definition: HcalUnpacker.h:38
void unpack(const FEDRawData &raw, const HcalElectronicsMap &emap, std::vector< HcalHistogramDigi > &histoDigis)
edm::SortedCollection< HcalCalibDataFrame > HcalCalibDigiCollection
std::vector< HcalCalibDataFrame > * calibCont
Definition: HcalUnpacker.h:27
std::vector< int > fedUnpackList_
Definition: HcalRawToDigi.h:37
const T & get() const
Definition: EventSetup.h:55
bool active() const
whether any filters are on
std::vector< HODataFrame > * hoCont
Definition: HcalUnpacker.h:25
HBHEDigiCollection filter(const HBHEDigiCollection &incol, HcalUnpackerReport &r)
filter HB/HE data frames
HcalUnpacker unpacker_
Definition: HcalRawToDigi.h:35
edm::SortedCollection< HcalTTPDigi > HcalTTPDigiCollection
edm::SortedCollection< HFDataFrame > HFDigiCollection
std::vector< ZDCDataFrame > * zdcCont
Definition: HcalUnpacker.h:28
edm::SortedCollection< HBHEDataFrame > HBHEDigiCollection
edm::InputTag dataTag_
Definition: HcalRawToDigi.h:34