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 &)
 
- 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_
 
unsigned int m_off
 
std::vector< Paramm_Params
 
const MagneticFieldmagfield_
 
float tanLorentzAnglePerTesla_
 
bool use_LorentzAngle_DB_
 

Detailed Description

Definition at line 13 of file Phase2StripCPE.h.

Member Typedef Documentation

Definition at line 16 of file Phase2StripCPE.h.

Definition at line 17 of file Phase2StripCPE.h.

Constructor & Destructor Documentation

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

Definition at line 5 of file Phase2StripCPE.cc.

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

6  : magfield_(magf), geom_(geom), tanLorentzAnglePerTesla_(conf.getParameter<double>("TanLorentzAnglePerTesla")) {
7  use_LorentzAngle_DB_ = conf.getParameter<bool>("LorentzAngle_DB");
9  throw cms::Exception("Lorentz Angle from DB not implemented yet");
11  // old code: LorentzAngleMap_.getLorentzAngle(det->geographicalId().rawId());
12  }
13  fillParam();
14 }
T getParameter(std::string const &) const
float tanLorentzAnglePerTesla_
bool use_LorentzAngle_DB_
const TrackerGeometry & geom_
const MagneticField & magfield_

Member Function Documentation

LocalVector Phase2StripCPE::driftDirection ( const Phase2TrackerGeomDetUnit det) const

Definition at line 28 of file Phase2StripCPE.cc.

References MagneticField::inTesla(), magfield_, GloballyPositioned< T >::position(), GeomDet::surface(), tanLorentzAnglePerTesla_, toLocal(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by fillParam().

28  {
29  LocalVector lbfield = (det.surface()).toLocal(magfield_.inTesla(det.surface().position()));
30 
31  float dir_x = -tanLorentzAnglePerTesla_ * lbfield.y();
32  float dir_y = tanLorentzAnglePerTesla_ * lbfield.x();
33  float dir_z = 1.f; // E field always in z direction
34 
35  return LocalVector(dir_x, dir_y, dir_z);
36 }
Local3DVector LocalVector
Definition: LocalVector.h:12
float tanLorentzAnglePerTesla_
T y() const
Definition: PV3DBase.h:60
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
T x() const
Definition: PV3DBase.h:59
const PositionType & position() const
const MagneticField & magfield_
void Phase2StripCPE::fillParam ( )
private

Definition at line 38 of file Phase2StripCPE.cc.

References Surface::bounds(), 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().

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

Implements ClusterParameterEstimator< Phase2TrackerCluster1D >.

Definition at line 16 of file Phase2StripCPE.cc.

References Phase2TrackerCluster1D::center(), Phase2TrackerCluster1D::column(), f, dqmMemoryStats::float, GeomDet::index(), m_off, m_Params, and AlCaHLTBitMon_ParallelJobs::p.

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

Member Data Documentation

const TrackerGeometry& Phase2StripCPE::geom_
private

Definition at line 36 of file Phase2StripCPE.h.

Referenced by fillParam().

unsigned int Phase2StripCPE::m_off
private

Definition at line 38 of file Phase2StripCPE.h.

Referenced by fillParam(), and localParameters().

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

Definition at line 33 of file Phase2StripCPE.h.

Referenced by fillParam(), and localParameters().

const MagneticField& Phase2StripCPE::magfield_
private

Definition at line 35 of file Phase2StripCPE.h.

Referenced by driftDirection().

float Phase2StripCPE::tanLorentzAnglePerTesla_
private

Definition at line 37 of file Phase2StripCPE.h.

Referenced by driftDirection(), and Phase2StripCPE().

bool Phase2StripCPE::use_LorentzAngle_DB_
private

Definition at line 40 of file Phase2StripCPE.h.

Referenced by Phase2StripCPE().