CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes
SiStripRecHitConverterAlgorithm Class Reference

#include <SiStripRecHitConverterAlgorithm.h>

Classes

struct  products
 

Public Member Functions

void initialize (const edm::EventSetup &)
 
void run (edm::Handle< edmNew::DetSetVector< SiStripCluster > > input, products &output)
 
void run (edm::Handle< edmNew::DetSetVector< SiStripCluster > > input, products &output, LocalVector trackdirection)
 
 SiStripRecHitConverterAlgorithm (const edm::ParameterSet &, edm::ConsumesCollector)
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &desc)
 

Private Types

typedef SiStripRecHit2DCollection::FastFiller Collector
 

Private Member Functions

void fillBad128StripBlocks (const uint32_t detid, bool bad128StripBlocks[6]) const
 
bool isMasked (const SiStripCluster &cluster, bool bad128StripBlocks[6]) const
 
void match (products &output, LocalVector trackdirection) const
 
bool useModule (const uint32_t id) const
 

Private Attributes

edm::ESGetToken< StripClusterParameterEstimator, TkStripCPERecordcpeToken
 
bool doMatching
 
bool maskBad128StripBlocks
 
const SiStripRecHitMatchermatcher = nullptr
 
edm::ESGetToken< SiStripRecHitMatcher, TkStripCPERecordmatcherToken
 
const StripClusterParameterEstimatorparameterestimator = nullptr
 
const SiStripQualityquality = nullptr
 
edm::ESGetToken< SiStripQuality, SiStripQualityRcdqualityToken
 
const TrackerGeometrytracker = nullptr
 
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtrackerToken
 
bool useQuality
 

Detailed Description

Definition at line 25 of file SiStripRecHitConverterAlgorithm.h.

Member Typedef Documentation

◆ Collector

Definition at line 69 of file SiStripRecHitConverterAlgorithm.h.

Constructor & Destructor Documentation

◆ SiStripRecHitConverterAlgorithm()

SiStripRecHitConverterAlgorithm::SiStripRecHitConverterAlgorithm ( const edm::ParameterSet conf,
edm::ConsumesCollector  iC 
)

Definition at line 16 of file SiStripRecHitConverterAlgorithm.cc.

References doMatching, edm::ConsumesCollector::esConsumes(), edm::ParameterSet::getParameter(), matcherToken, qualityToken, AlCaHarvesting_cff::SiStripQuality, and useQuality.

18  : useQuality(conf.getParameter<bool>("useSiStripQuality")),
19  maskBad128StripBlocks(conf.getParameter<bool>("MaskBadAPVFibers")),
20  doMatching(conf.getParameter<bool>("doMatching")),
23  conf.getParameter<edm::ESInputTag>("StripCPE"))) {
24  if (doMatching) {
26  }
27  if (useQuality) {
28  qualityToken =
29  iC.esConsumes<SiStripQuality, SiStripQualityRcd>(conf.getParameter<edm::ESInputTag>("siStripQualityLabel"));
30  }
31 }
edm::ESGetToken< SiStripRecHitMatcher, TkStripCPERecord > matcherToken
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< StripClusterParameterEstimator, TkStripCPERecord > cpeToken
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > qualityToken
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerToken

Member Function Documentation

◆ fillBad128StripBlocks()

void SiStripRecHitConverterAlgorithm::fillBad128StripBlocks ( const uint32_t  detid,
bool  bad128StripBlocks[6] 
) const
private

Definition at line 279 of file SiStripRecHitConverterAlgorithm.cc.

References SiStripQuality::getBadApvs(), SiStripQuality::getBadFibers(), dqmiolumiharvest::j, maskBad128StripBlocks, and quality.

Referenced by run().

279  {
280  if (maskBad128StripBlocks) {
281  short badApvs = quality->getBadApvs(detid);
282  short badFibers = quality->getBadFibers(detid);
283  for (int j = 0; j < 6; j++) {
284  bad128StripBlocks[j] = (badApvs & (1 << j));
285  }
286  for (int j = 0; j < 3; j++) {
287  if (badFibers & (1 << j)) {
288  bad128StripBlocks[2 * j + 0] = true;
289  bad128StripBlocks[2 * j + 1] = true;
290  }
291  }
292  }
293 }
short getBadFibers(uint32_t detid) const
short getBadApvs(uint32_t detid) const

