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 // 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), badPixelInfo_(nullptr), regions_(nullptr), hCPU(nullptr), hDigi(nullptr) {
47  includeErrors = config_.getParameter<bool>("IncludeErrors");
48  useQuality = config_.getParameter<bool>("UseQualityInfo");
49 
50  tkerrorlist = config_.getParameter<std::vector<int>>("ErrorList");
51  usererrorlist = config_.getParameter<std::vector<int>>("UserErrorList");
52 
53  tFEDRawDataCollection = consumes<FEDRawDataCollection>(config_.getParameter<edm::InputTag>("InputLabel"));
54 
55  //start counters
56  ndigis = 0;
57  nwords = 0;
58 
59  // Products
60  produces<edm::DetSetVector<PixelDigi>>();
61  if (includeErrors) {
62  produces<edm::DetSetVector<SiPixelRawDataError>>();
63  produces<DetIdCollection>();
64  produces<DetIdCollection>("UserErrorModules");
65  produces<edmNew::DetSetVector<PixelFEDChannel>>();
66  }
67 
68  // regions
69  if (!config_.getParameter<edm::ParameterSet>("Regions").getParameterNames().empty()) {
70  regions_ = new PixelUnpackingRegions(config_, consumesCollector());
71  }
72 
73  // Timing
74  bool timing = config_.getUntrackedParameter<bool>("Timing", false);
75  if (timing) {
76  theTimer.reset(new edm::CPUTimer);
77  hCPU = new TH1D("hCPU", "hCPU", 100, 0., 0.050);
78  hDigi = new TH1D("hDigi", "hDigi", 50, 0., 15000.);
79  }
80 
81  // Control the usage of pilot-blade data, FED=40
82  usePilotBlade = config_.getParameter<bool>("UsePilotBlade");
83  if (usePilotBlade)
84  edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
85 
86  // Control the usage of phase1
87  usePhase1 = config_.getParameter<bool>("UsePhase1");
88  if (usePhase1)
89  edm::LogInfo("SiPixelRawToDigi") << " Using phase1";
90 
91  //CablingMap could have a label //Tav
92  cablingMapLabel = config_.getParameter<std::string>("CablingMapLabel");
93 }
94 
95 // -----------------------------------------------------------------------------
97  edm::LogInfo("SiPixelRawToDigi") << " HERE ** SiPixelRawToDigi destructor!";
98 
99  if (regions_)
100  delete regions_;
101 
102  if (theTimer) {
103  TFile rootFile("analysis.root", "RECREATE", "my histograms");
104  hCPU->Write();
105  hDigi->Write();
106  }
107 }
108 
111  desc.add<bool>("IncludeErrors", true);
112  desc.add<bool>("UseQualityInfo", false);
113  {
114  std::vector<int> temp1;
115  temp1.reserve(1);
116  temp1.push_back(29);
117  desc.add<std::vector<int>>("ErrorList", temp1)
118  ->setComment("## ErrorList: list of error codes used by tracking to invalidate modules");
119  }
120  {
121  std::vector<int> temp1;
122  temp1.reserve(1);
123  temp1.push_back(40);
124  desc.add<std::vector<int>>("UserErrorList", temp1)
125  ->setComment("## UserErrorList: list of error codes used by Pixel experts for investigation");
126  }
127  desc.add<edm::InputTag>("InputLabel", edm::InputTag("siPixelRawData"));
128  {
130  psd0.addOptional<std::vector<edm::InputTag>>("inputs");
131  psd0.addOptional<std::vector<double>>("deltaPhi");
132  psd0.addOptional<std::vector<double>>("maxZ");
133  psd0.addOptional<edm::InputTag>("beamSpot");
134  desc.add<edm::ParameterSetDescription>("Regions", psd0)
135  ->setComment("## Empty Regions PSet means complete unpacking");
136  }
137  desc.addUntracked<bool>("Timing", false);
138  desc.add<bool>("UsePilotBlade", false)->setComment("## Use pilot blades");
139  desc.add<bool>("UsePhase1", false)->setComment("## Use phase1");
140  desc.add<std::string>("CablingMapLabel", "")->setComment("CablingMap label"); //Tav
141  desc.addOptional<bool>("CheckPixelOrder"); // never used, kept for back-compatibility
142  descriptions.add("siPixelRawToDigi", desc);
143 }
144 
145 // -----------------------------------------------------------------------------
146 
147 // -----------------------------------------------------------------------------
149  const uint32_t dummydetid = 0xffffffff;
151 
152  // initialize cabling map or update if necessary
153  if (recordWatcher.check(es)) {
154  // cabling map, which maps online address (fed->link->ROC->local pixel) to offline (DetId->global pixel)
156  es.get<SiPixelFedCablingMapRcd>().get(cablingMapLabel, cablingMap); //Tav
157  fedIds = cablingMap->fedIds();
158  cabling_ = cablingMap->cablingTree();
159  LogDebug("map version:") << cabling_->version();
160  }
161  // initialize quality record or update if necessary
162  if (qualityWatcher.check(es) && useQuality) {
163  // quality info for dead pixel modules or ROCs
164  edm::ESHandle<SiPixelQuality> qualityInfo;
165  es.get<SiPixelQualityRcd>().get(qualityInfo);
166  badPixelInfo_ = qualityInfo.product();
167  if (!badPixelInfo_) {
168  edm::LogError("SiPixelQualityNotPresent")
169  << " Configured to use SiPixelQuality, but SiPixelQuality not present" << endl;
170  }
171  }
172 
174  ev.getByToken(tFEDRawDataCollection, buffers);
175 
176  // create product (digis & errors)
177  auto collection = std::make_unique<edm::DetSetVector<PixelDigi>>();
178  // collection->reserve(8*1024);
179  auto errorcollection = std::make_unique<edm::DetSetVector<SiPixelRawDataError>>();
180  auto tkerror_detidcollection = std::make_unique<DetIdCollection>();
181  auto usererror_detidcollection = std::make_unique<DetIdCollection>();
182  auto disabled_channelcollection = std::make_unique<edmNew::DetSetVector<PixelFEDChannel>>();
183 
184  PixelDataFormatter formatter(cabling_.get(), usePhase1); // for phase 1 & 0
185 
186  formatter.setErrorStatus(includeErrors);
187 
188  if (useQuality)
189  formatter.setQualityStatus(useQuality, badPixelInfo_);
190 
191  if (theTimer)
192  theTimer->start();
193  bool errorsInEvent = false;
194  PixelDataFormatter::DetErrors nodeterrors;
195 
196  if (regions_) {
197  regions_->run(ev, es);
198  formatter.setModulesToUnpack(regions_->modulesToUnpack());
199  LogDebug("SiPixelRawToDigi") << "region2unpack #feds: " << regions_->nFEDs();
200  LogDebug("SiPixelRawToDigi") << "region2unpack #modules (BPIX,EPIX,total): " << regions_->nBarrelModules() << " "
201  << regions_->nForwardModules() << " " << regions_->nModules();
202  }
203 
204  for (auto aFed = fedIds.begin(); aFed != fedIds.end(); ++aFed) {
205  int fedId = *aFed;
206 
207  if (!usePilotBlade && (fedId == 40))
208  continue; // skip pilot blade data
209 
210  if (regions_ && !regions_->mayUnpackFED(fedId))
211  continue;
212 
213  if (debug)
214  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)
257  ch.roc_first = roc->idInDetUnit();
258  if (roc->idInDetUnit() > ch.roc_last)
259  ch.roc_last = roc->idInDetUnit();
260  }
261  disabledChannelsDetSet.push_back(ch);
262  }
263  }
264  } else {
265  // fill list of detIds to be turned off by tracking
266  if (!tkerrorlist.empty()) {
267  std::vector<int>::iterator it_find =
268  find(tkerrorlist.begin(), tkerrorlist.end(), aPixelError.getType());
269  if (it_find != tkerrorlist.end()) {
270  tkerror_detidcollection->push_back(errordetid);
271  }
272  }
273  }
274 
275  // fill list of detIds with errors to be studied
276  if (!usererrorlist.empty()) {
277  std::vector<int>::iterator it_find =
278  find(usererrorlist.begin(), usererrorlist.end(), aPixelError.getType());
279  if (it_find != usererrorlist.end()) {
280  usererror_detidcollection->push_back(errordetid);
281  }
282  }
283 
284  } // loop on DetSet of errors
285 
286  if (!disabledChannelsDetSet.empty()) {
287  disabled_channelcollection->insert(
288  errordetid, disabledChannelsDetSet.data(), disabledChannelsDetSet.size());
289  }
290  } // if error assigned to a real DetId
291  } // loop on errors in event for this FED
292  } // if errors to be included in the event
293  } // loop on FED data to be unpacked
294 
295  if (includeErrors) {
296  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(dummydetid);
297  errorDetSet.data = nodeterrors;
298  }
299  if (errorsInEvent)
300  LogDebug("SiPixelRawToDigi") << "Error words were stored in this event";
301 
302  if (theTimer) {
303  theTimer->stop();
304  LogDebug("SiPixelRawToDigi") << "TIMING IS: (real)" << theTimer->realTime();
305  ndigis += formatter.nDigis();
306  nwords += formatter.nWords();
307  LogDebug("SiPixelRawToDigi") << " (Words/Digis) this ev: " << formatter.nWords() << "/" << formatter.nDigis()
308  << "--- all :" << nwords << "/" << ndigis;
309  hCPU->Fill(theTimer->realTime());
310  hDigi->Fill(formatter.nDigis());
311  }
312 
313  ev.put(std::move(collection));
314  if (includeErrors) {
315  ev.put(std::move(errorcollection));
316  ev.put(std::move(tkerror_detidcollection));
317  ev.put(std::move(usererror_detidcollection), "UserErrorModules");
318  ev.put(std::move(disabled_channelcollection));
319  }
320 }
321 // 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:131
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:525
#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:19
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:37
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:81
unsigned int nFEDs() const
various informational accessors:
T get() const
Definition: EventSetup.h:73
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)