CMS 3D CMS Logo

Typedefs | Functions
shallow Namespace Reference

Typedefs

typedef std::map< std::pair< uint32_t, uint16_t >, unsigned int > CLUSTERMAP
 

Functions

LocalVector drift (const StripGeomDetUnit *, const MagneticField &, const SiStripLorentzAngle &)
 
int findTrackIndex (const edm::Handle< edm::View< reco::Track > > &h, const reco::Track *t)
 
CLUSTERMAP make_cluster_map (const edm::Event &, const edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > &)
 

Typedef Documentation

◆ CLUSTERMAP

typedef std::map<std::pair<uint32_t, uint16_t>, unsigned int> shallow::CLUSTERMAP

Definition at line 21 of file ShallowTools.h.

Function Documentation

◆ drift()

LocalVector shallow::drift ( const StripGeomDetUnit stripDet,
const MagneticField magfield,
const SiStripLorentzAngle lorentzAngle 
)

Definition at line 36 of file ShallowTools.cc.

References Surface::bounds(), GeomDet::geographicalId(), SiStripLorentzAngle::getLorentzAngle(), volumeBasedMagneticField_160812_cfi::magfield, GloballyPositioned< T >::position(), GeomDet::specificSurface(), GeomDet::surface(), Bounds::thickness(), toLocal(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by analyzer::SiPixelLorentzAngle::analyze(), SiPixelLorentzAnglePCLWorker::analyze(), DTLinearDriftAlgo::compute(), DTNoDriftAlgo::compute(), DTLinearDriftFromDBAlgo::compute(), DTParametrizedDriftAlgo::compute(), SiLinearChargeDivider::driftXPos(), Phase2StripCPE::fillParam(), ClusterShapeHitFilter::getSizes(), SiStripLorentzAngleRunInfoTableProducer::globalBeginRunProduce(), ShallowSimhitClustersProducer::produce(), ShallowTrackClustersProducer::produce(), SiStripTrackingRecHitsValid::rechitanalysis(), and SiStripTrackingRecHitsValid::rechitanalysis_matched().

38  {
39  LocalVector lbfield = (stripDet->surface()).toLocal(magfield.inTesla(stripDet->surface().position()));
40  float tanLorentzAnglePerTesla = lorentzAngle.getLorentzAngle(stripDet->geographicalId());
41  float driftz = stripDet->specificSurface().bounds().thickness();
42  float driftx = -tanLorentzAnglePerTesla * lbfield.y() * driftz;
43  float drifty = tanLorentzAnglePerTesla * lbfield.x() * driftz;
44  return LocalVector(driftx, drifty, driftz);
45  }
Local3DVector LocalVector
Definition: LocalVector.h:12
virtual float thickness() const =0
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
float getLorentzAngle(const uint32_t &) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const PositionType & position() const
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
const Bounds & bounds() const
Definition: Surface.h:87

◆ findTrackIndex()

int shallow::findTrackIndex ( const edm::Handle< edm::View< reco::Track > > &  h,
const reco::Track t 
)

Definition at line 25 of file ShallowTools.cc.

References mps_fire::end, HLT_2022v14_cff::track, and tracks.

Referenced by ShallowSimTracksProducer::produce(), SiStripOnTrackClusterTableProducerBase::produce(), ShallowTrackClustersProducer::produce(), and ShallowGainCalibration::produce().

25  {
27  //Compare addresses
28  for (; it != end; it++) {
29  if (&(*it) == track) {
30  return it - tracks->begin();
31  }
32  }
33  return -2;
34  }
auto const & tracks
cannot be loose
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86

◆ make_cluster_map()

CLUSTERMAP shallow::make_cluster_map ( const edm::Event iEvent,
const edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > &  cluster_token 
)

Definition at line 12 of file ShallowTools.cc.

References bsc_activity_cfg::clusters, and iEvent.

Referenced by ShallowRechitClustersProducer::produce(), ShallowSimhitClustersProducer::produce(), and ShallowTrackClustersProducer::produce().

13  {
14  CLUSTERMAP clustermap;
16  iEvent.getByToken(cluster_token, clusters);
17 
18  unsigned int clusterindex = 0;
19  for (auto const& ds : *clusters)
20  for (auto const& cluster : ds)
21  clustermap.insert(std::make_pair(std::make_pair(ds.detId(), cluster.firstStrip()), clusterindex++));
22  return clustermap;
23  }
int iEvent
Definition: GenABIO.cc:224
std::map< std::pair< uint32_t, uint16_t >, unsigned int > CLUSTERMAP
Definition: ShallowTools.h:21