CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
Phase2StripCPE Class Referencefinal

#include <Phase2StripCPE.h>

Inheritance diagram for Phase2StripCPE:
ClusterParameterEstimator< Phase2TrackerCluster1D >

Classes

struct  Param
 

Public Types

using Phase2TrackerGeomDetUnit = PixelGeomDetUnit
 
using Phase2TrackerTopology = PixelTopology
 
- Public Types inherited from ClusterParameterEstimator< Phase2TrackerCluster1D >
typedef std::pair< LocalPoint, LocalErrorLocalValues
 
typedef std::vector< LocalValuesVLocalValues
 

Public Member Functions

LocalVector driftDirection (const Phase2TrackerGeomDetUnit &det) const
 
LocalValues localParameters (const Phase2TrackerCluster1D &cluster, const GeomDetUnit &det) const override
 
 Phase2StripCPE (edm::ParameterSet &conf, const MagneticField &, const TrackerGeometry &, const SiPhase2OuterTrackerLorentzAngle &)
 
- Public Member Functions inherited from ClusterParameterEstimator< Phase2TrackerCluster1D >
virtual void clearParameters () const
 
virtual void enterLocalParameters (unsigned int id, std::pair< int, int > &row_col, LocalValues pos_err_info) const
 
virtual void enterLocalParameters (uint32_t id, uint16_t firstStrip, LocalValues pos_err_info) const
 
virtual LocalValues localParameters (const Phase2TrackerCluster1D &cluster, const GeomDetUnit &gd, const LocalTrajectoryParameters &) const
 
virtual LocalValues localParameters (const Phase2TrackerCluster1D &cluster, const GeomDetUnit &gd, const TrajectoryStateOnSurface &tsos) const
 
virtual VLocalValues localParametersV (const Phase2TrackerCluster1D &cluster, const GeomDetUnit &gd) const
 
virtual VLocalValues localParametersV (const Phase2TrackerCluster1D &cluster, const GeomDetUnit &gd, const TrajectoryStateOnSurface &tsos) const
 
virtual ~ClusterParameterEstimator ()
 

Private Member Functions

void fillParam ()
 

Private Attributes

const TrackerGeometrygeom_
 
const SiPhase2OuterTrackerLorentzAnglelorentzAngleMap_
 
unsigned int m_off
 
std::vector< Paramm_Params
 
const MagneticFieldmagfield_
 
float tanLorentzAnglePerTesla_
 
bool use_LorentzAngle_DB_
 

Detailed Description

Definition at line 14 of file Phase2StripCPE.h.

Member Typedef Documentation

◆ Phase2TrackerGeomDetUnit

Definition at line 17 of file Phase2StripCPE.h.

◆ Phase2TrackerTopology

Definition at line 18 of file Phase2StripCPE.h.

Constructor & Destructor Documentation

◆ Phase2StripCPE()

Phase2StripCPE::Phase2StripCPE ( edm::ParameterSet conf,
const MagneticField magf,
const TrackerGeometry geom,
const SiPhase2OuterTrackerLorentzAngle LorentzAngle 
)

Definition at line 5 of file Phase2StripCPE.cc.

References fillParam(), edm::ParameterSet::getParameter(), and use_LorentzAngle_DB_.

9  : magfield_(magf),
10  geom_(geom),
12  tanLorentzAnglePerTesla_(conf.getParameter<double>("TanLorentzAnglePerTesla")) {
13  use_LorentzAngle_DB_ = conf.getParameter<bool>("LorentzAngle_DB");
14  fillParam();
15 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
float tanLorentzAnglePerTesla_
const SiPhase2OuterTrackerLorentzAngle & lorentzAngleMap_
bool use_LorentzAngle_DB_
const TrackerGeometry & geom_
const MagneticField & magfield_

Member Function Documentation

◆ driftDirection()

LocalVector Phase2StripCPE::driftDirection ( const Phase2TrackerGeomDetUnit det) const

Definition at line 29 of file Phase2StripCPE.cc.

References GeomDet::geographicalId(), SiPhase2OuterTrackerLorentzAngle::getLorentzAngle(), MagneticField::inTesla(), lorentzAngleMap_, magfield_, GloballyPositioned< T >::position(), DetId::rawId(), GeomDet::surface(), tanLorentzAnglePerTesla_, toLocal(), use_LorentzAngle_DB_, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by fillParam().

29  {
30  LocalVector lbfield = (det.surface()).toLocal(magfield_.inTesla(det.surface().position()));
31 
32  float langle =
34 
35  float dir_x = -langle * lbfield.y();
36  float dir_y = langle * lbfield.x();
37  float dir_z = 1.f; // E field always in z direction
38 
39  return LocalVector(dir_x, dir_y, dir_z);
40 }
Local3DVector LocalVector
Definition: LocalVector.h:12
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
float tanLorentzAnglePerTesla_
const SiPhase2OuterTrackerLorentzAngle & lorentzAngleMap_
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)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const PositionType & position() const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
bool use_LorentzAngle_DB_
const MagneticField & magfield_

◆ fillParam()

void Phase2StripCPE::fillParam ( )
private

Definition at line 42 of file Phase2StripCPE.cc.

References cms::cuda::assert(), Surface::bounds(), ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), TrackerGeometry::detUnits(), shallow::drift(), driftDirection(), geom_, mps_fire::i, GeomDet::index(), createfilelist::int, LogDebug, m_off, m_Params, TrackerGeometry::offsetDU(), AlCaHLTBitMon_ParallelJobs::p, GeomDet::specificSurface(), PixelGeomDetUnit::specificTopology(), Bounds::thickness(), Calorimetry_cff::thickness, and GeomDetEnumerators::tkDetEnum.

