CMS 3D CMS Logo

SiStripFineDelayHit.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiStripFineDelayHit
4 // Class: SiStripFineDelayHit
5 //
13 //
14 // Original Author: Christophe DELAERE
15 // Created: Fri Nov 17 10:52:42 CET 2006
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <utility>
22 #include <vector>
23 #include <algorithm>
24 
25 // user include files
58 
59 // ROOT includes
60 #include "TMath.h"
61 
62 //
63 // constructors and destructor
64 //
66  //register your products
67  produces<edm::DetSetVector<SiStripRawDigi> >("FineDelaySelection");
68  //now do what ever other initialization is needed
70  cosmic_ = iConfig.getParameter<bool>("cosmic");
71  field_ = iConfig.getParameter<bool>("MagneticField");
72  maxAngle_ = iConfig.getParameter<double>("MaxTrackAngle");
73  minTrackP2_ = iConfig.getParameter<double>("MinTrackMomentum") * iConfig.getParameter<double>("MinTrackMomentum");
74  maxClusterDistance_ = iConfig.getParameter<double>("MaxClusterDistance");
75  /*
76  clusterLabel_ = iConfig.getParameter<edm::InputTag>("ClustersLabel");
77  trackLabel_ = iConfig.getParameter<edm::InputTag>("TracksLabel");
78  seedLabel_ = iConfig.getParameter<edm::InputTag>("SeedsLabel");
79  inputModuleLabel_ = iConfig.getParameter<edm::InputTag>( "InputModuleLabel" ) ;
80  digiLabel_ = iConfig.getParameter<edm::InputTag>("DigiLabel");
81  */
83  consumes<edmNew::DetSetVector<SiStripCluster> >(iConfig.getParameter<edm::InputTag>("ClustersLabel"));
84  trackToken_ = consumes<std::vector<Trajectory> >(iConfig.getParameter<edm::InputTag>("TracksLabel"));
85  trackCollectionToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("TracksLabel"));
86  seedcollToken_ = consumes<TrajectorySeedCollection>(iConfig.getParameter<edm::InputTag>("SeedsLabel"));
87  inputModuleToken_ = consumes<SiStripEventSummary>(iConfig.getParameter<edm::InputTag>("InputModuleLabel"));
88  digiToken_ = consumes<edm::DetSetVector<SiStripDigi> >(iConfig.getParameter<edm::InputTag>("DigiLabel"));
89 
90  homeMadeClusters_ = iConfig.getParameter<bool>("NoClustering");
91  explorationWindow_ = iConfig.getParameter<uint32_t>("ExplorationWindow");
92  noTracking_ = iConfig.getParameter<bool>("NoTracking");
93  mode_ = 0;
94 
97  fedCablingToken_ = esConsumes<edm::Transition::BeginRun>();
99 }
100 
102  // do anything here that needs to be done at desctruction time
103  // (e.g. close files, deallocate resources etc.)
104  delete anglefinder_;
105 }
106 
107 //
108 // member functions
109 //
111  const int substructure,
112  const TrackerTopology* tkrTopo) {
113  uint32_t rootDetId = 0;
114  uint32_t maskDetId = 0;
115 
116  switch (subdet) {
117  case StripSubdetector::TIB: {
118  rootDetId = tkrTopo->tibDetId(substructure, 0, 0, 0, 0, 0).rawId();
119  maskDetId = tkrTopo->tibDetId(15, 0, 0, 0, 0, 0).rawId();
120  break;
121  }
122  case StripSubdetector::TID: {
123  rootDetId = tkrTopo->tidDetId(substructure > 0 ? 2 : 1, abs(substructure), 0, 0, 0, 0).rawId();
124  maskDetId = tkrTopo->tidDetId(3, 15, 0, 0, 0, 0).rawId();
125  break;
126  }
127  case StripSubdetector::TOB: {
128  rootDetId = tkrTopo->tobDetId(substructure, 0, 0, 0, 0).rawId();
129  maskDetId = tkrTopo->tobDetId(15, 0, 0, 0, 0).rawId();
130  break;
131  }
132  case StripSubdetector::TEC: {
133  rootDetId = tkrTopo->tecDetId(substructure > 0 ? 2 : 1, abs(substructure), 0, 0, 0, 0, 0).rawId();
134  maskDetId = tkrTopo->tecDetId(3, 15, 0, 0, 0, 0, 0).rawId();
135  break;
136  }
137  default:
138  break;
139  }
140  return std::make_pair(maskDetId, rootDetId);
141 }
142 
143 std::vector<std::pair<uint32_t, std::pair<double, double> > > SiStripFineDelayHit::detId(
144  const TrackerGeometry& tracker,
145  const TrackerTopology* tkrTopo,
146  const reco::Track* tk,
147  const std::vector<Trajectory>& trajVec,
148  const StripSubdetector::SubDetector subdet,
149  const int substructure) {
150  if (substructure == 0xff)
151  return detId(tracker, tkrTopo, tk, trajVec, 0, 0);
152  // first determine the root detId we are looking for
153  DeviceMask mask = deviceMask(subdet, substructure, tkrTopo);
154  // then call the method that loops on recHits
155  return detId(tracker, tkrTopo, tk, trajVec, mask.first, mask.second);
156 }
157 
158 std::vector<std::pair<uint32_t, std::pair<double, double> > > SiStripFineDelayHit::detId(
159  const TrackerGeometry& tracker,
160  const TrackerTopology* tkrTopo,
161  const reco::Track* tk,
162  const std::vector<Trajectory>& trajVec,
163  const uint32_t& maskDetId,
164  const uint32_t& rootDetId) {
165  bool onDisk = ((maskDetId == tkrTopo->tidDetId(3, 15, 0, 0, 0, 0).rawId()) ||
166  (maskDetId == tkrTopo->tecDetId(3, 15, 0, 0, 0, 0, 0).rawId()));
167  std::vector<std::pair<uint32_t, std::pair<double, double> > > result;
168  std::vector<uint32_t> usedDetids;
169  // now loop on recHits to find the right detId plus the track local angle
170  std::vector<std::pair<std::pair<DetId, LocalPoint>, float> > hitangle;
171  if (!cosmic_) {
172  // use trajectories in event.
173  // we have first to find the right trajectory for the considered track.
174  for (std::vector<Trajectory>::const_iterator traj = trajVec.begin(); traj < trajVec.end(); ++traj) {
175  if (((traj->lastMeasurement().recHit()->geographicalId().rawId() ==
176  (*(tk->recHitsEnd() - 1))->geographicalId().rawId()) &&
177  (traj->lastMeasurement().recHit()->localPosition().x() == (*(tk->recHitsEnd() - 1))->localPosition().x())) ||
178  ((traj->firstMeasurement().recHit()->geographicalId().rawId() ==
179  (*(tk->recHitsEnd() - 1))->geographicalId().rawId()) &&
180  (traj->firstMeasurement().recHit()->localPosition().x() == (*(tk->recHitsEnd() - 1))->localPosition().x()))) {
181  hitangle = anglefinder_->findtrackangle(*traj);
182  break;
183  }
184  }
185  } else {
187  // event_->getByLabel(seedLabel_,seedcoll);
188  event_->getByToken(seedcollToken_, seedcoll);
189  // use trajectories in event.
190  hitangle = anglefinder_->findtrackangle(trajVec);
191  }
192  LogDebug("DetId") << "number of hits for the track: " << hitangle.size();
193  std::vector<std::pair<std::pair<DetId, LocalPoint>, float> >::iterator iter;
194  // select the interesting DetIds, based on the ID and TLA
195  for (iter = hitangle.begin(); iter != hitangle.end(); iter++) {
196  // check the detId.
197  // if substructure was 0xff, then maskDetId and rootDetId == 0
198  // this implies all detids are accepted. (also if maskDetId=rootDetId=0 explicitely).
199  // That "unusual" mode of operation allows to analyze also Latency scans
200  LogDebug("DetId") << "check the detid: " << std::hex << (iter->first.first.rawId()) << " vs " << rootDetId
201  << " with a mask of " << maskDetId << std::dec << std::endl;
202 
203  if (((iter->first.first.rawId() & maskDetId) != rootDetId))
204  continue;
205  if (std::find(usedDetids.begin(), usedDetids.end(), iter->first.first.rawId()) != usedDetids.end())
206  continue;
207  // check the local angle (extended to the equivalent angle correction)
208  LogDebug("DetId") << "check the angle: " << fabs((iter->second));
209  if (1 - fabs(fabs(iter->second) - 1) < cos(maxAngle_ / 180. * TMath::Pi()))
210  continue;
211  // returns the detid + the time of flight to there
212  std::pair<uint32_t, std::pair<double, double> > el;
213  std::pair<double, double> subel;
214  el.first = iter->first.first.rawId();
215  // here, we compute the TOF.
216  // For cosmics, some track parameters are missing. Parameters are recomputed.
217  // for our calculation, the track momemtum at any point is enough:
218  // only used without B field or for the sign of Pz.
219  double trackParameters[5];
220  for (int i = 0; i < 5; i++)
221  trackParameters[i] = tk->parameters()[i];
222  if (cosmic_)
223  SiStripFineDelayTOF::trackParameters(*tk, trackParameters);
224  double hit[3];
225  const GeomDetUnit* det(tracker.idToDetUnit(iter->first.first));
226  Surface::GlobalPoint gp = det->surface().toGlobal(iter->first.second);
227  hit[0] = gp.x();
228  hit[1] = gp.y();
229  hit[2] = gp.z();
230  double phit[3];
231  phit[0] = tk->momentum().x();
232  phit[1] = tk->momentum().y();
233  phit[2] = tk->momentum().z();
234  subel.first = SiStripFineDelayTOF::timeOfFlight(cosmic_, field_, trackParameters, hit, phit, onDisk);
235  subel.second = iter->second;
236  el.second = subel;
237  // returns the detid + TOF
238  result.push_back(el);
239  usedDetids.push_back(el.first);
240  }
241  return result;
242 }
243 
244 bool SiStripFineDelayHit::rechit(reco::Track* tk, uint32_t det_id) {
245  for (trackingRecHit_iterator it = tk->recHitsBegin(); it != tk->recHitsEnd(); it++)
246  if ((*it)->geographicalId().rawId() == det_id) {
247  return (*it)->isValid();
248  break;
249  }
250  return false;
251 }
252 
253 // VI January 2012: FIXME
254 // do not understand what is going on here: each hit has a cluster: by definition will be the closest!
255 std::pair<const SiStripCluster*, double> SiStripFineDelayHit::closestCluster(
256  const TrackerGeometry& tracker,
257  const reco::Track* tk,
258  const uint32_t& det_id,
261  std::pair<const SiStripCluster*, double> result(nullptr, 0.);
262  double hitStrip = -1;
263  int nstrips = -1;
264  // localize the crossing point of the track on the module
265  for (trackingRecHit_iterator it = tk->recHitsBegin(); it != tk->recHitsEnd(); it++) {
266  LogDebug("closestCluster") << "(*it)->geographicalId().rawId() vs det_id" << (*it)->geographicalId().rawId() << " "
267  << det_id;
268  //handle the mono rechits
269  if ((*it)->geographicalId().rawId() == det_id) {
270  if (!(*it)->isValid())
271  continue;
272  LogDebug("closestCluster") << " using the single mono hit";
273  LocalPoint lp = (*it)->localPosition();
274  const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet((*it)->geographicalId()));
276  hitStrip = p.x();
277  nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips();
278  break;
279  }
280  /* FIXME: local position is not there anymore...
281  //handle stereo part of matched hits
282  //one could try to cast to SiStripMatchedRecHit2D but it is faster to look at the detid
283  else if((det_id - (*it)->geographicalId().rawId())==1) {
284  const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it));
285  if(!hit2D) continue; // this is a security that should never trigger
286  const SiStripRecHit2D* stereo = hit2D->stereoHit();
287  if(!stereo) continue; // this is a security that should never trigger
288  if(!stereo->isValid()) continue;
289  LogDebug("closestCluster") << " using the stereo hit";
290  LocalPoint lp = stereo->localPosition();
291  const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet(stereo->geographicalId()));
292  MeasurementPoint p = gdu->topology().measurementPosition(lp);
293  hitStrip = p.x();
294  nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips();
295  break;
296  }
297  //handle mono part of matched hits
298  //one could try to cast to SiStripMatchedRecHit2D but it is faster to look at the detid
299  else if((det_id - (*it)->geographicalId().rawId())==2) {
300  const SiStripMatchedRecHit2D* hit2D = dynamic_cast<const SiStripMatchedRecHit2D*>(&(**it));
301  if(!hit2D) continue; // this is a security that should never trigger
302  const SiStripRecHit2D* mono = hit2D->monoHit();
303  if(!mono) continue; // this is a security that should never trigger
304  if(!mono->isValid()) continue;
305  LogDebug("closestCluster") << " using the mono hit";
306  LocalPoint lp = mono->localPosition();
307  const GeomDetUnit* gdu = static_cast<const GeomDetUnit*>(tracker.idToDet(mono->geographicalId()));
308  MeasurementPoint p = gdu->topology().measurementPosition(lp);
309  hitStrip = p.x();
310  nstrips = (dynamic_cast<const StripTopology*>(&(gdu->topology())))->nstrips();
311  break;
312  }
313  */
314  }
315  LogDebug("closestCluster") << " hit strip = " << hitStrip;
316  if (hitStrip < 0)
317  return result;
318  if (homeMadeClusters_) {
319  // take the list of digis on the module
320  for (edm::DetSetVector<SiStripDigi>::const_iterator DSViter = hits.begin(); DSViter != hits.end(); DSViter++) {
321  if (DSViter->id == det_id) {
322  // loop from hitstrip-n to hitstrip+n (explorationWindow_) and select the highest strip
323  int minStrip = int(round(hitStrip)) - explorationWindow_;
324  minStrip = minStrip < 0 ? 0 : minStrip;
325  int maxStrip = int(round(hitStrip)) + explorationWindow_ + 1;
326  maxStrip = maxStrip >= nstrips ? nstrips - 1 : maxStrip;
327  edm::DetSet<SiStripDigi>::const_iterator rangeStart = DSViter->end();
328  edm::DetSet<SiStripDigi>::const_iterator rangeStop = DSViter->end();
329  for (edm::DetSet<SiStripDigi>::const_iterator digiIt = DSViter->begin(); digiIt != DSViter->end(); ++digiIt) {
330  if (digiIt->strip() >= minStrip && rangeStart == DSViter->end())
331  rangeStart = digiIt;
332  if (digiIt->strip() <= maxStrip)
333  rangeStop = digiIt;
334  }
335  if (rangeStart != DSViter->end()) {
336  if (rangeStop != DSViter->end())
337  ++rangeStop;
338  // build a fake cluster
339  LogDebug("closestCluster") << "build a fake cluster ";
340  SiStripCluster* newCluster =
341  new SiStripCluster(SiStripCluster::SiStripDigiRange(rangeStart, rangeStop)); // /!\ ownership transfered
342  result.first = newCluster;
343  result.second = fabs(newCluster->barycenter() - hitStrip);
344  }
345  break;
346  }
347  }
348  } else {
349  // loop on the detsetvector<cluster> to find the right one
350  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = clusters.begin(); DSViter != clusters.end();
351  DSViter++)
352  if (DSViter->id() == det_id) {
353  LogDebug("closestCluster") << " detset with the right detid. ";
354  edmNew::DetSet<SiStripCluster>::const_iterator begin = DSViter->begin();
356  //find the cluster close to the hitStrip
357  result.second = 1000.;
358  for (edmNew::DetSet<SiStripCluster>::const_iterator iter = begin; iter != end; ++iter) {
359  double dist = fabs(iter->barycenter() - hitStrip);
360  if (dist < result.second) {
361  result.second = dist;
362  result.first = &(*iter);
363  }
364  }
365  break;
366  }
367  }
368  return result;
369 }
370 
371 // ------------ method called to produce the data ------------
373  using namespace edm;
374  // Retrieve commissioning information from "event summary"
376  // iEvent.getByLabel( inputModuleLabel_, runsummary );
377  iEvent.getByToken(inputModuleToken_, runsummary);
378  if (runsummary->runType() == sistrip::APV_LATENCY)
379  mode_ = 2; // LatencyScan
380  else if (runsummary->runType() == sistrip::FINE_DELAY)
381  mode_ = 1; // DelayScan
382  else {
383  mode_ = 0; //unknown
384  return;
385  }
386 
387  if (noTracking_) {
388  produceNoTracking(iEvent, iSetup);
389  return;
390  }
391  event_ = &iEvent;
392  // container for the selected hits
393  std::vector<edm::DetSet<SiStripRawDigi> > output;
394  output.reserve(100);
395  // access the tracks
397  // iEvent.getByLabel(trackLabel_,trackCollection);
399  const reco::TrackCollection* tracks = trackCollection.product();
400  const auto& tracker = iSetup.getData(tkGeomToken_);
401  if (!tracks->empty()) {
402  anglefinder_->init(iEvent, iSetup);
403  LogDebug("produce") << "Found " << tracks->size() << " tracks.";
404  // look at the hits if one needs them
406  const edm::DetSetVector<SiStripDigi>* hitSet = nullptr;
407  if (homeMadeClusters_) {
408  // iEvent.getByLabel(digiLabel_,hits);
409  iEvent.getByToken(digiToken_, hits);
410  hitSet = hits.product();
411  }
412  // look at the clusters
414  // iEvent.getByLabel(clusterLabel_, clusters);
415  iEvent.getByToken(clustersToken_, clusters);
416  const edmNew::DetSetVector<SiStripCluster>* clusterSet = clusters.product();
417  // look at the trajectories if they are in the event
418  std::vector<Trajectory> trajVec;
420  // iEvent.getByLabel(trackLabel_,TrajectoryCollection);
422  trajVec = *(TrajectoryCollection.product());
423  // Get TrackerTopology
424  const auto tTopo = &iSetup.getData(tTopoToken_);
425  // loop on tracks
426  for (reco::TrackCollection::const_iterator itrack = tracks->begin(); itrack < tracks->end(); itrack++) {
427  // first check the track Pt
428  if ((itrack->px() * itrack->px() + itrack->py() * itrack->py() + itrack->pz() * itrack->pz()) < minTrackP2_)
429  continue;
430  // check that we have something in the layer we are interested in
431  std::vector<std::pair<uint32_t, std::pair<double, double> > > intersections;
432  if (mode_ == 1) {
433  // Retrieve and decode commissioning information from "event summary"
435  // iEvent.getByLabel( inputModuleLabel_, summary );
436  iEvent.getByToken(inputModuleToken_, summary);
437  uint32_t layerCode = (const_cast<SiStripEventSummary*>(summary.product())->layerScanned()) >> 16;
439  if (((layerCode >> 6) & 0x3) == 0)
440  subdet = StripSubdetector::TIB;
441  else if (((layerCode >> 6) & 0x3) == 1)
442  subdet = StripSubdetector::TOB;
443  else if (((layerCode >> 6) & 0x3) == 2)
444  subdet = StripSubdetector::TID;
445  else if (((layerCode >> 6) & 0x3) == 3)
446  subdet = StripSubdetector::TEC;
447  int32_t layerIdx = (layerCode & 0xF) * (((layerCode >> 4) & 0x3) ? -1 : 1);
448  intersections = detId(tracker, tTopo, &(*itrack), trajVec, subdet, layerIdx);
449  } else {
450  // for latency scans, no layer is specified -> no cut on detid
451  intersections = detId(tracker, tTopo, &(*itrack), trajVec);
452  }
453  LogDebug("produce") << " Found " << intersections.size() << " interesting intersections." << std::endl;
454  for (std::vector<std::pair<uint32_t, std::pair<double, double> > >::iterator it = intersections.begin();
455  it < intersections.end();
456  it++) {
457  std::pair<const SiStripCluster*, double> candidateCluster =
458  closestCluster(tracker, &(*itrack), it->first, *clusterSet, *hitSet);
459  if (candidateCluster.first) {
460  LogDebug("produce") << " Found a cluster." << std::endl;
461  // cut on the distance
462  if (candidateCluster.second > maxClusterDistance_)
463  continue;
464  LogDebug("produce") << " The cluster is close enough." << std::endl;
465  // build the rawdigi corresponding to the leading strip and save it
466  // here, only the leading strip is retained. All other rawdigis in the module are set to 0.
467  const auto& amplitudes = candidateCluster.first->amplitudes();
468  uint8_t leadingCharge = 0;
469  uint8_t leadingStrip = candidateCluster.first->firstStrip();
470  uint8_t leadingPosition = 0;
471  for (auto amplit = amplitudes.begin(); amplit < amplitudes.end(); amplit++, leadingStrip++) {
472  if (leadingCharge < *amplit) {
473  leadingCharge = *amplit;
474  leadingPosition = leadingStrip;
475  }
476  }
477 
478  // look for an existing detset
479  std::vector<edm::DetSet<SiStripRawDigi> >::iterator newdsit = output.begin();
480  for (; newdsit != output.end() && newdsit->detId() != connectionMap_[it->first]; ++newdsit) {
481  }
482  // if there is no detset yet, create it.
483  if (newdsit == output.end()) {
485  output.push_back(newds);
486  newdsit = output.end() - 1;
487  }
488 
489  LogDebug("produce") << " New Hit... TOF:" << it->second.first << ", charge: " << int(leadingCharge)
490  << " at " << int(leadingPosition) << "." << std::endl
491  << "Angular correction: " << it->second.second << " giving a final value of "
492  << int(leadingCharge * fabs(it->second.second))
493  << " for fed key = " << connectionMap_[it->first] << " (detid=" << it->first << ")";
494  // apply corrections to the leading charge, but only if it has not saturated.
495  if (leadingCharge < 255) {
496  // correct the leading charge for the crossing angle
497  leadingCharge = uint8_t(leadingCharge * fabs(it->second.second));
498  // correct for module thickness for TEC and TOB
499  if ((((it->first >> 25) & 0x7f) == 0xd) ||
500  ((((it->first >> 25) & 0x7f) == 0xe) && (((it->first >> 5) & 0x7) > 4)))
501  leadingCharge = uint8_t((leadingCharge * 0.64));
502  }
503  //code the time of flight in the digi
504  unsigned int tof = abs(int(round(it->second.first * 10)));
505  tof = tof > 255 ? 255 : tof;
506  SiStripRawDigi newSiStrip(leadingCharge + (tof << 8));
507  newdsit->push_back(newSiStrip);
508  LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added.";
509  }
510  if (homeMadeClusters_)
511  delete candidateCluster.first; // we are owner of home-made clusters
512  }
513  }
514  }
515  // add the selected hits to the event.
516  LogDebug("produce") << "Putting " << output.size() << " new hits in the event.";
517  std::unique_ptr<edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output));
518  iEvent.put(std::move(formatedOutput), "FineDelaySelection");
519 }
520 
521 // Simple solution when tracking is not available/ not working
523  event_ = &iEvent;
524  // Get TrackerTopology
525  const auto tTopo = &iSetup.getData(tTopoToken_);
526  // container for the selected hits
527  std::vector<edm::DetSet<SiStripRawDigi> > output;
528  output.reserve(100);
529  // Retrieve and decode commissioning information from "event summary"
531  // iEvent.getByLabel( inputModuleLabel_, summary );
532  iEvent.getByToken(inputModuleToken_, summary);
533  uint32_t layerCode = (const_cast<SiStripEventSummary*>(summary.product())->layerScanned()) >> 16;
535  if (((layerCode >> 6) & 0x3) == 0)
536  subdet = StripSubdetector::TIB;
537  else if (((layerCode >> 6) & 0x3) == 1)
538  subdet = StripSubdetector::TOB;
539  else if (((layerCode >> 6) & 0x3) == 2)
540  subdet = StripSubdetector::TID;
541  else if (((layerCode >> 6) & 0x3) == 3)
542  subdet = StripSubdetector::TEC;
543  int32_t layerIdx = (layerCode & 0xF) * (((layerCode >> 4) & 0x3) ? -1 : 1);
544  DeviceMask mask = deviceMask(subdet, layerIdx, tTopo);
545  // look at the clusters
547  // iEvent.getByLabel(clusterLabel_,clusters);
548  iEvent.getByToken(clustersToken_, clusters);
549  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = clusters->begin(); DSViter != clusters->end();
550  DSViter++) {
551  // check that we are in the layer of interest
552  if (mode_ == 1 && ((DSViter->id() & mask.first) != mask.second))
553  continue;
554  // iterate over clusters
555  edmNew::DetSet<SiStripCluster>::const_iterator begin = DSViter->begin();
557  edm::DetSet<SiStripRawDigi> newds(connectionMap_[DSViter->id()]);
558  for (edmNew::DetSet<SiStripCluster>::const_iterator iter = begin; iter != end; ++iter) {
559  // build the rawdigi corresponding to the leading strip and save it
560  // here, only the leading strip is retained. All other rawdigis in the module are set to 0.
561  auto const& amplitudes = iter->amplitudes();
562  uint8_t leadingCharge = 0;
563  uint8_t leadingStrip = iter->firstStrip();
564  uint8_t leadingPosition = 0;
565  for (auto amplit = amplitudes.begin(); amplit < amplitudes.end(); amplit++, leadingStrip++) {
566  if (leadingCharge < *amplit) {
567  leadingCharge = *amplit;
568  leadingPosition = leadingStrip;
569  }
570  }
571  // apply some sanity cuts. This is needed since we don't use tracking to clean clusters
572  // 1.5< noise <8
573  // charge<250
574  // 50 > s/n > 10
575  const auto& noises = iSetup.getData(noiseToken_);
576  SiStripNoises::Range detNoiseRange = noises.getRange(DSViter->id());
577  float noise = noises.getNoise(leadingPosition, detNoiseRange);
578  if (noise < 1.5)
579  continue;
580  if (leadingCharge >= 250 || noise >= 8 || leadingCharge / noise > 50 || leadingCharge / noise < 10)
581  continue;
582  // apply some correction to the leading charge, but only if it has not saturated.
583  if (leadingCharge < 255) {
584  // correct for modulethickness for TEC and TOB
585  if ((((((DSViter->id()) >> 25) & 0x7f) == 0xd) || ((((DSViter->id()) >> 25) & 0x7f) == 0xe)) &&
586  ((((DSViter->id()) >> 5) & 0x7) > 4))
587  leadingCharge = uint8_t((leadingCharge * 0.64));
588  }
589  //code the time of flight == 0 in the digi
590  SiStripRawDigi newSiStrip(leadingCharge);
591  newds.push_back(newSiStrip);
592  }
593  //store into the detsetvector
594  output.push_back(newds);
595  LogDebug("produce") << "New edm::DetSet<SiStripRawDigi> added with fedkey = " << std::hex << std::setfill('0')
596  << std::setw(8) << connectionMap_[DSViter->id()] << std::dec;
597  }
598  // add the selected hits to the event.
599  LogDebug("produce") << "Putting " << output.size() << " new hits in the event.";
600  std::unique_ptr<edm::DetSetVector<SiStripRawDigi> > formatedOutput(new edm::DetSetVector<SiStripRawDigi>(output));
601  iEvent.put(std::move(formatedOutput), "FineDelaySelection");
602 }
603 
604 // ------------ method called once each job just before starting event loop ------------
606  // Retrieve FED cabling object
607  const auto& cabling = iSetup.getData(fedCablingToken_);
608  auto feds = cabling.fedIds();
609  for (auto fedid = feds.begin(); fedid < feds.end(); ++fedid) {
610  auto connections = cabling.fedConnections(*fedid);
611  for (std::vector<FedChannelConnection>::const_iterator conn = connections.begin(); conn < connections.end();
612  ++conn) {
613  /*
614  SiStripFedKey key(conn->fedId(),
615  SiStripFedKey::feUnit(conn->fedCh()),
616  SiStripFedKey::feChan(conn->fedCh()));
617  connectionMap_[conn->detId()] = key.key();
618  */
619  // the key is computed using an alternate formula for performance reasons.
620  connectionMap_[conn->detId()] = ((conn->fedId() & sistrip::invalid_) << 16) | (conn->fedCh() & sistrip::invalid_);
621  }
622  }
623 }
iterator end()
Definition: DetSet.h:58
const double Pi
static constexpr auto TEC
virtual void produceNoTracking(edm::Event &, const edm::EventSetup &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
DetId tibDetId(uint32_t layer, uint32_t str_fw_bw, uint32_t str_int_ext, uint32_t str, uint32_t module, uint32_t ster) const
virtual const Topology & topology() const
Definition: GeomDet.cc:67
std::pair< SiStripDigiIter, SiStripDigiIter > SiStripDigiRange
edm::EDGetTokenT< SiStripEventSummary > inputModuleToken_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
void init(const edm::Event &e, const edm::EventSetup &c)
const sistrip::RunType & runType() const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
DetId tobDetId(uint32_t layer, uint32_t rod_fw_bw, uint32_t rod, uint32_t module, uint32_t ster) const
Constants and enumerated type for sistrip::RunType.
data_type const * const_iterator
Definition: DetSetNew.h:31
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
constexpr uint32_t mask
Definition: gpuClustering.h:26
std::vector< std::pair< uint32_t, std::pair< double, double > > > detId(const TrackerGeometry &tracker, const TrackerTopology *tkrTopo, const reco::Track *tk, const std::vector< Trajectory > &trajVec, const StripSubdetector::SubDetector subdet=StripSubdetector::TIB, const int substructure=0xff)
DetId tidDetId(uint32_t side, uint32_t wheel, uint32_t ring, uint32_t module_fw_bw, uint32_t module, uint32_t ster) const
std::vector< std::pair< std::pair< DetId, LocalPoint >, float > > findtrackangle(const std::vector< Trajectory > &traj)
SiStripFineDelayHit(const edm::ParameterSet &)
std::pair< const SiStripCluster *, double > closestCluster(const TrackerGeometry &tracker, const reco::Track *tk, const uint32_t &detId, const edmNew::DetSetVector< SiStripCluster > &clusters, const edm::DetSetVector< SiStripDigi > &hits)
int iEvent
Definition: GenABIO.cc:224
edm::ESGetToken< SiStripNoises, SiStripNoisesRcd > noiseToken_
static double timeOfFlight(bool cosmics, bool field, double *trackParameters, double *hit, double *phit, bool onDisk)
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
SiStripFineDelayTLA * anglefinder_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > clustersToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:91
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
static void trackParameters(const reco::Track &tk, double *trackParameters)
static constexpr auto TOB
DeviceMask deviceMask(const StripSubdetector::SubDetector subdet, const int substructure, const TrackerTopology *tkrTopo)
trackCollection
Definition: JetHT_cfg.py:51
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:88
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
edm::EDGetTokenT< TrajectorySeedCollection > seedcollToken_
static constexpr auto TIB
edm::EDGetTokenT< edm::DetSetVector< SiStripDigi > > digiToken_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::pair< uint32_t, uint32_t > DeviceMask
static const uint16_t invalid_
Definition: Constants.h:16
bool rechit(reco::Track *tk, uint32_t detId)
std::map< uint32_t, uint32_t > connectionMap_
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
void beginRun(const edm::Run &, const edm::EventSetup &) override
HLT enums.
const edm::Event * event_
conn
Definition: getInfo.py:9
float barycenter() const
void produce(edm::Event &, const edm::EventSetup &) override
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
std::vector< Trajectory > TrajectoryCollection
edm::EDGetTokenT< std::vector< Trajectory > > trackToken_
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:711
Definition: output.py:1
edm::EDGetTokenT< reco::TrackCollection > trackCollectionToken_
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:47
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
edm::ESGetToken< SiStripFedCabling, SiStripFedCablingRcd > fedCablingToken_
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
static constexpr auto TID
DetId tecDetId(uint32_t side, uint32_t wheel, uint32_t petal_fw_bw, uint32_t petal, uint32_t ring, uint32_t module, uint32_t ster) const
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
#define LogDebug(id)