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 
12 
14 
16 
18 
22 
28 
30 
33 
34 #include "TH1D.h"
35 #include "TFile.h"
36 
37 using namespace std;
38 
39 // -----------------------------------------------------------------------------
41  : config_(conf),
42  badPixelInfo_(0),
43  regions_(0),
44  hCPU(0), hDigi(0)
45 {
46 
47  includeErrors = config_.getParameter<bool>("IncludeErrors");
48  useQuality = config_.getParameter<bool>("UseQualityInfo");
49  if (config_.exists("ErrorList")) {
50  tkerrorlist = config_.getParameter<std::vector<int> > ("ErrorList");
51  }
52  if (config_.exists("UserErrorList")) {
53  usererrorlist = config_.getParameter<std::vector<int> > ("UserErrorList");
54  }
55  tFEDRawDataCollection = consumes <FEDRawDataCollection> (config_.getParameter<edm::InputTag>("InputLabel"));
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.reset( new edm::CPUTimer );
81  hCPU = new TH1D ("hCPU","hCPU",100,0.,0.050);
82  hDigi = new TH1D("hDigi","hDigi",50,0.,15000.);
83  }
84 
85  // Control the usage of pilot-blade data, FED=40
86  usePilotBlade = false;
87  if (config_.exists("UsePilotBlade")) {
88  usePilotBlade = config_.getParameter<bool> ("UsePilotBlade");
89  if(usePilotBlade) edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
90  }
91 
92  // Control the usage of phase1
93  usePhase1 = false;
94  if (config_.exists("UsePhase1")) {
95  usePhase1 = config_.getParameter<bool> ("UsePhase1");
96  if(usePhase1) edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
97  }
98 
99 }
100 
101 
102 // -----------------------------------------------------------------------------
104  edm::LogInfo("SiPixelRawToDigi") << " HERE ** SiPixelRawToDigi destructor!";
105 
106  if (regions_) delete regions_;
107 
108  if (theTimer) {
109  TFile rootFile("analysis.root", "RECREATE", "my histograms");
110  hCPU->Write();
111  hDigi->Write();
112  }
113 
114 }
115 
116 
117 // -----------------------------------------------------------------------------
118 
119 
120 // -----------------------------------------------------------------------------
122  const edm::EventSetup& es)
123 {
124  const uint32_t dummydetid = 0xffffffff;
126 
127 // initialize cabling map or update if necessary
128  if (recordWatcher.check( es )) {
129  // cabling map, which maps online address (fed->link->ROC->local pixel) to offline (DetId->global pixel)
131  es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
132  fedIds = cablingMap->fedIds();
133  cabling_ = cablingMap->cablingTree();
134  LogDebug("map version:")<< cabling_->version();
135  }
136 // initialize quality record or update if necessary
137  if (qualityWatcher.check( es )&&useQuality) {
138  // quality info for dead pixel modules or ROCs
139  edm::ESHandle<SiPixelQuality> qualityInfo;
140  es.get<SiPixelQualityRcd>().get( qualityInfo );
141  badPixelInfo_ = qualityInfo.product();
142  if (!badPixelInfo_) {
143  edm::LogError("SiPixelQualityNotPresent")<<" Configured to use SiPixelQuality, but SiPixelQuality not present"<<endl;
144  }
145  }
146 
148  ev.getByToken(tFEDRawDataCollection, buffers);
149 
150 // create product (digis & errors)
151  std::auto_ptr< edm::DetSetVector<PixelDigi> > collection( new edm::DetSetVector<PixelDigi> );
152  // collection->reserve(8*1024);
153  std::auto_ptr< edm::DetSetVector<SiPixelRawDataError> > errorcollection( new edm::DetSetVector<SiPixelRawDataError> );
154  std::auto_ptr< DetIdCollection > tkerror_detidcollection(new DetIdCollection());
155  std::auto_ptr< DetIdCollection > usererror_detidcollection(new DetIdCollection());
156 
157  //PixelDataFormatter formatter(cabling_.get()); // phase 0 only
158  PixelDataFormatter formatter(cabling_.get(), usePhase1); // for phase 1 & 0
159 
160  formatter.setErrorStatus(includeErrors);
161 
162  if (useQuality) formatter.setQualityStatus(useQuality, badPixelInfo_);
163 
164  if (theTimer) theTimer->start();
165  bool errorsInEvent = false;
166  PixelDataFormatter::DetErrors nodeterrors;
167 
168  if (regions_) {
169  regions_->run(ev, es);
170  formatter.setModulesToUnpack(regions_->modulesToUnpack());
171  LogDebug("SiPixelRawToDigi") << "region2unpack #feds (BPIX,EPIX,total): "<<regions_->nBarrelFEDs()<<" "<<regions_->nForwardFEDs()<<" "<<regions_->nFEDs();
172  LogDebug("SiPixelRawToDigi") << "region2unpack #modules (BPIX,EPIX,total): "<<regions_->nBarrelModules()<<" "<<regions_->nForwardModules()<<" "<<regions_->nModules();
173  }
174 
175  for (auto aFed = fedIds.begin(); aFed != fedIds.end(); ++aFed) {
176  int fedId = *aFed;
177 
178  if(!usePilotBlade && (fedId==40) ) continue; // skip pilot blade data
179 
180  if (regions_ && !regions_->mayUnpackFED(fedId)) continue;
181 
182  if(debug) LogDebug("SiPixelRawToDigi")<< " PRODUCE DIGI FOR FED: " << fedId << endl;
183 
185 
186  //get event data for this fed
187  const FEDRawData& fedRawData = buffers->FEDData( fedId );
188 
189  //convert data to digi and strip off errors
190  formatter.interpretRawData( errorsInEvent, fedId, fedRawData, *collection, errors);
191 
192  //pack errors into collection
193  if(includeErrors) {
194  typedef PixelDataFormatter::Errors::iterator IE;
195  for (IE is = errors.begin(); is != errors.end(); is++) {
196  uint32_t errordetid = is->first;
197  if (errordetid==dummydetid) { // errors given dummy detId must be sorted by Fed
198  nodeterrors.insert( nodeterrors.end(), errors[errordetid].begin(), errors[errordetid].end() );
199  } else {
200  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(errordetid);
201  errorDetSet.data.insert(errorDetSet.data.end(), is->second.begin(), is->second.end());
202  // Fill detid of the detectors where there is error AND the error number is listed
203  // in the configurable error list in the job option cfi.
204  // Code needs to be here, because there can be a set of errors for each
205  // entry in the for loop over PixelDataFormatter::Errors
206  if(!tkerrorlist.empty() || !usererrorlist.empty()){
207  DetId errorDetId(errordetid);
208  edm::DetSet<SiPixelRawDataError>::const_iterator itPixelError=errorDetSet.begin();
209  for(; itPixelError!=errorDetSet.end(); ++itPixelError){
210  // fill list of detIds to be turned off by tracking
211  if(!tkerrorlist.empty()) {
212  std::vector<int>::iterator it_find = find(tkerrorlist.begin(), tkerrorlist.end(), itPixelError->getType());
213  if(it_find != tkerrorlist.end()){
214  tkerror_detidcollection->push_back(errordetid);
215  }
216  }
217  // fill list of detIds with errors to be studied
218  if(!usererrorlist.empty()) {
219  std::vector<int>::iterator it_find = find(usererrorlist.begin(), usererrorlist.end(), itPixelError->getType());
220  if(it_find != usererrorlist.end()){
221  usererror_detidcollection->push_back(errordetid);
222  }
223  }
224  }
225  }
226  }
227  }
228  }
229  }
230 
231  if(includeErrors) {
232  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(dummydetid);
233  errorDetSet.data = nodeterrors;
234  }
235  if (errorsInEvent) LogDebug("SiPixelRawToDigi") << "Error words were stored in this event";
236 
237  if (theTimer) {
238  theTimer->stop();
239  LogDebug("SiPixelRawToDigi") << "TIMING IS: (real)" << theTimer->realTime() ;
240  ndigis += formatter.nDigis();
241  nwords += formatter.nWords();
242  LogDebug("SiPixelRawToDigi") << " (Words/Digis) this ev: "
243  <<formatter.nWords()<<"/"<<formatter.nDigis() << "--- all :"<<nwords<<"/"<<ndigis;
244  hCPU->Fill( theTimer->realTime() );
245  hDigi->Fill(formatter.nDigis());
246  }
247 
248  //send digis and errors back to framework
249  ev.put( collection );
250  if(includeErrors){
251  ev.put( errorcollection );
252  ev.put( tkerror_detidcollection );
253  ev.put( usererror_detidcollection, "UserErrorModules" );
254  }
255 }
#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
unsigned int nForwardFEDs() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
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:113
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
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:55
T const * product() const
Definition: ESHandle.h:86
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_