◆ fillPSetDescription()

void SiStripRecHitConverterAlgorithm::fillPSetDescription ( edm::ParameterSetDescription desc)
static

Definition at line 33 of file SiStripRecHitConverterAlgorithm.cc.

References submitPVResolutionJobs::desc.

Referenced by SiStripRecHitConverter::fillDescriptions().

33  {
34  desc.add<bool>("useSiStripQuality", false);
35  desc.add<bool>("MaskBadAPVFibers", false);
36  desc.add<bool>("doMatching", true);
37  desc.add<edm::ESInputTag>("StripCPE", edm::ESInputTag("StripCPEfromTrackAngleESProducer", "StripCPEfromTrackAngle"));
38  desc.add<edm::ESInputTag>("Matcher", edm::ESInputTag("SiStripRecHitMatcherESProducer", "StandardMatcher"));
39  desc.add<edm::ESInputTag>("siStripQualityLabel", edm::ESInputTag());
40 }

◆ initialize()

void SiStripRecHitConverterAlgorithm::initialize ( const edm::EventSetup es)

Definition at line 42 of file SiStripRecHitConverterAlgorithm.cc.

References cpeToken, doMatching, edm::EventSetup::getData(), matcher, matcherToken, parameterestimator, quality, qualityToken, tracker, trackerToken, and useQuality.

Referenced by SiStripRecHitConverter::produce().

42  {
45  if (doMatching) {
47  }
48  if (useQuality) {
50  }
51 }
edm::ESGetToken< SiStripRecHitMatcher, TkStripCPERecord > matcherToken
edm::ESGetToken< StripClusterParameterEstimator, TkStripCPERecord > cpeToken
const StripClusterParameterEstimator * parameterestimator
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > qualityToken
bool getData(T &iHolder) const
Definition: EventSetup.h:122
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerToken

◆ isMasked()

bool SiStripRecHitConverterAlgorithm::isMasked ( const SiStripCluster cluster,
bool  bad128StripBlocks[6] 
) const
inlineprivate

Definition at line 295 of file SiStripRecHitConverterAlgorithm.cc.

References SiStripCluster::amplitudes(), SiStripCluster::barycenter(), SiStripCluster::firstStrip(), maskBad128StripBlocks, and SiStripCluster::size().

Referenced by run().

295  {
296  if (maskBad128StripBlocks) {
297  if (bad128StripBlocks[cluster.firstStrip() >> 7]) {
298  if (bad128StripBlocks[(cluster.firstStrip() + cluster.amplitudes().size()) >> 7] ||
299  bad128StripBlocks[static_cast<int32_t>(cluster.barycenter() - 0.499999) >> 7]) {
300  return true;
301  }
302  } else {
303  if (bad128StripBlocks[(cluster.firstStrip() + cluster.amplitudes().size()) >> 7] &&
304  bad128StripBlocks[static_cast<int32_t>(cluster.barycenter() - 0.499999) >> 7]) {
305  return true;
306  }
307  }
308  }
309  return false;
310 }
uint16_t firstStrip() const
SiStripCluster const & amplitudes() const
auto size() const
float barycenter() const

◆ match()

void SiStripRecHitConverterAlgorithm::match ( products output,
LocalVector  trackdirection 
) const
private

Definition at line 147 of file SiStripRecHitConverterAlgorithm.cc.

References edmNew::DetSet< T >::begin(), edm::OwnVector< T, P >::begin(), edm::OwnVector< T, P >::clear(), filterCSVwithJSON::copy, edmNew::DetSet< T >::detId(), SiStripRecHitMatcher::doubleMatch(), edmNew::DetSet< T >::empty(), edmNew::DetSet< T >::end(), edm::OwnVector< T, P >::end(), trigObjTnPSource_cfi::filler, TrackerGeometry::idToDet(), SiStripRecHitMatcher::match(), matcher, convertSQLitetoXML_cfg::output, StripSubdetector::partnerDetId(), edmNew::DetSet< T >::size(), edm::OwnVector< T, P >::size(), jetUpdater_cfi::sort, and tracker.

Referenced by run().

