CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelRawToDigi.cc
Go to the documentation of this file.
1 // Include parameter driven interface to SiPixelQuality for study purposes
2 // exclude ROC(raw) based on bad ROC list in SiPixelQuality
3 // enabled by: process.siPixelDigis.UseQualityInfo = True
4 // 20-10-2010 Andrew York (Tennessee)
5 
6 #include "SiPixelRawToDigi.h"
7 
11 
13 
15 
17 
21 
27 
29 
33 
34 #include "TH1D.h"
35 #include "TFile.h"
36 
37 using namespace std;
38 
39 // -----------------------------------------------------------------------------
41  : config_(conf),
42  cabling_(0),
43  badPixelInfo_(0),
44  regions_(0),
45  hCPU(0), hDigi(0), theTimer(0)
46 {
47 
48  includeErrors = config_.getParameter<bool>("IncludeErrors");
49  useQuality = config_.getParameter<bool>("UseQualityInfo");
50  useCablingTree_ = config_.getUntrackedParameter<bool>("UseCablingTree",true);
51  if (config_.exists("ErrorList")) {
52  tkerrorlist = config_.getParameter<std::vector<int> > ("ErrorList");
53  }
54  if (config_.exists("UserErrorList")) {
55  usererrorlist = config_.getParameter<std::vector<int> > ("UserErrorList");
56  }
57  tFEDRawDataCollection = consumes <FEDRawDataCollection> (config_.getParameter<edm::InputTag>("InputLabel"));
58 
59  //start counters
60  ndigis = 0;
61  nwords = 0;
62 
63  // Products
64  produces< edm::DetSetVector<PixelDigi> >();
65  if(includeErrors){
66  produces< edm::DetSetVector<SiPixelRawDataError> >();
67  produces<DetIdCollection>();
68  produces<DetIdCollection>("UserErrorModules");
69  }
70 
71  // regions
72  if (config_.exists("Regions")) {
73  if(config_.getParameter<edm::ParameterSet>("Regions").getParameterNames().size() > 0)
74  {
76  }
77  }
78 
79  // Timing
80  bool timing = config_.getUntrackedParameter<bool>("Timing",false);
81  if (timing) {
82  theTimer = new R2DTimerObserver("**** MY TIMING REPORT ***");
83  hCPU = new TH1D ("hCPU","hCPU",100,0.,0.050);
84  hDigi = new TH1D("hDigi","hDigi",50,0.,15000.);
85  }
86 }
87 
88 
89 // -----------------------------------------------------------------------------
91  edm::LogInfo("SiPixelRawToDigi") << " HERE ** SiPixelRawToDigi destructor!";
92 
93  if (useCablingTree_) delete cabling_;
94  if (regions_) delete regions_;
95 
96  if (theTimer) {
97  TFile rootFile("analysis.root", "RECREATE", "my histograms");
98  hCPU->Write();
99  hDigi->Write();
100  delete theTimer;
101  }
102 
103 }
104 
105 
106 // -----------------------------------------------------------------------------
107 
108 
109 // -----------------------------------------------------------------------------
111  const edm::EventSetup& es)
112 {
113  const uint32_t dummydetid = 0xffffffff;
115 
116 // initialize cabling map or update if necessary
117  if (recordWatcher.check( es )) {
118  // cabling map, which maps online address (fed->link->ROC->local pixel) to offline (DetId->global pixel)
119  if (useCablingTree_) {
120  delete cabling_;
121  // we are going to make our own copy so safe to let the map be deleted early
123  es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
124  fedIds = cablingMap->fedIds();
125  cabling_ = cablingMap->cablingTree();
126  } else {
127  // we are going to hold the pointer so we need the map to stick around
129  es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
130  fedIds = cablingMap->fedIds();
131  cabling_ = cablingMap.product();
132  }
133  LogDebug("map version:")<< cabling_->version();
134  }
135 // initialize quality record or update if necessary
136  if (qualityWatcher.check( es )&&useQuality) {
137  // quality info for dead pixel modules or ROCs
138  edm::ESHandle<SiPixelQuality> qualityInfo;
139  es.get<SiPixelQualityRcd>().get( qualityInfo );
140  badPixelInfo_ = qualityInfo.product();
141  if (!badPixelInfo_) {
142  edm::LogError("SiPixelQualityNotPresent")<<" Configured to use SiPixelQuality, but SiPixelQuality not present"<<endl;
143  }
144  }
145 
147  label = config_.getParameter<edm::InputTag>("InputLabel");
148  ev.getByToken(tFEDRawDataCollection, buffers);
149 
150 // create product (digis & errors)
151  std::auto_ptr< edm::DetSetVector<PixelDigi> > collection( new edm::DetSetVector<PixelDigi> );
152  std::auto_ptr< edm::DetSetVector<SiPixelRawDataError> > errorcollection( new edm::DetSetVector<SiPixelRawDataError> );
153  std::auto_ptr< DetIdCollection > tkerror_detidcollection(new DetIdCollection());
154  std::auto_ptr< DetIdCollection > usererror_detidcollection(new DetIdCollection());
155 
156  PixelDataFormatter formatter(cabling_);
157  formatter.setErrorStatus(includeErrors);
158 
160 
161  if (theTimer) theTimer->start();
162  bool errorsInEvent = false;
163  PixelDataFormatter::DetErrors nodeterrors;
164 
165  if (regions_) {
166  regions_->run(ev, es);
168  LogDebug("SiPixelRawToDigi") << "region2unpack #feds (BPIX,EPIX,total): "<<regions_->nBarrelFEDs()<<" "<<regions_->nForwardFEDs()<<" "<<regions_->nFEDs();
169  LogDebug("SiPixelRawToDigi") << "region2unpack #modules (BPIX,EPIX,total): "<<regions_->nBarrelModules()<<" "<<regions_->nForwardModules()<<" "<<regions_->nModules();
170  }
171 
172  typedef std::vector<unsigned int>::const_iterator IF;
173  for (IF aFed = fedIds.begin(); aFed != fedIds.end(); ++aFed) {
174  int fedId = *aFed;
175 
176  if (regions_ && !regions_->mayUnpackFED(fedId)) continue;
177 
178  if(debug) LogDebug("SiPixelRawToDigi")<< " PRODUCE DIGI FOR FED: " << fedId << endl;
181 
182  //get event data for this fed
183  const FEDRawData& fedRawData = buffers->FEDData( fedId );
184 
185  //convert data to digi and strip off errors
186  formatter.interpretRawData( errorsInEvent, fedId, fedRawData, digis, errors);
187 
188  //pack digi into collection
189  typedef PixelDataFormatter::Digis::iterator ID;
190  for (ID it = digis.begin(); it != digis.end(); it++) {
191  uint32_t detid = it->first;
192  edm::DetSet<PixelDigi>& detSet = collection->find_or_insert(detid);
193  detSet.data.insert(detSet.data.end(), it->second.begin(), it->second.end());
194  }
195 
196  //pack errors into collection
197  if(includeErrors) {
198  typedef PixelDataFormatter::Errors::iterator IE;
199  for (IE is = errors.begin(); is != errors.end(); is++) {
200  uint32_t errordetid = is->first;
201  if (errordetid==dummydetid) { // errors given dummy detId must be sorted by Fed
202  nodeterrors.insert( nodeterrors.end(), errors[errordetid].begin(), errors[errordetid].end() );
203  } else {
204  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(errordetid);
205  errorDetSet.data.insert(errorDetSet.data.end(), is->second.begin(), is->second.end());
206  // Fill detid of the detectors where there is error AND the error number is listed
207  // in the configurable error list in the job option cfi.
208  // Code needs to be here, because there can be a set of errors for each
209  // entry in the for loop over PixelDataFormatter::Errors
210  if(!tkerrorlist.empty() || !usererrorlist.empty()){
211  DetId errorDetId(errordetid);
212  edm::DetSet<SiPixelRawDataError>::const_iterator itPixelError=errorDetSet.begin();
213  for(; itPixelError!=errorDetSet.end(); ++itPixelError){
214  // fill list of detIds to be turned off by tracking
215  if(!tkerrorlist.empty()) {
216  std::vector<int>::iterator it_find = find(tkerrorlist.begin(), tkerrorlist.end(), itPixelError->getType());
217  if(it_find != tkerrorlist.end()){
218  tkerror_detidcollection->push_back(errordetid);
219  }
220  }
221  // fill list of detIds with errors to be studied
222  if(!usererrorlist.empty()) {
223  std::vector<int>::iterator it_find = find(usererrorlist.begin(), usererrorlist.end(), itPixelError->getType());
224  if(it_find != usererrorlist.end()){
225  usererror_detidcollection->push_back(errordetid);
226  }
227  }
228  }
229  }
230  }
231  }
232  }
233  }
234 
235  if(includeErrors) {
236  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(dummydetid);
237  errorDetSet.data = nodeterrors;
238  }
239  if (errorsInEvent) LogDebug("SiPixelRawToDigi") << "Error words were stored in this event";
240 
241  if (theTimer) {
242  theTimer->stop();
243  LogDebug("SiPixelRawToDigi") << "TIMING IS: (real)" << theTimer->lastMeasurement().real() ;
244  ndigis += formatter.nDigis();
245  nwords += formatter.nWords();
246  LogDebug("SiPixelRawToDigi") << " (Words/Digis) this ev: "
247  <<formatter.nWords()<<"/"<<formatter.nDigis() << "--- all :"<<nwords<<"/"<<ndigis;
248  hCPU->Fill( theTimer->lastMeasurement().real() );
249  hDigi->Fill(formatter.nDigis());
250  }
251 
252  //send digis and errors back to framework
253  ev.put( collection );
254  if(includeErrors){
255  ev.put( errorcollection );
256  ev.put( tkerror_detidcollection );
257  ev.put( usererror_detidcollection, "UserErrorModules" );
258  }
259 }
#define LogDebug(id)
iterator end()
Definition: DetSet.h:60
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::map< cms_uint32_t, DetDigis > Digis
unsigned int nForwardFEDs() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
uint32_t ID
Definition: Definitions.h:26
static MessageDrop * instance()
Definition: MessageDrop.cc:60
bool exists(std::string const &parameterName) const
checks if a parameter exists
void interpretRawData(bool &errorsInEvent, int fedId, const FEDRawData &data, Digis &digis, Errors &errors)
std::vector< int > usererrorlist
const std::set< unsigned int > * modulesToUnpack() const
full set of module ids to unpack
void setErrorStatus(bool ErrorStatus)
virtual std::string version() const =0
edm::InputTag label
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
PixelUnpackingRegions * regions_
std::vector< SiPixelRawDataError > DetErrors
std::vector< unsigned int > fedIds
SiPixelRawToDigi(const edm::ParameterSet &)
ctor
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
unsigned int nBarrelModules() const
edm::EDGetTokenT< FEDRawDataCollection > tFEDRawDataCollection
void setQualityStatus(bool QualityStatus, const SiPixelQuality *QualityInfo)
std::vector< std::string > getParameterNames() const
tuple conf
Definition: dbtoconf.py:185
std::map< cms_uint32_t, DetErrors > Errors
virtual ~SiPixelRawToDigi()
dtor
iterator begin()
Definition: DetSet.h:59
edm::ESWatcher< SiPixelFedCablingMapRcd > recordWatcher
void run(const edm::Event &e, const edm::EventSetup &es)
has to be run during each event
Definition: DetId.h:18
unsigned int nModules() const
void setModulesToUnpack(const std::set< unsigned int > *moduleIds)
const LastMeasurement & lastMeasurement()
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::vector< int > tkerrorlist
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:58
edm::ParameterSet config_
unsigned int nForwardModules() const
collection_type data
Definition: DetSet.h:78
unsigned int nFEDs() const
various informational accessors:
unsigned int nBarrelFEDs() const
bool mayUnpackFED(unsigned int fed_n) const
check whether a FED has to be unpacked
edm::EDCollection< DetId > DetIdCollection
R2DTimerObserver * theTimer
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
edm::ESWatcher< SiPixelQualityRcd > qualityWatcher
const SiPixelQuality * badPixelInfo_
const SiPixelFedCabling * cabling_
virtual void produce(edm::Event &, const edm::EventSetup &)
get data, convert to digis attach againe to Event