CMS 3D CMS Logo

SiPixelStatusProducer.cc
Go to the documentation of this file.
1 // C++ standard
2 #include <string>
3 // ROOT
4 #include "TMath.h"
5 // CMSSW FW
9 // CMSSW DataFormats
17 // "FED error 25"
19 // CMSSW CondFormats
26 
27 // Header file
29 
31  //NOTE: Token for all stream replicas are identical and constructors for the replicas are called
32  // sequentially so there is no race condition.
33  iCache->trackerGeometryToken_ = esConsumes<edm::Transition::BeginRun>();
34  iCache->trackerTopologyToken_ = esConsumes<edm::Transition::BeginRun>();
35  iCache->siPixelFedCablingMapToken_ = esConsumes<edm::Transition::BeginRun>();
36 
37  /* badPixelFEDChannelCollections */
38  std::vector<edm::InputTag> badPixelFEDChannelCollectionLabels =
39  iConfig.getParameter<edm::ParameterSet>("SiPixelStatusProducerParameters")
40  .getParameter<std::vector<edm::InputTag>>("badPixelFEDChannelCollections");
42  theBadPixelFEDChannelsTokens_.push_back(consumes<PixelFEDChannelCollection>(t));
43 
44  /* pixel clusters */
45  fPixelClusterLabel_ = iConfig.getParameter<edm::ParameterSet>("SiPixelStatusProducerParameters")
46  .getUntrackedParameter<edm::InputTag>("pixelClusterLabel");
47  fSiPixelClusterToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(fPixelClusterLabel_);
48 
49  //debug_ = iConfig.getUntrackedParameter<bool>("debug");
50 
51  /* register products */
52  produces<SiPixelDetectorStatus, edm::Transition::EndLuminosityBlock>("siPixelStatus");
53 }
54 
56 
57 //--------------------------------------------------------------------------------------------------
58 
59 std::shared_ptr<SiPixelTopoFinder> SiPixelStatusProducer::globalBeginRun(edm::Run const& iRun,
60  edm::EventSetup const& iSetup,
61  GlobalCache const* iCache) {
62  const TrackerGeometry* trackerGeometry = &iSetup.getData(iCache->trackerGeometryToken_);
63  const TrackerTopology* trackerTopology = &iSetup.getData(iCache->trackerTopologyToken_);
64  const SiPixelFedCablingMap* cablingMap = &iSetup.getData(iCache->siPixelFedCablingMapToken_);
65 
66  auto returnValue = std::make_shared<SiPixelTopoFinder>();
67 
68  returnValue->init(trackerGeometry, trackerTopology, cablingMap);
69  return returnValue;
70 }
71 
73  /*Is it good to pass the objects stored in runCache to set class private members values *
74  or just call runCahche every time in the calss function?*/
75 
76  edm::LogInfo("SiPixelStatusProducer") << "beginRun: update the std::map for pixel geo/topo " << std::endl;
77  /* update the std::map for pixel geo/topo */
78  /* vector of all <int> detIds */
79  fDetIds_ = runCache()->getDetIds(); //getDetIds();
80  /* ROC size (number of row, number of columns for each det id) */
81  fSensors_ = runCache()->getSensors();
82  /* the roc layout on a module */
83  fSensorLayout_ = runCache()->getSensorLayout();
84  /* fedId as a function of detId */
85  fFedIds_ = runCache()->getFedIds();
86  /* map the index ROC to rocId */
87  fRocIds_ = runCache()->getRocIds();
88 }
89 
91  edm::LogInfo("SiPixelStatusProducer") << "beginlumi instance" << std::endl;
92 
93  /* initialize fDet_ with a set of modules(detIds) and clean the fFEDerror25_ */
95  for (unsigned int itDetId = 0; itDetId < fDetIds_.size(); ++itDetId) {
96  int detid = fDetIds_[itDetId];
97  int nrocs = fSensorLayout_[detid].first * fSensorLayout_[detid].second;
98 
99  fDet_.addModule(detid, nrocs);
100  }
101 
102  fFEDerror25_.clear();
103  ftotalevents_ = 0;
104 }
105 
107  edm::LogInfo("SiPixelStatusProducer") << "start cluster analyzer " << std::endl;
108 
109  /* count number of events for the current module instance in the luminosityBlock */
110  ftotalevents_++;
111 
112  /* ----------------------------------------------------------------------
113  -- Pixel cluster analysis
114  ----------------------------------------------------------------------*/
115 
117  if (!iEvent.getByToken(fSiPixelClusterToken_, hClusterColl)) {
118  edm::LogWarning("SiPixelStatusProducer")
119  << " edmNew::DetSetVector<SiPixelCluster> " << fPixelClusterLabel_ << " does not exist!" << std::endl;
120  return;
121  }
122 
123  iEvent.getByToken(fSiPixelClusterToken_, hClusterColl);
124 
125  if (hClusterColl.isValid()) {
126  for (const auto& clusters : *hClusterColl) { /*loop over different clusters in a clusters vector (module)*/
127  for (const auto& clu : clusters) { /*loop over cluster in a given detId (module)*/
128  int detid = clusters.detId();
129  int rowsperroc = fSensors_[detid].first;
130  int colsperroc = fSensors_[detid].second;
131 
132  //int nROCrows = fSensorLayout_[detid].first;
133  int nROCcolumns = fSensorLayout_[detid].second;
134 
135  int roc(-1);
136  std::map<int, int> rocIds_detid;
137  if (fRocIds_.find(detid) != fRocIds_.end()) {
138  rocIds_detid = fRocIds_[detid];
139  }
140 
141  /* A module is made with a few ROCs
142  Need to convert global row/column (on a module) to local row/column (on a ROC) */
143  const std::vector<SiPixelCluster::Pixel>& pixvector = clu.pixels();
144  for (unsigned int i = 0; i < pixvector.size(); ++i) {
145  int mr0 = pixvector[i].x; /* constant column direction is along x-axis */
146  int mc0 = pixvector[i].y; /* constant row direction is along y-axis */
147 
148  int irow = mr0 / rowsperroc;
149  int icol = mc0 / colsperroc;
150 
151  int key = indexROC(irow, icol, nROCcolumns);
152  if (rocIds_detid.find(key) != rocIds_detid.end()) {
153  roc = rocIds_detid[key];
154  }
155 
156  fDet_.fillDIGI(detid, roc);
157 
158  } /* loop over pixels in a cluster */
159 
160  } /* loop over cluster in a detId (module) */
161 
162  } /* loop over detId-grouped clusters in cluster detId-grouped clusters-vector* */
163 
164  } /* hClusterColl.isValid() */
165  else {
166  edm::LogWarning("SiPixelStatusProducer")
167  << " edmNew::DetSetVector<SiPixelCluster> " << fPixelClusterLabel_ << " is NOT Valid!" << std::endl;
168  }
169 
170  /*|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||*/
171 
172  /* FEDerror25 per-ROC per-event */
173  edm::Handle<PixelFEDChannelCollection> pixelFEDChannelCollectionHandle;
174 
175  /* look over different resouces of tokens */
177  if (!iEvent.getByToken(tk, pixelFEDChannelCollectionHandle)) {
178  edm::LogWarning("SiPixelStatusProducer")
179  << " PixelFEDChannelCollection with index " << tk.index() << " does NOT exist!" << std::endl;
180  continue;
181  }
182 
183  iEvent.getByToken(tk, pixelFEDChannelCollectionHandle);
184  if (!pixelFEDChannelCollectionHandle.isValid()) {
185  edm::LogWarning("SiPixelStatusProducer")
186  << " PixelFEDChannelCollection with index " << tk.index() << " is NOT valid!" << std::endl;
187  continue;
188  }
189 
190  /* FEDerror channels for the current events */
191  std::map<int, std::vector<PixelFEDChannel>> tmpFEDerror25;
192  for (const auto& disabledChannels : *pixelFEDChannelCollectionHandle) {
193  /* loop over different PixelFED in a PixelFED vector (module) */
194  for (const auto& ch : disabledChannels) {
195  DetId detId = disabledChannels.detId();
196  int detid = detId.rawId();
197 
198  if (ftotalevents_ == 1) {
199  /* FEDerror25 channels for the "first" event in the lumi section (first for the current instance of the module) */
200  fFEDerror25_[detid].push_back(ch);
201  } else
202  tmpFEDerror25[detid].push_back(ch);
203 
204  } /* loop over different PixelFED in a PixelFED vector (different channel for a module) */
205 
206  } /* loop over different (different DetId) PixelFED vectors in PixelFEDChannelCollection */
207 
208  /* Compare the current FEDerror list with the first event's FEDerror list
209  * and save the common channels */
210  if (!tmpFEDerror25.empty() && !fFEDerror25_.empty()) {
211  std::map<int, std::vector<PixelFEDChannel>>::iterator itFEDerror25;
212  for (itFEDerror25 = fFEDerror25_.begin(); itFEDerror25 != fFEDerror25_.end(); itFEDerror25++) {
213  int detid = itFEDerror25->first;
214  if (tmpFEDerror25.find(detid) != tmpFEDerror25.end()) {
215  std::vector<PixelFEDChannel> chs = itFEDerror25->second;
216  std::vector<PixelFEDChannel> chs_tmp = tmpFEDerror25[detid];
217 
218  std::vector<PixelFEDChannel> chs_common;
219  for (unsigned int ich = 0; ich < chs.size(); ich++) {
220  PixelFEDChannel ch = chs[ich];
221  /* loop over the current FEDerror25 channels, save the common FED channels */
222  for (unsigned int ich_tmp = 0; ich_tmp < chs_tmp.size(); ich_tmp++) {
223  PixelFEDChannel ch_tmp = chs_tmp[ich_tmp];
224  if ((ch.fed == ch_tmp.fed) && (ch.link == ch_tmp.link)) { /* the same FED channel */
225  chs_common.push_back(ch);
226  break;
227  }
228  }
229  }
230  /* remove the full module from FEDerror25 list if no common channels are left */
231  if (chs_common.empty())
232  fFEDerror25_.erase(itFEDerror25);
233  /* otherwise replace with the common channels */
234  else {
235  fFEDerror25_[detid].clear();
236  fFEDerror25_[detid] = chs_common;
237  }
238  } else { /* remove the full module from FEDerror25 list if the module doesn't appear in the current event's FEDerror25 list */
239  fFEDerror25_.erase(itFEDerror25);
240  }
241 
242  } /* loop over modules that have FEDerror25 in the first event in the lumi section */
243 
244  } /* non-empty FEDerror lists */
245 
246  } /* look over different resouces of takens */
247 
248  /* Caveat
249  no per-event collection put into iEvent
250  If use produce() function and no collection is put into iEvent, produce() will not run in unScheduled mode
251  Now since CMSSW_10_1_X, the accumulate() function will run whatsoever in the unScheduled mode
252  Accumulate() is NOT available for releases BEFORE CMSSW_10_1_X */
253 }
254 
256  /* set total number of events through ftotalevents_ */
258 
259  if (ftotalevents_ > 0) {
260  /* Add FEDerror25 information into SiPixelDetectorStatus fDet_ for FED channels stored in fFEDerror25_ */
261  if (!fFEDerror25_.empty()) { // non-empty FEDerror25
262  std::map<int, std::vector<PixelFEDChannel>>::iterator itFEDerror25;
263  for (itFEDerror25 = fFEDerror25_.begin(); itFEDerror25 != fFEDerror25_.end();
264  itFEDerror25++) { // loop over detIds
265  int detid = itFEDerror25->first;
266  std::vector<PixelFEDChannel> chs = itFEDerror25->second;
267  for (unsigned int ich = 0; ich < chs.size(); ich++) {
268  fDet_.fillFEDerror25(detid, chs[ich]);
269  }
270  } // loop over detIds
271  } // if non-empty FEDerror25
272 
273  } // only for non-zero events
274 }
275 
277  edm::LuminosityBlock const& iLumi,
278  edm::EventSetup const&,
279  std::vector<SiPixelDetectorStatus>* siPixelDetectorStatusVtr) const {
280  /*add the Stream's partial information to the full information*/
281 
282  /* only save for the lumi sections with NON-ZERO events */
283  if (ftotalevents_ > 0)
284  siPixelDetectorStatusVtr->push_back(fDet_);
285 }
286 
287 /* helper function */
288 int SiPixelStatusProducer::indexROC(int irow, int icol, int nROCcolumns) {
289  return int(icol + irow * nROCcolumns);
290 
291  /* generate the folling roc index that is going to map with ROC id as
292  8 9 10 11 12 13 14 15
293  0 1 2 3 4 5 6 7 */
294 }
295 
static std::shared_ptr< SiPixelTopoFinder > globalBeginRun(edm::Run const &iRun, edm::EventSetup const &iSetup, GlobalCache const *iCache)
void addModule(int detid, int nrocs)
std::map< int, std::pair< int, int > > fSensors_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned int fed
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::vector< int > fDetIds_
edm::ESGetToken< SiPixelFedCablingMap, SiPixelFedCablingMapRcd > siPixelFedCablingMapToken_
void fillFEDerror25(int detid, PixelFEDChannel ch)
unsigned long int ftotalevents_
unsigned int link
void fillDIGI(int detid, int roc)
std::unordered_map< uint32_t, unsigned int > fFedIds_
int iEvent
Definition: GenABIO.cc:224
void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) final
SiPixelDetectorStatus fDet_
void setNevents(unsigned long int N)
std::vector< edm::EDGetTokenT< PixelFEDChannelCollection > > theBadPixelFEDChannelsTokens_
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > fSiPixelClusterToken_
key
prepare the HTCondor submission files and eventually submit them
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopologyToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeometryToken_
Log< level::Info, false > LogInfo
std::map< int, std::vector< PixelFEDChannel > > fFEDerror25_
Definition: DetId.h:17
std::map< int, std::pair< int, int > > fSensorLayout_
void accumulate(edm::Event const &iEvent, edm::EventSetup const &iSetup) final
bool isValid() const
Definition: HandleBase.h:70
SiPixelStatusProducer(edm::ParameterSet const &iPSet, SiPixelStatusCache const *)
void beginRun(edm::Run const &, edm::EventSetup const &) final
void endLuminosityBlockSummary(edm::LuminosityBlock const &iLumi, edm::EventSetup const &, std::vector< SiPixelDetectorStatus > *siPixelDetectorStatusVtr) const final
std::map< int, std::map< int, int > > fRocIds_
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) final
virtual int indexROC(int irow, int icol, int nROCcolumns) final
Log< level::Warning, false > LogWarning
Definition: Run.h:45