147  {
148  int nmatch = 0;
149  edm::OwnVector<SiStripMatchedRecHit2D> collectorMatched; // gp/FIXME: avoid this
150 
151  // Remember the ends of the collections, as we will use them a lot
152  SiStripRecHit2DCollection::const_iterator edStereoDet = output.stereo->end();
153  SiStripRecHit2DCollection::const_iterator edRPhiDet = output.rphi->end();
154 
155  // two work vectors for bookeeping clusters used by the stereo part of the matched hits
156  std::vector<SiStripRecHit2D::ClusterRef::key_type> matchedSteroClusters;
157 
158  for (SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->begin(); itRPhiDet != edRPhiDet;
159  ++itRPhiDet) {
160  edmNew::DetSet<SiStripRecHit2D> rphiHits = *itRPhiDet;
161  StripSubdetector specDetId(rphiHits.detId());
162  uint32_t partnerId = specDetId.partnerDetId();
163 
164  // if not part of a glued pair
165  if (partnerId == 0) {
166  // I must copy these as unmatched
167  if (!rphiHits.empty()) {
168  SiStripRecHit2DCollection::FastFiller filler(*output.rphiUnmatched, rphiHits.detId());
169  filler.resize(rphiHits.size());
170  std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
171  }
172  continue;
173  }
174 
175  SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->find(partnerId);
176 
177  // if the partner is not found (which probably can happen if it's empty)
178  if (itStereoDet == edStereoDet) {
179  // I must copy these as unmatched
180  if (!rphiHits.empty()) {
181  SiStripRecHit2DCollection::FastFiller filler(*output.rphiUnmatched, rphiHits.detId());
182  filler.resize(rphiHits.size());
183  std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
184  }
185  continue;
186  }
187 
188  edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
189 
190  // Get ready for making glued hits
191  const GluedGeomDet* gluedDet = (const GluedGeomDet*)tracker->idToDet(DetId(specDetId.glued()));
193  Collector collector(*output.matched, specDetId.glued());
194 
195  // Prepare also the list for unmatched rphi hits
196  SiStripRecHit2DCollection::FastFiller fillerRphiUnm(*output.rphiUnmatched, rphiHits.detId());
197 
198  // a list of clusters used by the matched part of the stereo hits in this detector
199  matchedSteroClusters.clear(); // at the beginning, empty
200 
201 #ifdef DOUBLE_MATCH
202  CollectorHelper chelper(collector, collectorMatched, fillerRphiUnm, matchedSteroClusters);
204  rphiHits.begin(), rphiHits.end(), stereoHits.begin(), stereoHits.end(), gluedDet, trackdirection, chelper);
205  nmatch += chelper.nmatch;
206 #else
207  // Make simple collection of this (gp:FIXME: why do we need it?)
209  // gp:FIXME: use std::transform
210  stereoSimpleHits.reserve(stereoHits.size());
211  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed;
212  ++it) {
213  stereoSimpleHits.push_back(&*it);
214  }
215 
216  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = rphiHits.begin(), ed = rphiHits.end(); it != ed; ++it) {
217  matcher->match(
218  &(*it), stereoSimpleHits.begin(), stereoSimpleHits.end(), collectorMatched, gluedDet, trackdirection);
219  if (collectorMatched.size() > 0) {
220  nmatch += collectorMatched.size();
222  edm = collectorMatched.end();
223  itm != edm;
224  ++itm) {
225  collector.push_back(*itm);
226  // mark the stereo hit cluster as used, so that the hit won't go in the unmatched stereo ones
227  matchedSteroClusters.push_back(itm->stereoClusterRef().key());
228  }
229  collectorMatched.clear();
230  } else {
231  // store a copy of this rphi hit as an unmatched rphi hit
232  fillerRphiUnm.push_back(*it);
233  }
234  }
235 
236 #endif
237 
238  // discard matched hits if the collection is empty
239  if (collector.empty())
240  collector.abort();
241 
242  // discard unmatched rphi hits if there are none
243  if (fillerRphiUnm.empty())
244  fillerRphiUnm.abort();
245 
246  // now look for unmatched stereo hits
247  SiStripRecHit2DCollection::FastFiller fillerStereoUnm(*output.stereoUnmatched, stereoHits.detId());
248  std::sort(matchedSteroClusters.begin(), matchedSteroClusters.end());
249  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed;
250  ++it) {
251  if (!std::binary_search(matchedSteroClusters.begin(), matchedSteroClusters.end(), it->cluster().key())) {
252  fillerStereoUnm.push_back(*it);
253  }
254  }
255  if (fillerStereoUnm.empty())
256  fillerStereoUnm.abort();
257  }
258 
259  for (SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->begin(); itStereoDet != edStereoDet;
260  ++itStereoDet) {
261  edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
262  StripSubdetector specDetId(stereoHits.detId());
263  uint32_t partnerId = specDetId.partnerDetId();
264  if (partnerId == 0)
265  continue;
266  SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->find(partnerId);
267  if (itRPhiDet == edRPhiDet) {
268  if (!stereoHits.empty()) {
269  SiStripRecHit2DCollection::FastFiller filler(*output.stereoUnmatched, stereoHits.detId());
270  filler.resize(stereoHits.size());
271  std::copy(stereoHits.begin(), stereoHits.end(), filler.begin());
272  }
273  }
274  }
275 
276  edm::LogInfo("SiStripRecHitConverter") << "found\n" << nmatch << " matched RecHits\n";
277 }
unsigned int partnerDetId() const
std::unique_ptr< SiStripMatchedRecHit2D > match(const SiStripRecHit2D *monoRH, const SiStripRecHit2D *stereoRH, const GluedGeomDet *gluedDet, LocalVector trackdirection, bool force) const
data_type const * const_iterator
Definition: DetSetNew.h:31
bool empty() const
Definition: DetSetNew.h:70
id_type detId() const
Definition: DetSetNew.h:66
iterator begin()
Definition: OwnVector.h:280
size_type size() const
Definition: OwnVector.h:300
void clear()
Definition: OwnVector.h:481
const TrackerGeomDet * idToDet(DetId) const override
Log< level::Info, false > LogInfo
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
iterator end()
Definition: OwnVector.h:285
std::vector< const SiStripRecHit2D * > SimpleHitCollection
SiStripRecHit2DCollection::FastFiller Collector
size_type size() const
Definition: DetSetNew.h:68
void doubleMatch(MonoIterator monoRHiter, MonoIterator monoRHend, StereoIterator seconditer, StereoIterator seconditerend, const GluedGeomDet *gluedDet, LocalVector trdir, CollectorHelper &collectorHelper) const
HLT enums.
iterator end()
Definition: DetSetNew.h:56
iterator begin()
Definition: DetSetNew.h:54

