test
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 
8 #include "SiPixelRawToDigi.h"
9 
15 
17 
19 
21 
25 
31 
33 
36 
37 #include "TH1D.h"
38 #include "TFile.h"
39 
40 using namespace std;
41 
42 // -----------------------------------------------------------------------------
44  : config_(conf),
45  badPixelInfo_(0),
46  regions_(0),
47  hCPU(0), hDigi(0)
48 {
49 
50  includeErrors = config_.getParameter<bool>("IncludeErrors");
51  useQuality = config_.getParameter<bool>("UseQualityInfo");
52  if (config_.exists("ErrorList")) {
53  tkerrorlist = config_.getParameter<std::vector<int> > ("ErrorList");
54  }
55  if (config_.exists("UserErrorList")) {
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  }
71 
72  // regions
73  if (config_.exists("Regions")) {
74  if(config_.getParameter<edm::ParameterSet>("Regions").getParameterNames().size() > 0)
75  {
77  }
78  }
79 
80  // Timing
81  bool timing = config_.getUntrackedParameter<bool>("Timing",false);
82  if (timing) {
83  theTimer.reset( new edm::CPUTimer );
84  hCPU = new TH1D ("hCPU","hCPU",100,0.,0.050);
85  hDigi = new TH1D("hDigi","hDigi",50,0.,15000.);
86  }
87 
88  // Control the usage of pilot-blade data, FED=40
89  usePilotBlade = false;
90  if (config_.exists("UsePilotBlade")) {
91  usePilotBlade = config_.getParameter<bool> ("UsePilotBlade");
92  if(usePilotBlade) edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
93  }
94 
95  // Control the usage of phase1
96  usePhase1 = false;
97  if (config_.exists("UsePhase1")) {
98  usePhase1 = config_.getParameter<bool> ("UsePhase1");
99  if(usePhase1) edm::LogInfo("SiPixelRawToDigi") << " Use pilot blade data (FED 40)";
100  }
101  //CablingMap could have a label //Tav
102  cablingMapLabel = config_.getParameter<std::string> ("CablingMapLabel");
103 
104 }
105 
106 
107 // -----------------------------------------------------------------------------
109  edm::LogInfo("SiPixelRawToDigi") << " HERE ** SiPixelRawToDigi destructor!";
110 
111  if (regions_) delete regions_;
112 
113  if (theTimer) {
114  TFile rootFile("analysis.root", "RECREATE", "my histograms");
115  hCPU->Write();
116  hDigi->Write();
117  }
118 
119 }
120 
121 void
124  desc.add<bool>("IncludeErrors",true);
125  desc.add<bool>("UseQualityInfo",false);
126  {
127  std::vector<int> temp1;
128  temp1.reserve(1);
129  temp1.push_back(29);
130  desc.add<std::vector<int> >("ErrorList",temp1)->setComment("## ErrorList: list of error codes used by tracking to invalidate modules");
131  }
132  {
133  std::vector<int> temp1;
134  temp1.reserve(1);
135  temp1.push_back(40);
136  desc.add<std::vector<int> >("UserErrorList",temp1)->setComment("## UserErrorList: list of error codes used by Pixel experts for investigation");
137  }
138  desc.add<edm::InputTag>("InputLabel",edm::InputTag("siPixelRawData"));
139  {
141  psd0.addOptional<std::vector<edm::InputTag>>("inputs");
142  psd0.addOptional<std::vector<double>>("deltaPhi");
143  psd0.addOptional<std::vector<double>>("maxZ");
144  psd0.addOptional<edm::InputTag>("beamSpot");
145  desc.add<edm::ParameterSetDescription>("Regions",psd0)->setComment("## Empty Regions PSet means complete unpacking");
146  }
147  desc.addUntracked<bool>("Timing",false);
148  desc.add<bool>("UsePilotBlade",false)->setComment("## Use pilot blades");
149  desc.add<bool>("UsePhase1",false)->setComment("## Use phase1");
150  desc.add<std::string>("CablingMapLabel","")->setComment("CablingMap label"); //Tav
151  desc.addOptional<bool>("CheckPixelOrder"); // never used, kept for back-compatibility
152  descriptions.add("siPixelRawToDigi",desc);
153 }
154 
155 // -----------------------------------------------------------------------------
156 
157 
158 // -----------------------------------------------------------------------------
160  const edm::EventSetup& es)
161 {
162  const uint32_t dummydetid = 0xffffffff;
164 
165 // initialize cabling map or update if necessary
166  if (recordWatcher.check( es )) {
167  // cabling map, which maps online address (fed->link->ROC->local pixel) to offline (DetId->global pixel)
169  es.get<SiPixelFedCablingMapRcd>().get( cablingMapLabel, cablingMap ); //Tav
170  fedIds = cablingMap->fedIds();
171  cabling_ = cablingMap->cablingTree();
172  LogDebug("map version:")<< cabling_->version();
173  }
174 // initialize quality record or update if necessary
175  if (qualityWatcher.check( es )&&useQuality) {
176  // quality info for dead pixel modules or ROCs
177  edm::ESHandle<SiPixelQuality> qualityInfo;
178  es.get<SiPixelQualityRcd>().get( qualityInfo );
179  badPixelInfo_ = qualityInfo.product();
180  if (!badPixelInfo_) {
181  edm::LogError("SiPixelQualityNotPresent")<<" Configured to use SiPixelQuality, but SiPixelQuality not present"<<endl;
182  }
183  }
184 
186  ev.getByToken(tFEDRawDataCollection, buffers);
187 
188 // create product (digis & errors)
189  std::auto_ptr< edm::DetSetVector<PixelDigi> > collection( new edm::DetSetVector<PixelDigi> );
190  // collection->reserve(8*1024);
191  std::auto_ptr< edm::DetSetVector<SiPixelRawDataError> > errorcollection( new edm::DetSetVector<SiPixelRawDataError> );
192  std::auto_ptr< DetIdCollection > tkerror_detidcollection(new DetIdCollection());
193  std::auto_ptr< DetIdCollection > usererror_detidcollection(new DetIdCollection());
194 
195  //PixelDataFormatter formatter(cabling_.get()); // phase 0 only
196  PixelDataFormatter formatter(cabling_.get(), usePhase1); // for phase 1 & 0
197 
198  formatter.setErrorStatus(includeErrors);
199 
200  if (useQuality) formatter.setQualityStatus(useQuality, badPixelInfo_);
201 
202  if (theTimer) theTimer->start();
203  bool errorsInEvent = false;
204  PixelDataFormatter::DetErrors nodeterrors;
205 
206  if (regions_) {
207  regions_->run(ev, es);
208  formatter.setModulesToUnpack(regions_->modulesToUnpack());
209  LogDebug("SiPixelRawToDigi") << "region2unpack #feds (BPIX,EPIX,total): "<<regions_->nBarrelFEDs()<<" "<<regions_->nForwardFEDs()<<" "<<regions_->nFEDs();
210  LogDebug("SiPixelRawToDigi") << "region2unpack #modules (BPIX,EPIX,total): "<<regions_->nBarrelModules()<<" "<<regions_->nForwardModules()<<" "<<regions_->nModules();
211  }
212 
213  for (auto aFed = fedIds.begin(); aFed != fedIds.end(); ++aFed) {
214  int fedId = *aFed;
215 
216  if(!usePilotBlade && (fedId==40) ) continue; // skip pilot blade data
217 
218  if (regions_ && !regions_->mayUnpackFED(fedId)) continue;
219 
220  if(debug) LogDebug("SiPixelRawToDigi")<< " PRODUCE DIGI FOR FED: " << fedId << endl;
221 
223 
224  //get event data for this fed
225  const FEDRawData& fedRawData = buffers->FEDData( fedId );
226 
227  //convert data to digi and strip off errors
228  formatter.interpretRawData( errorsInEvent, fedId, fedRawData, *collection, errors);
229 
230  //pack errors into collection
231  if(includeErrors) {
232  typedef PixelDataFormatter::Errors::iterator IE;
233  for (IE is = errors.begin(); is != errors.end(); is++) {
234  uint32_t errordetid = is->first;
235  if (errordetid==dummydetid) { // errors given dummy detId must be sorted by Fed
236  nodeterrors.insert( nodeterrors.end(), errors[errordetid].begin(), errors[errordetid].end() );
237  } else {
238  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(errordetid);
239  errorDetSet.data.insert(errorDetSet.data.end(), is->second.begin(), is->second.end());
240  // Fill detid of the detectors where there is error AND the error number is listed
241  // in the configurable error list in the job option cfi.
242  // Code needs to be here, because there can be a set of errors for each
243  // entry in the for loop over PixelDataFormatter::Errors
244  if(!tkerrorlist.empty() || !usererrorlist.empty()){
245  DetId errorDetId(errordetid);
246  edm::DetSet<SiPixelRawDataError>::const_iterator itPixelError=errorDetSet.begin();
247  for(; itPixelError!=errorDetSet.end(); ++itPixelError){
248  // fill list of detIds to be turned off by tracking
249  if(!tkerrorlist.empty()) {
250  std::vector<int>::iterator it_find = find(tkerrorlist.begin(), tkerrorlist.end(), itPixelError->getType());
251  if(it_find != tkerrorlist.end()){
252  tkerror_detidcollection->push_back(errordetid);
253  }
254  }
255  // fill list of detIds with errors to be studied
256  if(!usererrorlist.empty()) {
257  std::vector<int>::iterator it_find = find(usererrorlist.begin(), usererrorlist.end(), itPixelError->getType());
258  if(it_find != usererrorlist.end()){
259  usererror_detidcollection->push_back(errordetid);
260  }
261  }
262  }
263  }
264  }
265  }
266  }
267  }
268 
269  if(includeErrors) {
270  edm::DetSet<SiPixelRawDataError>& errorDetSet = errorcollection->find_or_insert(dummydetid);
271  errorDetSet.data = nodeterrors;
272  }
273  if (errorsInEvent) LogDebug("SiPixelRawToDigi") << "Error words were stored in this event";
274 
275  if (theTimer) {
276  theTimer->stop();
277  LogDebug("SiPixelRawToDigi") << "TIMING IS: (real)" << theTimer->realTime() ;
278  ndigis += formatter.nDigis();
279  nwords += formatter.nWords();
280  LogDebug("SiPixelRawToDigi") << " (Words/Digis) this ev: "
281  <<formatter.nWords()<<"/"<<formatter.nDigis() << "--- all :"<<nwords<<"/"<<ndigis;
282  hCPU->Fill( theTimer->realTime() );
283  hDigi->Fill(formatter.nDigis());
284  }
285 
286  //send digis and errors back to framework
287  ev.put( collection );
288  if(includeErrors){
289  ev.put( errorcollection );
290  ev.put( tkerror_detidcollection );
291  ev.put( usererror_detidcollection, "UserErrorModules" );
292  }
293 }
#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:462
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)
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:121
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)
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)