CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/EventFilter/HcalRawToDigi/plugins/HcalRawToDigi.cc

Go to the documentation of this file.
00001 using namespace std;
00002 #include "EventFilter/HcalRawToDigi/plugins/HcalRawToDigi.h"
00003 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00004 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00005 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00008 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include <iostream>
00011 
00012 HcalRawToDigi::HcalRawToDigi(edm::ParameterSet const& conf):
00013   dataTag_(conf.getParameter<edm::InputTag>("InputLabel")),
00014   unpacker_(conf.getUntrackedParameter<int>("HcalFirstFED",int(FEDNumbering::MINHCALFEDID)),conf.getParameter<int>("firstSample"),conf.getParameter<int>("lastSample")),
00015   filter_(conf.getParameter<bool>("FilterDataQuality"),conf.getParameter<bool>("FilterDataQuality"),
00016           false,
00017           0, 0, 
00018           -1),
00019   fedUnpackList_(conf.getUntrackedParameter<std::vector<int> >("FEDs", std::vector<int>())),
00020   firstFED_(conf.getUntrackedParameter<int>("HcalFirstFED",FEDNumbering::MINHCALFEDID)),
00021   unpackCalib_(conf.getUntrackedParameter<bool>("UnpackCalib",false)),
00022   unpackZDC_(conf.getUntrackedParameter<bool>("UnpackZDC",false)),
00023   unpackTTP_(conf.getUntrackedParameter<bool>("UnpackTTP",false)),
00024   silent_(conf.getUntrackedParameter<bool>("silent",true)),
00025   complainEmptyData_(conf.getUntrackedParameter<bool>("ComplainEmptyData",false)),
00026   expectedOrbitMessageTime_(conf.getUntrackedParameter<int>("ExpectedOrbitMessageTime",-1))
00027 {
00028   if (fedUnpackList_.empty()) {
00029     for (int i=FEDNumbering::MINHCALFEDID; i<=FEDNumbering::MAXHCALFEDID; i++)
00030       fedUnpackList_.push_back(i);
00031   } 
00032   
00033   unpacker_.setExpectedOrbitMessageTime(expectedOrbitMessageTime_);
00034   std::ostringstream ss;
00035   for (unsigned int i=0; i<fedUnpackList_.size(); i++) 
00036     ss << fedUnpackList_[i] << " ";
00037   edm::LogInfo("HCAL") << "HcalRawToDigi will unpack FEDs ( " << ss.str() << ")";
00038     
00039   // products produced...
00040   produces<HBHEDigiCollection>();
00041   produces<HFDigiCollection>();
00042   produces<HODigiCollection>();
00043   produces<HcalTrigPrimDigiCollection>();
00044   produces<HOTrigPrimDigiCollection>();
00045   produces<HcalUnpackerReport>();
00046   if (unpackCalib_)
00047     produces<HcalCalibDigiCollection>();
00048   if (unpackZDC_)
00049     produces<ZDCDigiCollection>();
00050   if (unpackTTP_)
00051     produces<HcalTTPDigiCollection>();
00052 
00053   memset(&stats_,0,sizeof(stats_));
00054 
00055 }
00056 
00057 // Virtual destructor needed.
00058 HcalRawToDigi::~HcalRawToDigi() { }  
00059 
00060 // Functions that gets called by framework every event
00061 void HcalRawToDigi::produce(edm::Event& e, const edm::EventSetup& es)
00062 {
00063   // Step A: Get Inputs 
00064   edm::Handle<FEDRawDataCollection> rawraw;  
00065   e.getByLabel(dataTag_,rawraw);
00066   // get the mapping
00067   edm::ESHandle<HcalDbService> pSetup;
00068   es.get<HcalDbRecord>().get( pSetup );
00069   const HcalElectronicsMap* readoutMap=pSetup->getHcalMapping();
00070   
00071   // Step B: Create empty output  : three vectors for three classes...
00072   std::vector<HBHEDataFrame> hbhe;
00073   std::vector<HODataFrame> ho;
00074   std::vector<HFDataFrame> hf;
00075   std::vector<HcalTriggerPrimitiveDigi> htp;
00076   std::vector<HcalCalibDataFrame> hc;
00077   std::vector<ZDCDataFrame> zdc;
00078   std::vector<HcalTTPDigi> ttp;
00079   std::vector<HOTriggerPrimitiveDigi> hotp;
00080   std::auto_ptr<HcalUnpackerReport> report(new HcalUnpackerReport);
00081 
00082   // Heuristics: use ave+(max-ave)/8
00083   if (stats_.max_hbhe>0) hbhe.reserve(stats_.ave_hbhe+(stats_.max_hbhe-stats_.ave_hbhe)/8);
00084   if (stats_.max_ho>0) ho.reserve(stats_.ave_ho+(stats_.max_ho-stats_.ave_ho)/8);
00085   if (stats_.max_hf>0) hf.reserve(stats_.ave_hf+(stats_.max_hf-stats_.ave_hf)/8);
00086   if (stats_.max_calib>0) hc.reserve(stats_.ave_calib+(stats_.max_calib-stats_.ave_calib)/8);
00087   if (stats_.max_tp>0) htp.reserve(stats_.ave_tp+(stats_.max_tp-stats_.ave_tp)/8);
00088   if (stats_.max_tpho>0) hotp.reserve(stats_.ave_tpho+(stats_.max_tpho-stats_.ave_tpho)/8);
00089 
00090   if (unpackZDC_) zdc.reserve(24);
00091 
00092 
00093   HcalUnpacker::Collections colls;
00094   colls.hbheCont=&hbhe;
00095   colls.hoCont=&ho;
00096   colls.hfCont=&hf;
00097   colls.tpCont=&htp;
00098   colls.tphoCont=&hotp;
00099   colls.calibCont=&hc;
00100   colls.zdcCont=&zdc;
00101   if (unpackTTP_) colls.ttp=&ttp;
00102  
00103   // Step C: unpack all requested FEDs
00104   for (std::vector<int>::const_iterator i=fedUnpackList_.begin(); i!=fedUnpackList_.end(); i++) {
00105     const FEDRawData& fed = rawraw->FEDData(*i);
00106     if (fed.size()==0) {
00107       if (complainEmptyData_) {
00108         if (!silent_) edm::LogWarning("EmptyData") << "No data for FED " << *i;
00109         report->addError(*i);
00110       }
00111     } else if (fed.size()<8*3) {
00112       if (!silent_) edm::LogWarning("EmptyData") << "Tiny data " << fed.size() << " for FED " << *i;
00113       report->addError(*i);
00114     } else {
00115       try {
00116         unpacker_.unpack(fed,*readoutMap,colls, *report,silent_);
00117         report->addUnpacked(*i);
00118       } catch (cms::Exception& e) {
00119         if (!silent_) edm::LogWarning("Unpacking error") << e.what();
00120         report->addError(*i);
00121       } catch (...) {
00122         if (!silent_) edm::LogWarning("Unpacking exception");
00123         report->addError(*i);
00124       }
00125     }
00126   }
00127 
00128 
00129   // gather statistics
00130   stats_.max_hbhe=std::max(stats_.max_hbhe,(int)hbhe.size());
00131   stats_.ave_hbhe=(stats_.ave_hbhe*stats_.n+hbhe.size())/(stats_.n+1);
00132   stats_.max_ho=std::max(stats_.max_ho,(int)ho.size());
00133   stats_.ave_ho=(stats_.ave_ho*stats_.n+ho.size())/(stats_.n+1);
00134   stats_.max_hf=std::max(stats_.max_hf,(int)hf.size());
00135   stats_.ave_hf=(stats_.ave_hf*stats_.n+hf.size())/(stats_.n+1);
00136   stats_.max_tp=std::max(stats_.max_tp,(int)htp.size());
00137   stats_.ave_tp=(stats_.ave_tp*stats_.n+htp.size())/(stats_.n+1);
00138   stats_.max_tpho=std::max(stats_.max_tpho,(int)hotp.size());
00139   stats_.ave_tpho=(stats_.ave_tpho*stats_.n+hotp.size())/(stats_.n+1);
00140   stats_.max_calib=std::max(stats_.max_calib,(int)hc.size());
00141   stats_.ave_calib=(stats_.ave_calib*stats_.n+hc.size())/(stats_.n+1);
00142 
00143 
00144   stats_.n++;
00145 
00146   // Step B: encapsulate vectors in actual collections
00147   std::auto_ptr<HBHEDigiCollection> hbhe_prod(new HBHEDigiCollection()); 
00148   std::auto_ptr<HFDigiCollection> hf_prod(new HFDigiCollection());
00149   std::auto_ptr<HODigiCollection> ho_prod(new HODigiCollection());
00150   std::auto_ptr<HcalTrigPrimDigiCollection> htp_prod(new HcalTrigPrimDigiCollection());  
00151   std::auto_ptr<HOTrigPrimDigiCollection> hotp_prod(new HOTrigPrimDigiCollection());  
00152 
00153   hbhe_prod->swap_contents(hbhe);
00154   hf_prod->swap_contents(hf);
00155   ho_prod->swap_contents(ho);
00156   htp_prod->swap_contents(htp);
00157   hotp_prod->swap_contents(hotp);
00158 
00159   // Step C2: filter FEDs, if required
00160   if (filter_.active()) {
00161     HBHEDigiCollection filtered_hbhe=filter_.filter(*hbhe_prod,*report);
00162     HODigiCollection filtered_ho=filter_.filter(*ho_prod,*report);
00163     HFDigiCollection filtered_hf=filter_.filter(*hf_prod,*report);
00164     
00165     hbhe_prod->swap(filtered_hbhe);
00166     ho_prod->swap(filtered_ho);
00167     hf_prod->swap(filtered_hf);    
00168   }
00169 
00170 
00171   // Step D: Put outputs into event
00172   // just until the sorting is proven
00173   hbhe_prod->sort();
00174   ho_prod->sort();
00175   hf_prod->sort();
00176   htp_prod->sort();
00177   hotp_prod->sort();
00178 
00179   e.put(hbhe_prod);
00180   e.put(ho_prod);
00181   e.put(hf_prod);
00182   e.put(htp_prod);
00183   e.put(hotp_prod);
00184 
00186   if (unpackCalib_) {
00187     std::auto_ptr<HcalCalibDigiCollection> hc_prod(new HcalCalibDigiCollection());
00188     hc_prod->swap_contents(hc);
00189 
00190     if (filter_.active()) {
00191       HcalCalibDigiCollection filtered_calib=filter_.filter(*hc_prod,*report);
00192       hc_prod->swap(filtered_calib);
00193     }
00194 
00195     hc_prod->sort();
00196     e.put(hc_prod);
00197   }
00198 
00200   if (unpackZDC_) {
00201     std::auto_ptr<ZDCDigiCollection> prod(new ZDCDigiCollection());
00202     prod->swap_contents(zdc);
00203     
00204     if (filter_.active()) {
00205       ZDCDigiCollection filtered_zdc=filter_.filter(*prod,*report);
00206       prod->swap(filtered_zdc);
00207     }
00208 
00209     prod->sort();
00210     e.put(prod);
00211   }
00212 
00213   if (unpackTTP_) {
00214     std::auto_ptr<HcalTTPDigiCollection> prod(new HcalTTPDigiCollection());
00215     prod->swap_contents(ttp);
00216     
00217     prod->sort();
00218     e.put(prod);
00219   }
00220   e.put(report);
00221 
00222 
00223 }
00224 
00225