◆ run() [1/2]

void SiStripRecHitConverterAlgorithm::run ( edm::Handle< edmNew::DetSetVector< SiStripCluster > >  input,
products output 
)

Definition at line 53 of file SiStripRecHitConverterAlgorithm.cc.

References input, and convertSQLitetoXML_cfg::output.

Referenced by SiStripRecHitConverter::produce().

53  {
54  run(input, output, LocalVector(0., 0., 0.));
55 }
Local3DVector LocalVector
Definition: LocalVector.h:12
void run(edm::Handle< edmNew::DetSetVector< SiStripCluster > > input, products &output)
static std::string const input
Definition: EdmProvDump.cc:47

◆ run() [2/2]

void SiStripRecHitConverterAlgorithm::run ( edm::Handle< edmNew::DetSetVector< SiStripCluster > >  input,
products output,
LocalVector  trackdirection 
)

Definition at line 57 of file SiStripRecHitConverterAlgorithm.cc.

References edmNew::DetSetVector< T >::FastFiller::abort(), doMatching, edmNew::DetSetVector< T >::FastFiller::empty(), fillBad128StripBlocks(), TrackerGeometry::idToDetUnit(), isMasked(), StripClusterParameterEstimator::localParameters(), match(), convertSQLitetoXML_cfg::output, parameterestimator, edmNew::DetSetVector< T >::FastFiller::push_back(), StripSubdetector::stereo(), tracker, and useModule().

