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