CMS 3D CMS Logo

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 // Jan 2016 Tamas Almos Vami (Tav) (Wigner RCP) -- Cabling Map label option
7 // Jul 2017 Viktor Veszpremi -- added PixelFEDChannel
8 
9 #include "SiPixelRawToDigi.h"
10 
16 
18 
20 
22 
26 
32 
34 
38 
39 #include "TH1D.h"
40 #include "TFile.h"
41 
42 using namespace std;
43 
44 // -----------------------------------------------------------------------------
46  : config_(conf),
47  badPixelInfo_(nullptr),
48  regions_(nullptr),
49  hCPU(nullptr), hDigi(nullptr)
50 {
51 
52  includeErrors = config_.getParameter<bool>("IncludeErrors");
53  useQuality = config_.getParameter<bool>("UseQualityInfo");
54  if (config_.exists("ErrorList")) {
55  tkerrorlist = config_.getParameter<std::vector<int> > ("ErrorList");
56  }
57  if (config_.exists("UserErrorList")) {
58  usererrorlist = config_.getParameter<std::vector<int> > ("UserErrorList");
59  }
60  tFEDRawDataCollection = consumes <FEDRawDataCollection> (config_.getParameter<edm::InputTag>("InputLabel"));
61  theBadPixelFEDChannelsLabel = consumes<PixelFEDChannelCollection>(config_.getParameter<edm::InputTag>("BadPixelFEDChannelsInputLabel"));
62 
63  //start counters
64  ndigis = 0;
65  nwords = 0;
66 
67  // Products
68  produces< edm::DetSetVector<PixelDigi> >();
69  if(includeErrors){
70  produces< edm::DetSetVector<SiPixelRawDataError> >();
71  produces<DetIdCollection>();
72  produces<DetIdCollection>("UserErrorModules");
73  produces<edmNew::DetSetVector<PixelFEDChannel> >();
74  }
75 
76  // regions
77  if (config_.exists("Regions")) {
78  if(!config_.getParameter<edm::ParameterSet>("Regions").getParameterNames().empty())
79  {
80  regions_ = new PixelUnpackingRegions(config_, consumesCollector());
81  }
82  }
83 
84  // Timing
85  bool timing = config_.getUntrackedParameter<bool>("Timing",false);
86  if (timing) {
87  theTimer.reset( new edm::CPUTimer );
88  hCPU = new TH1D ("hCPU","hCPU",100,0.,0.050);
89  hDigi = new TH1D("hDigi","hDigi",50,0.,15000.);
90  }
91 
92  // Control the usage of pilot-blade data, FED=40
93  usePilotBlade = false;
94  if (config_.exists("UsePilotBlade")) {
95  usePilotBlade = config_.getParameter<bool> ("UsePilotBlade");
96  if(usePilotBlade) edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
97  }
98 
99  // Control the usage of phase1
100  usePhase1 = false;
101  if (config_.exists("UsePhase1")) {
102  usePhase1 = config_.getParameter<bool> ("UsePhase1");
103  if(usePhase1) edm::LogInfo("SiPixelRawToDigi") << " Using phase1";
104  }
105  //CablingMap could have a label //Tav
106  cablingMapLabel = config_.getParameter<std::string> ("CablingMapLabel");
107 
108 }
109 
110 
111 // -----------------------------------------------------------------------------
113  edm::LogInfo("SiPixelRawToDigi") << " HERE ** SiPixelRawToDigi destructor!";
114 
115  if (regions_) delete regions_;
116 
117  if (theTimer) {
118  TFile rootFile("analysis.root", "RECREATE", "my histograms");
119  hCPU->Write();
120  hDigi->Write();
121  }
122 
123 }
124 
125 void
128  desc.add<bool>("IncludeErrors",true);
129  desc.add<bool>("UseQualityInfo",false);
130  {
131  std::vector<int> temp1;
132  temp1.reserve(1);
133  temp1.push_back(29);
134  desc.add<std::vector<int> >("ErrorList",temp1)->setComment("## ErrorList: list of error codes used by tracking to invalidate modules");
135  }
136  {
137  std::vector<int> temp1;
138  temp1.reserve(1);
139  temp1.push_back(40);
140  desc.add<std::vector<int> >("UserErrorList",temp1)->setComment("## UserErrorList: list of error codes used by Pixel experts for investigation");
141  }
142  desc.add<edm::InputTag>("InputLabel",edm::InputTag("siPixelRawData"));
143  {
145  psd0.addOptional<std::vector<edm::InputTag>>("inputs");
146  psd0.addOptional<std::vector<double>>("deltaPhi");
147  psd0.addOptional<std::vector<double>>("maxZ");
148  psd0.addOptional<edm::InputTag>("beamSpot");
149  desc.add<edm::ParameterSetDescription>("Regions",psd0)->setComment("## Empty Regions PSet means complete unpacking");
150  }
151  desc.addUntracked<bool>("Timing",false);
152  desc.add<bool>("UsePilotBlade",false)->setComment("## Use pilot blades");
153  desc.add<bool>("UsePhase1",false)->setComment("## Use phase1");
154  desc.add<std::string>("CablingMapLabel","")->setComment("CablingMap label"); //Tav
155  desc.addOptional<bool>("CheckPixelOrder"); // never used, kept for back-compatibility
156  desc.add<edm::InputTag>("BadPixelFEDChannelsInputLabel",edm::InputTag("mix"));
157  descriptions.add("siPixelRawToDigi",desc);
158 
159 }
160 
161 // -----------------------------------------------------------------------------
162 
163 
164 // -----------------------------------------------------------------------------
166  const edm::EventSetup& es)
167 {
168  const uint32_t dummydetid = 0xffffffff;
170 
171 // initialize cabling map or update if necessary
172  if (recordWatcher.check( es )) {
173  // cabling map, which maps online address (fed->link->ROC->local pixel) to offline (DetId->global pixel)
175  es.get<SiPixelFedCablingMapRcd>().get( cablingMapLabel, cablingMap ); //Tav
176  fedIds = cablingMap->fedIds();
177  cabling_ = cablingMap->cablingTree();
178  LogDebug("map version:")<< cabling_->version();
179  }
180 // initialize quality record or update if necessary
181  if (qualityWatcher.check( es )&&useQuality) {
182  // quality info for dead pixel modules or ROCs
183  edm::ESHandle<SiPixelQuality> qualityInfo;
184  es.get<SiPixelQualityRcd>().get( qualityInfo );
185  badPixelInfo_ = qualityInfo.product();
186  if (!badPixelInfo_) {
187  edm::LogError("SiPixelQualityNotPresent")<<" Configured to use SiPixelQuality, but SiPixelQuality not present"<<endl;
188  }
189  }
190 
192  ev.getByToken(tFEDRawDataCollection, buffers);
193 
194 // create product (digis & errors)
195  auto collection = std::make_unique<edm::DetSetVector<PixelDigi>>();
196  // collection->reserve(8*1024);
197  auto errorcollection = std::make_unique<edm::DetSetVector<SiPixelRawDataError>>();
198  auto tkerror_detidcollection = std::make_unique<DetIdCollection>();
199  auto usererror_detidcollection = std::make_unique<DetIdCollection>();
200  auto disabled_channelcollection = std::make_unique<edmNew::DetSetVector<PixelFEDChannel> >();
201 
202  //PixelDataFormatter formatter(cabling_.get()); // phase 0 only
203  PixelDataFormatter formatter(cabling_.get(), usePhase1); // for phase 1 & 0
204 
205  formatter.setErrorStatus(includeErrors);
206 
207  if (useQuality) formatter.setQualityStatus(useQuality, badPixelInfo_);
208 
209  if (theTimer) theTimer->start();
210  bool errorsInEvent = false;
211  PixelDataFormatter::DetErrors nodeterrors;
212 
213  if (regions_) {
214  regions_->run(ev, es);
215  formatter.setModulesToUnpack(regions_->modulesToUnpack());
216  LogDebug("SiPixelRawToDigi") << "region2unpack #feds: "<<regions_->nFEDs();
217  LogDebug("SiPixelRawToDigi") << "region2unpack #modules (BPIX,EPIX,total): "<<regions_->nBarrelModules()<<" "<<regions_->nForwardModules()<<" "<<regions_->nModules();
218  }
219 
220  for (auto aFed = fedIds.begin(); aFed != fedIds.end(); ++aFed) {
221  int fedId = *aFed;
222 
223  if(!usePilotBlade && (fedId==40) ) continue; // skip pilot blade data
224 
225  if (regions_ && !regions_->mayUnpackFED(fedId)) continue;
226 
227  if(debug) LogDebug("SiPixelRawToDigi")<< " PRODUCE DIGI FOR FED: " << fedId << endl;
228 
230 
231  //get event data for this fed
232  const FEDRawData& fedRawData = buffers->FEDData( fedId );
233 
234  //convert data to digi and strip off errors
235  formatter.interpretRawData( errorsInEvent, fedId, fedRawData, *collection, errors);
236 
237  //pack errors into collection
238  if(includeErrors) {
239  typedef PixelDataFormatter::Errors::iterator IE;
240  for (IE is = errors.begin(); is != errors.end(); is++) {
241  uint32_t errordetid = is->first;
242  if (errordetid==dummydetid) { // errors given dummy detId must be sorted by Fed
243  nodeterrors.insert( nodeterrors.end(), errors[errordetid].begin(), errors[errordetid].end() );
244  } else {
245  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(errordetid);
246  errorDetSet.data.insert(errorDetSet.data.end(), is->second.begin(), is->second.end());
247  // Fill detid of the detectors where there is error AND the error number is listed
248  // in the configurable error list in the job option cfi.
249  // Code needs to be here, because there can be a set of errors for each
250  // entry in the for loop over PixelDataFormatter::Errors
251 
252  std::vector<PixelFEDChannel> disabledChannelsDetSet;
253 
254  for (auto const& aPixelError : errorDetSet) {
255  // For the time being, we extend the error handling functionality with ErrorType 25
256  // In the future, we should sort out how the usage of tkerrorlist can be generalized
257  if (aPixelError.getType()==25) {
258  assert(aPixelError.getFedId()==fedId);
259  const sipixelobjects::PixelFEDCabling* fed = cabling_->fed(fedId);
260  if (fed) {
261  cms_uint32_t linkId = formatter.linkId(aPixelError.getWord32());
262  const sipixelobjects::PixelFEDLink* link = fed->link(linkId);
263  if (link) {
264  // The "offline" 0..15 numbering is fixed by definition, also, the FrameConversion depends on it
265  // in contrast, the ROC-in-channel numbering is determined by hardware --> better to use the "offline" scheme
266  PixelFEDChannel ch = {fed->id(), linkId, 25, 0};
267  for (unsigned int iRoc=1; iRoc<=link->numberOfROCs(); iRoc++) {
268  const sipixelobjects::PixelROC * roc = link->roc(iRoc);
269  if (roc->idInDetUnit()<ch.roc_first) ch.roc_first=roc->idInDetUnit();
270  if (roc->idInDetUnit()>ch.roc_last) ch.roc_last=roc->idInDetUnit();
271  }
272  disabledChannelsDetSet.push_back(ch);
273  }
274  }
275  } else {
276  // fill list of detIds to be turned off by tracking
277  if(!tkerrorlist.empty()) {
278  std::vector<int>::iterator it_find = find(tkerrorlist.begin(), tkerrorlist.end(), aPixelError.getType());
279  if(it_find != tkerrorlist.end()){
280  tkerror_detidcollection->push_back(errordetid);
281  }
282  }
283  }
284 
285  // fill list of detIds with errors to be studied
286  if(!usererrorlist.empty()) {
287  std::vector<int>::iterator it_find = find(usererrorlist.begin(), usererrorlist.end(), aPixelError.getType());
288  if(it_find != usererrorlist.end()){
289  usererror_detidcollection->push_back(errordetid);
290  }
291  }
292 
293  } // loop on DetSet of errors
294 
295  if (!disabledChannelsDetSet.empty()) {
296  disabled_channelcollection->insert(errordetid, disabledChannelsDetSet.data(), disabledChannelsDetSet.size());
297  }
298  } // if error assigned to a real DetId
299  } // loop on errors in event for this FED
300  } // if errors to be included in the event
301  } // loop on FED data to be unpacked
302 
303  if(includeErrors) {
304  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(dummydetid);
305  errorDetSet.data = nodeterrors;
306  }
307  if (errorsInEvent) LogDebug("SiPixelRawToDigi") << "Error words were stored in this event";
308 
309  if (theTimer) {
310  theTimer->stop();
311  LogDebug("SiPixelRawToDigi") << "TIMING IS: (real)" << theTimer->realTime() ;
312  ndigis += formatter.nDigis();
313  nwords += formatter.nWords();
314  LogDebug("SiPixelRawToDigi") << " (Words/Digis) this ev: "
315  <<formatter.nWords()<<"/"<<formatter.nDigis() << "--- all :"<<nwords<<"/"<<ndigis;
316  hCPU->Fill( theTimer->realTime() );
317  hDigi->Fill(formatter.nDigis());
318  }
319 
320 
321  //------------------------------------
322  //send digis and errors back to framework
323 
324  edm::Handle<PixelFEDChannelCollection> pixelFEDChannelCollectionHandle;
325  std::unique_ptr<PixelFEDChannelCollection> PixelFEDChannelCollection_ = nullptr;
326  if (ev.getByToken(theBadPixelFEDChannelsLabel, pixelFEDChannelCollectionHandle)){
327 
328  const PixelFEDChannelCollection * pfcc= pixelFEDChannelCollectionHandle.product();
329  PixelFEDChannelCollection_ = std::make_unique<PixelFEDChannelCollection>(*pfcc);
330  }
331 
332  ev.put(std::move(collection));
333  if(includeErrors){
334  ev.put(std::move(errorcollection));
335  ev.put(std::move(tkerror_detidcollection));
336  ev.put(std::move(usererror_detidcollection), "UserErrorModules");
337  if (PixelFEDChannelCollection_ == nullptr){
338  ev.put(std::move(disabled_channelcollection));
339  }
340  else{
341  ev.put(std::move(PixelFEDChannelCollection_));
342  }
343  }
344 }
345 // declare this as a framework plugin
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< edm::CPUTimer > theTimer
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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
std::string cablingMapLabel
void setErrorStatus(bool ErrorStatus)
void produce(edm::Event &, const edm::EventSetup &) override
get data, convert to digis attach againe to Event
#define nullptr
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
PixelUnpackingRegions * regions_
unsigned int roc_last
std::vector< SiPixelRawDataError > DetErrors
std::vector< unsigned int > fedIds
SiPixelRawToDigi(const edm::ParameterSet &)
ctor
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
unsigned int idInDetUnit() const
id of this ROC in DetUnit etermined by token path
Definition: PixelROC.h:40
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)
unsigned int cms_uint32_t
Definition: typedefs.h:15
~SiPixelRawToDigi() override
dtor
std::map< cms_uint32_t, DetErrors > Errors
edm::ESWatcher< SiPixelFedCablingMapRcd > recordWatcher
void run(const edm::Event &e, const edm::EventSetup &es)
has to be run during each event
std::unique_ptr< SiPixelFedCablingTree > cablingTree() const
unsigned int nModules() const
T const * product() const
Definition: Handle.h:81
unsigned int roc_first
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:
T get() const
Definition: EventSetup.h:68
bool mayUnpackFED(unsigned int fed_n) const
check whether a FED has to be unpacked
edm::EDGetTokenT< PixelFEDChannelCollection > theBadPixelFEDChannelsLabel
std::vector< unsigned int > fedIds() const
edm::ESWatcher< SiPixelQualityRcd > qualityWatcher
T const * product() const
Definition: ESHandle.h:84
def move(src, dest)
Definition: eostools.py:511
const SiPixelQuality * badPixelInfo_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)