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 // Skip FED40 pilot-blade
2 // Include parameter driven interface to SiPixelQuality for study purposes
3 // exclude ROC(raw) based on bad ROC list in SiPixelQuality
4 // enabled by: process.siPixelDigis.UseQualityInfo = True (BY DEFAULT NOT USED)
5 // 20-10-2010 Andrew York (Tennessee)
6 
7 #include "SiPixelRawToDigi.h"
8 
14 
16 
18 
20 
24 
30 
32 
35 
36 #include "TH1D.h"
37 #include "TFile.h"
38 
39 using namespace std;
40 
41 // -----------------------------------------------------------------------------
43  : config_(conf),
44  badPixelInfo_(0),
45  regions_(0),
46  hCPU(0), hDigi(0)
47 {
48 
49  includeErrors = config_.getParameter<bool>("IncludeErrors");
50  useQuality = config_.getParameter<bool>("UseQualityInfo");
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.reset( new edm::CPUTimer );
83  hCPU = new TH1D ("hCPU","hCPU",100,0.,0.050);
84  hDigi = new TH1D("hDigi","hDigi",50,0.,15000.);
85  }
86 
87  // Control the usage of pilot-blade data, FED=40
88  usePilotBlade = false;
89  if (config_.exists("UsePilotBlade")) {
90  usePilotBlade = config_.getParameter<bool> ("UsePilotBlade");
91  if(usePilotBlade) edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
92  }
93 
94  // Control the usage of phase1
95  usePhase1 = false;
96  if (config_.exists("UsePhase1")) {
97  usePhase1 = config_.getParameter<bool> ("UsePhase1");
98  if(usePhase1) edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
99  }
100 
101 }
102 
103 
104 // -----------------------------------------------------------------------------
106  edm::LogInfo("SiPixelRawToDigi") << " HERE ** SiPixelRawToDigi destructor!";
107 
108  if (regions_) delete regions_;
109 
110  if (theTimer) {
111  TFile rootFile("analysis.root", "RECREATE", "my histograms");
112  hCPU->Write();
113  hDigi->Write();
114  }
115 
116 }
117 
118 void
121  desc.add<bool>("IncludeErrors",true);
122  desc.add<bool>("UseQualityInfo",false);
123  {
124  std::vector<int> temp1;
125  temp1.reserve(1);
126  temp1.push_back(29);
127  desc.add<std::vector<int> >("ErrorList",temp1)->setComment("## ErrorList: list of error codes used by tracking to invalidate modules");
128  }
129  {
130  std::vector<int> temp1;
131  temp1.reserve(1);
132  temp1.push_back(40);
133  desc.add<std::vector<int> >("UserErrorList",temp1)->setComment("## UserErrorList: list of error codes used by Pixel experts for investigation");
134  }
135  desc.add<edm::InputTag>("InputLabel",edm::InputTag("siPixelRawData"));
136  {
138  psd0.addOptional<std::vector<edm::InputTag>>("inputs");
139  psd0.addOptional<std::vector<double>>("deltaPhi");
140  psd0.addOptional<std::vector<double>>("maxZ");
141  psd0.addOptional<edm::InputTag>("beamSpot");
142  desc.add<edm::ParameterSetDescription>("Regions",psd0)->setComment("## Empty Regions PSet means complete unpacking");
143  }
144  desc.addUntracked<bool>("Timing",false);
145  desc.add<bool>("UsePilotBlade",false)->setComment("## Use pilot blades");
146  desc.add<bool>("UsePhase1",false)->setComment("## Use phase1");
147  desc.addOptional<bool>("CheckPixelOrder"); // never used, kept for back-compatibility
148  descriptions.add("siPixelRawToDigi",desc);
149 }
150 
151 // -----------------------------------------------------------------------------
152 
153 
154 // -----------------------------------------------------------------------------
156  const edm::EventSetup& es)
157 {
158  const uint32_t dummydetid = 0xffffffff;
160 
161 // initialize cabling map or update if necessary
162  if (recordWatcher.check( es )) {
163  // cabling map, which maps online address (fed->link->ROC->local pixel) to offline (DetId->global pixel)
165  es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
166  fedIds = cablingMap->fedIds();
167  cabling_ = cablingMap->cablingTree();
168  LogDebug("map version:")<< cabling_->version();
169  }
170 // initialize quality record or update if necessary
171  if (qualityWatcher.check( es )&&useQuality) {
172  // quality info for dead pixel modules or ROCs
173  edm::ESHandle<SiPixelQuality> qualityInfo;
174  es.get<SiPixelQualityRcd>().get( qualityInfo );
175  badPixelInfo_ = qualityInfo.product();
176  if (!badPixelInfo_) {
177  edm::LogError("SiPixelQualityNotPresent")<<" Configured to use SiPixelQuality, but SiPixelQuality not present"<<endl;
178  }
179  }
180 
182  ev.getByToken(tFEDRawDataCollection, buffers);
183 
184 // create product (digis & errors)
185  std::auto_ptr< edm::DetSetVector<PixelDigi> > collection( new edm::DetSetVector<PixelDigi> );
186  // collection->reserve(8*1024);
187  std::auto_ptr< edm::DetSetVector<SiPixelRawDataError> > errorcollection( new edm::DetSetVector<SiPixelRawDataError> );
188  std::auto_ptr< DetIdCollection > tkerror_detidcollection(new DetIdCollection());
189  std::auto_ptr< DetIdCollection > usererror_detidcollection(new DetIdCollection());
190 
191  //PixelDataFormatter formatter(cabling_.get()); // phase 0 only
192  PixelDataFormatter formatter(cabling_.get(), usePhase1); // for phase 1 & 0
193 
194  formatter.setErrorStatus(includeErrors);
195 
196  if (useQuality) formatter.setQualityStatus(useQuality, badPixelInfo_);
197 
198  if (theTimer) theTimer->start();
199  bool errorsInEvent = false;
200  PixelDataFormatter::DetErrors nodeterrors;
201 
202  if (regions_) {
203  regions_->run(ev, es);
204  formatter.setModulesToUnpack(regions_->modulesToUnpack());
205  LogDebug("SiPixelRawToDigi") << "region2unpack #feds (BPIX,EPIX,total): "<<regions_->nBarrelFEDs()<<" "<<regions_->nForwardFEDs()<<" "<<regions_->nFEDs();
206  LogDebug("SiPixelRawToDigi") << "region2unpack #modules (BPIX,EPIX,total): "<<regions_->nBarrelModules()<<" "<<regions_->nForwardModules()<<" "<<regions_->nModules();
207  }
208 
209  for (auto aFed = fedIds.begin(); aFed != fedIds.end(); ++aFed) {
210  int fedId = *aFed;
211 
212  if(!usePilotBlade && (fedId==40) ) continue; // skip pilot blade data
213 
214  if (regions_ && !regions_->mayUnpackFED(fedId)) continue;
215 
216  if(debug) LogDebug("SiPixelRawToDigi")<< " PRODUCE DIGI FOR FED: " << fedId << endl;
217 
219 
220  //get event data for this fed
221  const FEDRawData& fedRawData = buffers->FEDData( fedId );
222 
223  //convert data to digi and strip off errors
224  formatter.interpretRawData( errorsInEvent, fedId, fedRawData, *collection, errors);
225 
226  //pack errors into collection
227  if(includeErrors) {
228  typedef PixelDataFormatter::Errors::iterator IE;
229  for (IE is = errors.begin(); is != errors.end(); is++) {
230  uint32_t errordetid = is->first;
231  if (errordetid==dummydetid) { // errors given dummy detId must be sorted by Fed
232  nodeterrors.insert( nodeterrors.end(), errors[errordetid].begin(), errors[errordetid].end() );
233  } else {
234  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(errordetid);
235  errorDetSet.data.insert(errorDetSet.data.end(), is->second.begin(), is->second.end());
236  // Fill detid of the detectors where there is error AND the error number is listed
237  // in the configurable error list in the job option cfi.
238  // Code needs to be here, because there can be a set of errors for each
239  // entry in the for loop over PixelDataFormatter::Errors
240  if(!tkerrorlist.empty() || !usererrorlist.empty()){
241  DetId errorDetId(errordetid);
242  edm::DetSet<SiPixelRawDataError>::const_iterator itPixelError=errorDetSet.begin();
243  for(; itPixelError!=errorDetSet.end(); ++itPixelError){
244  // fill list of detIds to be turned off by tracking
245  if(!tkerrorlist.empty()) {
246  std::vector<int>::iterator it_find = find(tkerrorlist.begin(), tkerrorlist.end(), itPixelError->getType());
247  if(it_find != tkerrorlist.end()){
248  tkerror_detidcollection->push_back(errordetid);
249  }
250  }
251  // fill list of detIds with errors to be studied
252  if(!usererrorlist.empty()) {
253  std::vector<int>::iterator it_find = find(usererrorlist.begin(), usererrorlist.end(), itPixelError->getType());
254  if(it_find != usererrorlist.end()){
255  usererror_detidcollection->push_back(errordetid);
256  }
257  }
258  }
259  }
260  }
261  }
262  }
263  }
264 
265  if(includeErrors) {
266  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(dummydetid);
267  errorDetSet.data = nodeterrors;
268  }
269  if (errorsInEvent) LogDebug("SiPixelRawToDigi") << "Error words were stored in this event";
270 
271  if (theTimer) {
272  theTimer->stop();
273  LogDebug("SiPixelRawToDigi") << "TIMING IS: (real)" << theTimer->realTime() ;
274  ndigis += formatter.nDigis();
275  nwords += formatter.nWords();
276  LogDebug("SiPixelRawToDigi") << " (Words/Digis) this ev: "
277  <<formatter.nWords()<<"/"<<formatter.nDigis() << "--- all :"<<nwords<<"/"<<ndigis;
278  hCPU->Fill( theTimer->realTime() );
279  hDigi->Fill(formatter.nDigis());
280  }
281 
282  //send digis and errors back to framework
283  ev.put( collection );
284  if(includeErrors){
285  ev.put( errorcollection );
286  ev.put( tkerror_detidcollection );
287  ev.put( usererror_detidcollection, "UserErrorModules" );
288  }
289 }
#define LogDebug(id)
iterator end()
Definition: DetSet.h:60
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< edm::CPUTimer > theTimer
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
unsigned int nForwardFEDs() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
static MessageDrop * instance()
Definition: MessageDrop.cc:60
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< int > usererrorlist
const std::set< unsigned int > * modulesToUnpack() const
full set of module ids to unpack
bool ev
void setErrorStatus(bool ErrorStatus)
virtual void produce(edm::Event &, const edm::EventSetup &) override
get data, convert to digis attach againe to Event
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:120
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
std::unique_ptr< SiPixelFedCablingTree > cabling_
unsigned int nBarrelModules() const
edm::EDGetTokenT< FEDRawDataCollection > tFEDRawDataCollection
std::vector< std::string > getParameterNames() const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< int > tkerrorlist
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
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
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
edm::ESWatcher< SiPixelQualityRcd > qualityWatcher
const SiPixelQuality * badPixelInfo_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)