Referenced by Phase2StripCPE().

42  {
43  // in phase 2 they are all pixel topologies...
44  auto const& dus = geom_.detUnits();
45  m_off = dus.size();
46  // skip Barrel and Foward pixels...
47  for (unsigned int i = 3; i < 7; ++i) {
48  LogDebug("LookingForFirstPhase2OT") << " Subdetector " << i << " GeomDetEnumerator "
49  << GeomDetEnumerators::tkDetEnum[i] << " offset "
51  if (geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i]) != dus.size()) {
54  }
55  }
56  LogDebug("LookingForFirstPhase2OT") << " Chosen offset: " << m_off;
57 
58  m_Params.resize(dus.size() - m_off);
59  // very very minimal, for sure it will need to expand...
60  for (auto i = m_off; i != dus.size(); ++i) {
61  auto& p = m_Params[i - m_off];
62 
63  const Phase2TrackerGeomDetUnit& det = (const Phase2TrackerGeomDetUnit&)(*dus[i]);
64  assert(det.index() == int(i));
65  p.topology = &det.specificTopology();
66 
67  auto pitch_x = p.topology->pitch().first;
68  auto pitch_y = p.topology->pitch().second;
69 
70  // see https://github.com/cms-sw/cmssw/blob/CMSSW_8_1_X/RecoLocalTracker/SiStripRecHitConverter/src/StripCPE.cc
71  auto thickness = det.specificSurface().bounds().thickness();
72  auto drift = driftDirection(det) * thickness;
73  auto lvec = drift + LocalVector(0, 0, -thickness);
74  p.coveredStrips = lvec.x() / pitch_x; // simplifies wrt Phase0 tracker because only rectangular modules
75 
76  constexpr float o12 = 1. / 12;
77  p.localErr = LocalError(o12 * pitch_x * pitch_x, 0, o12 * pitch_y * pitch_y); // e2_xx, e2_xy, e2_yy
78  }
79 }
std::vector< Param > m_Params
Local3DVector LocalVector
Definition: LocalVector.h:12
int index() const
Definition: GeomDet.h:83
unsigned int offsetDU(SubDetector sid) const
LocalVector drift(const StripGeomDetUnit *, const MagneticField &, const SiStripLorentzAngle &)
Definition: ShallowTools.cc:36
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
assert(be >=bs)
virtual float thickness() const =0
LocalVector driftDirection(const Phase2TrackerGeomDetUnit &det) const
unsigned int m_off
const TrackerGeometry & geom_
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
constexpr SubDetector tkDetEnum[8]
#define LogDebug(id)
const Bounds & bounds() const
Definition: Surface.h:87

◆ localParameters()

Phase2StripCPE::LocalValues Phase2StripCPE::localParameters ( const Phase2TrackerCluster1D cluster,
const GeomDetUnit det 
) const
overridevirtual

Implements ClusterParameterEstimator< Phase2TrackerCluster1D >.

Definition at line 17 of file Phase2StripCPE.cc.

References Phase2TrackerCluster1D::center(), Phase2TrackerCluster1D::column(), f, nano_mu_digi_cff::float, GeomDet::index(), ALPAKA_ACCELERATOR_NAMESPACE::ecal::reconstruction::internal::endcap::ix(), ALPAKA_ACCELERATOR_NAMESPACE::ecal::reconstruction::internal::endcap::iy(), m_off, m_Params, and AlCaHLTBitMon_ParallelJobs::p.

18  {
19  auto const& p = m_Params[detunit.index() - m_off];
20  auto const& topo = *p.topology;
21  float ix = cluster.center() - 0.5f * p.coveredStrips;
22  float iy = float(cluster.column()) + 0.5f; // halfway the column
23 
24  LocalPoint lp(topo.localX(ix), topo.localY(iy), 0); // x, y, z
25 
26  return std::make_pair(lp, p.localErr);
27 }
std::vector< Param > m_Params
double f[11][100]
unsigned int m_off
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
unsigned int column() const
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)

Member Data Documentation

◆ geom_

const TrackerGeometry& Phase2StripCPE::geom_
private

Definition at line 40 of file Phase2StripCPE.h.

Referenced by fillParam().

◆ lorentzAngleMap_

const SiPhase2OuterTrackerLorentzAngle& Phase2StripCPE::lorentzAngleMap_
private

Definition at line 41 of file Phase2StripCPE.h.

Referenced by driftDirection().

◆ m_off

unsigned int Phase2StripCPE::m_off
private

Definition at line 44 of file Phase2StripCPE.h.

Referenced by fillParam(), and localParameters().

◆ m_Params

std::vector<Param> Phase2StripCPE::m_Params
private

Definition at line 37 of file Phase2StripCPE.h.

Referenced by fillParam(), and localParameters().

◆ magfield_

const MagneticField& Phase2StripCPE::magfield_
private

Definition at line 39 of file Phase2StripCPE.h.

Referenced by driftDirection().

◆ tanLorentzAnglePerTesla_

float Phase2StripCPE::tanLorentzAnglePerTesla_
private

Definition at line 43 of file Phase2StripCPE.h.

Referenced by driftDirection().

◆ use_LorentzAngle_DB_

bool Phase2StripCPE::use_LorentzAngle_DB_
private

Definition at line 46 of file Phase2StripCPE.h.

Referenced by driftDirection(), and Phase2StripCPE().