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