Go to the documentation of this file.00001
00002
00003
00004
00005
00006 #include "SiPixelRawToDigi.h"
00007
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "FWCore/Framework/interface/ESTransientHandle.h"
00011
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013
00014 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
00015
00016 #include "DataFormats/SiPixelRawData/interface/SiPixelRawDataError.h"
00017
00018 #include "DataFormats/Common/interface/DetSetVector.h"
00019 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00020 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
00021
00022 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00023
00024 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
00025 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingTree.h"
00026 #include "EventFilter/SiPixelRawToDigi/interface/PixelDataFormatter.h"
00027
00028 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
00029
00030 #include "EventFilter/SiPixelRawToDigi/interface/R2DTimerObserver.h"
00031
00032 #include "TH1D.h"
00033 #include "TFile.h"
00034
00035 using namespace std;
00036
00037
00038 SiPixelRawToDigi::SiPixelRawToDigi( const edm::ParameterSet& conf )
00039 : config_(conf),
00040 cabling_(0),
00041 badPixelInfo_(0),
00042 hCPU(0), hDigi(0), theTimer(0)
00043 {
00044
00045 includeErrors = config_.getParameter<bool>("IncludeErrors");
00046 useQuality = config_.getParameter<bool>("UseQualityInfo");
00047 useCablingTree_ = config_.getUntrackedParameter<bool>("UseCablingTree",true);
00048
00049
00050 ndigis = 0;
00051 nwords = 0;
00052
00053
00054 produces< edm::DetSetVector<PixelDigi> >();
00055 if(includeErrors) produces< edm::DetSetVector<SiPixelRawDataError> >();
00056
00057
00058 bool timing = config_.getUntrackedParameter<bool>("Timing",false);
00059 if (timing) {
00060 theTimer = new R2DTimerObserver("**** MY TIMING REPORT ***");
00061 hCPU = new TH1D ("hCPU","hCPU",100,0.,0.050);
00062 hDigi = new TH1D("hDigi","hDigi",50,0.,15000.);
00063 }
00064 }
00065
00066
00067
00068 SiPixelRawToDigi::~SiPixelRawToDigi() {
00069 edm::LogInfo("SiPixelRawToDigi") << " HERE ** SiPixelRawToDigi destructor!";
00070
00071 if(useCablingTree_) delete cabling_;
00072
00073 if (theTimer) {
00074 TFile rootFile("analysis.root", "RECREATE", "my histograms");
00075 hCPU->Write();
00076 hDigi->Write();
00077 delete theTimer;
00078 }
00079
00080 }
00081
00082
00083
00084
00085
00086
00087 void SiPixelRawToDigi::produce( edm::Event& ev,
00088 const edm::EventSetup& es)
00089 {
00090 const uint32_t dummydetid = 0xffffffff;
00091 debug = edm::MessageDrop::instance()->debugEnabled;
00092
00093
00094 if (recordWatcher.check( es )) {
00095
00096 edm::ESTransientHandle<SiPixelFedCablingMap> cablingMap;
00097 es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
00098 fedIds = cablingMap->fedIds();
00099 if (useCablingTree_ && cabling_) delete cabling_;
00100 if (useCablingTree_) cabling_ = cablingMap->cablingTree();
00101 else cabling_ = cablingMap.product();
00102 LogDebug("map version:")<< cabling_->version();
00103 }
00104
00105 if (qualityWatcher.check( es )&&useQuality) {
00106
00107 edm::ESHandle<SiPixelQuality> qualityInfo;
00108 es.get<SiPixelQualityRcd>().get( qualityInfo );
00109 badPixelInfo_ = qualityInfo.product();
00110 if (!badPixelInfo_) {
00111 edm::LogError("**SiPixelRawToDigi**")<<" Configured to use SiPixelQuality, but SiPixelQuality not present"<<endl;
00112 }
00113 }
00114
00115 edm::Handle<FEDRawDataCollection> buffers;
00116 label = config_.getParameter<edm::InputTag>("InputLabel");
00117 ev.getByLabel( label, buffers);
00118
00119
00120 std::auto_ptr< edm::DetSetVector<PixelDigi> > collection( new edm::DetSetVector<PixelDigi> );
00121 std::auto_ptr< edm::DetSetVector<SiPixelRawDataError> > errorcollection( new edm::DetSetVector<SiPixelRawDataError> );
00122
00123 PixelDataFormatter formatter(cabling_);
00124 formatter.setErrorStatus(includeErrors);
00125 if (useQuality) formatter.setQualityStatus(useQuality, badPixelInfo_);
00126
00127 if (theTimer) theTimer->start();
00128 bool errorsInEvent = false;
00129 PixelDataFormatter::DetErrors nodeterrors;
00130
00131 typedef std::vector<unsigned int>::const_iterator IF;
00132 for (IF aFed = fedIds.begin(); aFed != fedIds.end(); ++aFed) {
00133 int fedId = *aFed;
00134 if(debug) LogDebug("SiPixelRawToDigi")<< " PRODUCE DIGI FOR FED: " << fedId << endl;
00135 PixelDataFormatter::Digis digis;
00136 PixelDataFormatter::Errors errors;
00137
00138
00139 const FEDRawData& fedRawData = buffers->FEDData( fedId );
00140
00141
00142 formatter.interpretRawData( errorsInEvent, fedId, fedRawData, digis, errors);
00143
00144
00145 typedef PixelDataFormatter::Digis::iterator ID;
00146 for (ID it = digis.begin(); it != digis.end(); it++) {
00147 uint32_t detid = it->first;
00148 edm::DetSet<PixelDigi>& detSet = collection->find_or_insert(detid);
00149 detSet.data = it->second;
00150 }
00151
00152
00153 if(includeErrors) {
00154 typedef PixelDataFormatter::Errors::iterator IE;
00155 for (IE is = errors.begin(); is != errors.end(); is++) {
00156 uint32_t errordetid = is->first;
00157 if (errordetid==dummydetid) {
00158 nodeterrors.insert( nodeterrors.end(), errors[errordetid].begin(), errors[errordetid].end() );
00159 } else {
00160 edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(errordetid);
00161 errorDetSet.data = is->second;
00162 }
00163 }
00164 }
00165 }
00166
00167 if(includeErrors) {
00168 edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(dummydetid);
00169 errorDetSet.data = nodeterrors;
00170 }
00171 if (errorsInEvent) LogDebug("SiPixelRawToDigi") << "Error words were stored in this event";
00172
00173 if (theTimer) {
00174 theTimer->stop();
00175 LogDebug("SiPixelRawToDigi") << "TIMING IS: (real)" << theTimer->lastMeasurement().real() ;
00176 ndigis += formatter.nDigis();
00177 nwords += formatter.nWords();
00178 LogDebug("SiPixelRawToDigi") << " (Words/Digis) this ev: "
00179 <<formatter.nWords()<<"/"<<formatter.nDigis() << "--- all :"<<nwords<<"/"<<ndigis;
00180 hCPU->Fill( theTimer->lastMeasurement().real() );
00181 hDigi->Fill(formatter.nDigis());
00182 }
00183
00184
00185 ev.put( collection );
00186 if(includeErrors) ev.put( errorcollection );
00187 }