CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 &, edm::InputTag &)
 

Typedef Documentation

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

Definition at line 16 of file ShallowTools.h.

Function Documentation

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

Definition at line 39 of file ShallowTools.cc.

References Surface::bounds(), GeomDet::geographicalId(), SiStripLorentzAngle::getLorentzAngle(), MagneticField::inTesla(), 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(), DTLinearDriftAlgo::compute(), DTLinearDriftFromDBAlgo::compute(), DTNoDriftAlgo::compute(), ChargeDrifterFP420::drift(), ClusterShapeHitFilter::getSizes(), ShallowTrackClustersProducer::produce(), ShallowSimhitClustersProducer::produce(), SiStripTrackingRecHitsValid::rechitanalysis(), and SiStripTrackingRecHitsValid::rechitanalysis_matched().

39  {
40  LocalVector lbfield=( stripDet->surface()).toLocal( magfield.inTesla(stripDet->surface().position()));
41  float tanLorentzAnglePerTesla = lorentzAngle.getLorentzAngle(stripDet->geographicalId());
42  float driftz = stripDet->specificSurface().bounds().thickness();
43  float driftx =-tanLorentzAnglePerTesla * lbfield.y() * driftz;
44  float drifty = tanLorentzAnglePerTesla * lbfield.x() * driftz;
45  return LocalVector(driftx,drifty,driftz);
46 }
Local3DVector LocalVector
Definition: LocalVector.h:12
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
T y() const
Definition: PV3DBase.h:63
const Bounds & bounds() const
Definition: Surface.h:128
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
virtual float thickness() const =0
float getLorentzAngle(const uint32_t &) const
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:77
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:43
int shallow::findTrackIndex ( const edm::Handle< edm::View< reco::Track > > &  h,
const reco::Track t 
)

Definition at line 29 of file ShallowTools.cc.

References end, and testEve_cfg::tracks.

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

29  {
31  it = tracks->begin(),
32  end = tracks->end();
33  //Compare addresses
34  for(; it!=end; it++) { if (&(*it)==track) { return it - tracks->begin(); } }
35  return -2;
36 }
#define end
Definition: vmac.h:37
tuple tracks
Definition: testEve_cfg.py:39
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
CLUSTERMAP shallow::make_cluster_map ( const edm::Event iEvent,
edm::InputTag clusterLabel 
)

Definition at line 15 of file ShallowTools.cc.

References HLT_25ns14e33_v1_cff::clusters, edmNew::DetSet< T >::detId(), SiStripCluster::firstStrip(), and edm::Event::getByLabel().

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

15  {
16  CLUSTERMAP clustermap;
18  iEvent.getByLabel(clusterLabel, clusters);
19 
20  unsigned int clusterindex = 0;
21  BOOST_FOREACH(const edmNew::DetSet<SiStripCluster>& ds, *clusters)
22  BOOST_FOREACH(const SiStripCluster& cluster, ds)
23  clustermap.insert( std::make_pair( std::make_pair(ds.detId(),cluster.firstStrip()),
24  clusterindex++));
25  return clustermap;
26 }
return((rh^lh)&mask)
bool insert(Storage &iStorage, ItemType *iItem, const IdTag &iIdTag)
Definition: HCMethods.h:49
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
string const
Definition: compareJSON.py:14
std::map< std::pair< uint32_t, uint16_t >, unsigned int > CLUSTERMAP
Definition: ShallowTools.h:16