59  {
60  for (auto const& DS : *inputhandle) {
61  auto id = DS.id();
62  if (!useModule(id))
63  continue;
64 
65  Collector collector = StripSubdetector(id).stereo() ? Collector(*output.stereo, id) : Collector(*output.rphi, id);
66 
67  bool bad128StripBlocks[6];
68  fillBad128StripBlocks(id, bad128StripBlocks);
69 
70  GeomDetUnit const& du = *(tracker->idToDetUnit(id));
71  for (auto const& cluster : DS) {
72  if (isMasked(cluster, bad128StripBlocks))
73  continue;
74 
76  collector.push_back(
77  SiStripRecHit2D(parameters.first, parameters.second, du, DS.makeRefTo(inputhandle, &cluster)));
78  }
79 
80  if (collector.empty())
81  collector.abort();
82  }
83  if (doMatching) {
84  match(output, trackdirection);
85  }
86 }
std::pair< LocalPoint, LocalError > LocalValues
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
unsigned int stereo() const
stereo
virtual void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const
bool isMasked(const SiStripCluster &cluster, bool bad128StripBlocks[6]) const
const StripClusterParameterEstimator * parameterestimator
void fillBad128StripBlocks(const uint32_t detid, bool bad128StripBlocks[6]) const
void match(products &output, LocalVector trackdirection) const
SiStripRecHit2DCollection::FastFiller Collector

◆ useModule()

bool SiStripRecHitConverterAlgorithm::useModule ( const uint32_t  id) const
inlineprivate

Definition at line 312 of file SiStripRecHitConverterAlgorithm.cc.

References TrackerGeometry::idToDetUnit(), SiStripQuality::IsModuleUsable(), quality, tracker, and useQuality.

Referenced by run().

312  {
313  const StripGeomDetUnit* stripdet = (const StripGeomDetUnit*)tracker->idToDetUnit(id);
314  if (stripdet == nullptr)
315  edm::LogWarning("SiStripRecHitConverter") << "Detid=" << id << " not found";
316  return stripdet != nullptr && (!useQuality || quality->IsModuleUsable(id));
317 }
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Log< level::Warning, false > LogWarning
bool IsModuleUsable(uint32_t detid) const

Member Data Documentation

◆ cpeToken

edm::ESGetToken<StripClusterParameterEstimator, TkStripCPERecord> SiStripRecHitConverterAlgorithm::cpeToken
private

Definition at line 61 of file SiStripRecHitConverterAlgorithm.h.

Referenced by initialize().

◆ doMatching

bool SiStripRecHitConverterAlgorithm::doMatching
private

◆ maskBad128StripBlocks

bool SiStripRecHitConverterAlgorithm::maskBad128StripBlocks
private

Definition at line 59 of file SiStripRecHitConverterAlgorithm.h.

Referenced by fillBad128StripBlocks(), and isMasked().

◆ matcher

const SiStripRecHitMatcher* SiStripRecHitConverterAlgorithm::matcher = nullptr
private

Definition at line 66 of file SiStripRecHitConverterAlgorithm.h.

Referenced by initialize(), and match().

◆ matcherToken

edm::ESGetToken<SiStripRecHitMatcher, TkStripCPERecord> SiStripRecHitConverterAlgorithm::matcherToken
private

◆ parameterestimator

const StripClusterParameterEstimator* SiStripRecHitConverterAlgorithm::parameterestimator = nullptr
private

Definition at line 65 of file SiStripRecHitConverterAlgorithm.h.

Referenced by initialize(), and run().

◆ quality

const SiStripQuality* SiStripRecHitConverterAlgorithm::quality = nullptr
private

Definition at line 67 of file SiStripRecHitConverterAlgorithm.h.

Referenced by fillBad128StripBlocks(), initialize(), and useModule().

◆ qualityToken

edm::ESGetToken<SiStripQuality, SiStripQualityRcd> SiStripRecHitConverterAlgorithm::qualityToken
private

◆ tracker

const TrackerGeometry* SiStripRecHitConverterAlgorithm::tracker = nullptr
private

Definition at line 64 of file SiStripRecHitConverterAlgorithm.h.

Referenced by initialize(), match(), run(), and useModule().

◆ trackerToken

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> SiStripRecHitConverterAlgorithm::trackerToken
private

Definition at line 60 of file SiStripRecHitConverterAlgorithm.h.

Referenced by initialize().

◆ useQuality

bool SiStripRecHitConverterAlgorithm::useQuality
private