CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes
StripCPEfromTrackAngle Class Reference

#include <StripCPEfromTrackAngle.h>

Inheritance diagram for StripCPEfromTrackAngle:
StripCPE StripClusterParameterEstimator

Public Types

using AClusters = StripClusterParameterEstimator::AClusters
 
using AlgoParam = StripCPE::AlgoParam
 
using ALocalValues = StripClusterParameterEstimator::ALocalValues
 
- Public Types inherited from StripClusterParameterEstimator
using AClusters = DynArray< SiStripCluster const * >
 
using ALocalValues = DynArray< LocalValues >
 
using LocalValues = std::pair< LocalPoint, LocalError >
 
typedef std::vector< LocalValuesVLocalValues
 

Public Member Functions

float legacyStripErrorSquared (const unsigned N, const float uProj) const
 
void localParameters (AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const override
 
StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &cl, AlgoParam const &ap) const override
 
StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &, const GeomDetUnit &, const LocalTrajectoryParameters &) const override
 
 StripCPEfromTrackAngle (edm::ParameterSet &conf, const MagneticField &mag, const TrackerGeometry &geom, const SiStripLorentzAngle &lorentz, const SiStripBackPlaneCorrection &backPlaneCorrection, const SiStripConfObject &confObj, const SiStripLatency &latency)
 
float stripErrorSquared (const unsigned N, const float uProj, const SiStripDetId::SubDetector loc) const
 
- Public Member Functions inherited from StripCPE
LocalVector driftDirection (const StripGeomDetUnit *det) const override
 
AlgoParam getAlgoParam (const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
 
virtual LocalValues localParameters (const SiStripCluster &cluster, const GeomDetUnit &gd, const TrajectoryStateOnSurface &tsos) const
 
virtual LocalValues localParameters (const SiStripCluster &, const GeomDetUnit &) const
 
virtual LocalValues localParameters (const SiStripCluster &cluster, const GeomDetUnit &gd, const LocalTrajectoryParameters &) const
 
virtual void localParameters (AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const
 
StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &cl, const GeomDetUnit &) const override
 
 StripCPE (edm::ParameterSet &conf, const MagneticField &, const TrackerGeometry &, const SiStripLorentzAngle &, const SiStripBackPlaneCorrection &, const SiStripConfObject &, const SiStripLatency &)
 
- Public Member Functions inherited from StripClusterParameterEstimator
virtual LocalValues localParameters (const SiStripCluster &cluster, const GeomDetUnit &gd, const TrajectoryStateOnSurface &tsos) const
 
virtual VLocalValues localParametersV (const SiStripCluster &cluster, const GeomDetUnit &gd) const
 
virtual VLocalValues localParametersV (const SiStripCluster &cluster, const GeomDetUnit &gd, const TrajectoryStateOnSurface &tsos) const
 
virtual ~StripClusterParameterEstimator ()
 

Private Types

enum  Algo { Algo::legacy, Algo::mergeCK, Algo::chargeCK }
 

Private Member Functions

StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &cl, const GeomDetUnit &) const override
 
virtual LocalValues localParameters (const SiStripCluster &cluster, const GeomDetUnit &gd, const TrajectoryStateOnSurface &tsos) const
 
virtual LocalValues localParameters (const SiStripCluster &cluster, const GeomDetUnit &gd, const LocalTrajectoryParameters &) const
 
virtual LocalValues localParameters (const SiStripCluster &, const GeomDetUnit &) const
 
virtual void localParameters (AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const
 
virtual StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &cl, AlgoParam const &ap) const
 

Private Attributes

Algo m_algo
 
const float maxChgOneMIP
 
float mHC_P [4][2]
 
std::map< SiStripDetId::SubDetector, float > mHC_P0
 
std::map< SiStripDetId::SubDetector, float > mHC_P1
 
float mLC_P [3]
 
const bool useLegacyError
 

Additional Inherited Members

- Protected Member Functions inherited from StripCPE
Param const & param (const GeomDetUnit &det) const
 
- Protected Attributes inherited from StripCPE
const SiStripBackPlaneCorrectionBackPlaneCorrectionMap_
 
const TrackerGeometrygeom_
 
const SiStripLorentzAngleLorentzAngleMap_
 
const MagneticFieldmagfield_
 
const bool peakMode_
 
std::vector< float > xtalk1
 
std::vector< float > xtalk2
 

Detailed Description

Definition at line 6 of file StripCPEfromTrackAngle.h.

Member Typedef Documentation

◆ AClusters

Definition at line 31 of file StripCPEfromTrackAngle.h.

◆ AlgoParam

Definition at line 30 of file StripCPEfromTrackAngle.h.

◆ ALocalValues

Definition at line 32 of file StripCPEfromTrackAngle.h.

Member Enumeration Documentation

◆ Algo

enum StripCPEfromTrackAngle::Algo
strongprivate
Enumerator
legacy 
mergeCK 
chargeCK 

Definition at line 25 of file StripCPEfromTrackAngle.h.

25 { legacy, mergeCK, chargeCK };

Constructor & Destructor Documentation

◆ StripCPEfromTrackAngle()

StripCPEfromTrackAngle::StripCPEfromTrackAngle ( edm::ParameterSet conf,
const MagneticField mag,
const TrackerGeometry geom,
const SiStripLorentzAngle lorentz,
const SiStripBackPlaneCorrection backPlaneCorrection,
const SiStripConfObject confObj,
const SiStripLatency latency 
)

Definition at line 7 of file StripCPEfromTrackAngle.cc.

References edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), mHC_P, mLC_P, SiStripDetId::TEC, SiStripDetId::TIB, SiStripDetId::TID, and SiStripDetId::TOB.

14  : StripCPE(conf, mag, geom, lorentz, backPlaneCorrection, confObj, latency),
15  useLegacyError(conf.existsAs<bool>("useLegacyError") ? conf.getParameter<bool>("useLegacyError") : true),
16  maxChgOneMIP(conf.existsAs<float>("maxChgOneMIP") ? conf.getParameter<double>("maxChgOneMIP") : -6000.),
18  mLC_P[0] = conf.existsAs<double>("mLC_P0") ? conf.getParameter<double>("mLC_P0") : -.326;
19  mLC_P[1] = conf.existsAs<double>("mLC_P1") ? conf.getParameter<double>("mLC_P1") : .618;
20  mLC_P[2] = conf.existsAs<double>("mLC_P2") ? conf.getParameter<double>("mLC_P2") : .300;
21 
22  mHC_P[SiStripDetId::TIB - 3][0] = conf.existsAs<double>("mTIB_P0") ? conf.getParameter<double>("mTIB_P0") : -.742;
23  mHC_P[SiStripDetId::TIB - 3][1] = conf.existsAs<double>("mTIB_P1") ? conf.getParameter<double>("mTIB_P1") : .202;
24  mHC_P[SiStripDetId::TID - 3][0] = conf.existsAs<double>("mTID_P0") ? conf.getParameter<double>("mTID_P0") : -1.026;
25  mHC_P[SiStripDetId::TID - 3][1] = conf.existsAs<double>("mTID_P1") ? conf.getParameter<double>("mTID_P1") : .253;
26  mHC_P[SiStripDetId::TOB - 3][0] = conf.existsAs<double>("mTOB_P0") ? conf.getParameter<double>("mTOB_P0") : -1.427;
27  mHC_P[SiStripDetId::TOB - 3][1] = conf.existsAs<double>("mTOB_P1") ? conf.getParameter<double>("mTOB_P1") : .433;
28  mHC_P[SiStripDetId::TEC - 3][0] = conf.existsAs<double>("mTEC_P0") ? conf.getParameter<double>("mTEC_P0") : -1.885;
29  mHC_P[SiStripDetId::TEC - 3][1] = conf.existsAs<double>("mTEC_P1") ? conf.getParameter<double>("mTEC_P1") : .471;
30 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
StripCPE(edm::ParameterSet &conf, const MagneticField &, const TrackerGeometry &, const SiStripLorentzAngle &, const SiStripBackPlaneCorrection &, const SiStripConfObject &, const SiStripLatency &)
Definition: StripCPE.cc:9
static constexpr auto TID
Definition: SiStripDetId.h:39
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:172
static constexpr auto TOB
Definition: SiStripDetId.h:40
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
static constexpr auto TIB
Definition: SiStripDetId.h:38
static constexpr auto TEC
Definition: SiStripDetId.h:41

Member Function Documentation

◆ legacyStripErrorSquared()

float StripCPEfromTrackAngle::legacyStripErrorSquared ( const unsigned  N,
const float  uProj 
) const

Definition at line 40 of file StripCPEfromTrackAngle.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), f, myMath::fast_expf(), ALCARECOEcalPhiSym_cff::float, N, and UNLIKELY.

Referenced by localParameters().

40  {
41  if UNLIKELY ((float(N) - uProj) > 3.5f)
42  return float(N * N) / 12.f;
43  else {
44  static constexpr float P1 = -0.339;
45  static constexpr float P2 = 0.90;
46  static constexpr float P3 = 0.279;
47  const float uerr = P1 * uProj * vdt::fast_expf(-uProj * P2) + P3;
48  return uerr * uerr;
49  }
50 }
double f[11][100]
#define N
Definition: blowfish.cc:9
float fast_expf(float x)
#define UNLIKELY(x)
Definition: Likely.h:21

◆ localParameters() [1/9]

virtual LocalValues StripClusterParameterEstimator::localParameters
inlineprivate

Definition at line 45 of file StripClusterParameterEstimator.h.

47  {
48  return localParameters(cluster, gd, tsos.localParameters());
49  }
void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const override

◆ localParameters() [2/9]

virtual void StripClusterParameterEstimator::localParameters
inlineprivate

Definition at line 32 of file StripClusterParameterEstimator.h.

35  {}

◆ localParameters() [3/9]

virtual LocalValues StripClusterParameterEstimator::localParameters
inlineprivate

Definition at line 37 of file StripClusterParameterEstimator.h.

37  {
38  return std::make_pair(LocalPoint(), LocalError());
39  }
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30

◆ localParameters() [4/9]

virtual LocalValues StripClusterParameterEstimator::localParameters
inlineprivate

Definition at line 40 of file StripClusterParameterEstimator.h.

42  {
43  return localParameters(cluster, gd);
44  }
void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const override

◆ localParameters() [5/9]

virtual StripClusterParameterEstimator::LocalValues StripCPE::localParameters
inlineprivate

Definition at line 52 of file StripCPE.h.

53  {
54  return std::make_pair(LocalPoint(), LocalError());
55  }
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30

◆ localParameters() [6/9]

StripClusterParameterEstimator::LocalValues StripCPE::localParameters
overrideprivate

Definition at line 62 of file StripCPE.cc.

63  {
64  StripCPE::Param const& p = param(det);
65  const float barycenter = cluster.barycenter();
66  const float fullProjection =
67  p.coveredStrips(p.drift + LocalVector(0, 0, -p.thickness), p.topology->localPosition(barycenter));
68  const float strip = barycenter - 0.5f * (1.f - p.backplanecorrection) * fullProjection;
69 
70  return std::make_pair(p.topology->localPosition(strip), p.topology->localError(strip, 1.f / 12.f));
71 }
Param const & param(const GeomDetUnit &det) const
Definition: StripCPE.h:81
double f[11][100]

◆ localParameters() [7/9]

void StripCPEfromTrackAngle::localParameters ( AClusters const &  clusters,
ALocalValues retValues,
const GeomDetUnit gd,
const LocalTrajectoryParameters ltp 
) const
overridevirtual

Reimplemented from StripClusterParameterEstimator.

Definition at line 52 of file StripCPEfromTrackAngle.cc.

References chargeCK, siStripClusterTools::chargePerCM(), bsc_activity_cfg::clusters, alignCSCRings::corr, EcalPhiSymFlatTableProducers_cfi::fill, StripCPE::getAlgoParam(), mps_fire::i, legacy, legacyStripErrorSquared(), isotrackTrainRegressor::loc, m_algo, maxChgOneMIP, mergeCK, N, AlCaHLTBitMon_ParallelJobs::p, nano_mu_digi_cff::strip, stripErrorSquared(), mitigatedMETSequence_cff::U, and LocalTrajectoryParameters::vector().

Referenced by localParameters().

55  {
56  auto const& par = getAlgoParam(det, ltp);
57  auto const& p = par.p;
58  auto loc = par.loc;
59  auto corr = par.corr;
60  auto afp = par.afullProjection;
61 
62  auto fill = [&](unsigned int i, float uerr2) {
63  const float strip = clusters[i]->barycenter() + corr;
64  retValues[i].first = p.topology->localPosition(strip, ltp.vector());
65  retValues[i].second = p.topology->localError(strip, uerr2, ltp.vector());
66  };
67 
68  switch (m_algo) {
69  case Algo::chargeCK:
70  for (auto i = 0U; i < clusters.size(); ++i) {
71  auto dQdx = siStripClusterTools::chargePerCM(*clusters[i], ltp, p.invThickness);
72  auto N = clusters[i]->amplitudes().size();
73  auto uerr2 = dQdx > maxChgOneMIP ? legacyStripErrorSquared(N, afp) : stripErrorSquared(N, afp, loc);
74  fill(i, uerr2);
75  }
76  break;
77  case Algo::legacy:
78  for (auto i = 0U; i < clusters.size(); ++i) {
79  auto N = clusters[i]->amplitudes().size();
80  auto uerr2 = legacyStripErrorSquared(N, afp);
81  fill(i, uerr2);
82  }
83  break;
84  case Algo::mergeCK:
85  for (auto i = 0U; i < clusters.size(); ++i) {
86  auto N = clusters[i]->amplitudes().size();
87  auto uerr2 = clusters[i]->isMerged() ? legacyStripErrorSquared(N, afp) : stripErrorSquared(N, afp, loc);
88  fill(i, uerr2);
89  }
90  break;
91  }
92 }
float chargePerCM(DetId detid, Iter a, Iter b)
float legacyStripErrorSquared(const unsigned N, const float uProj) const
AlgoParam getAlgoParam(const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
Definition: StripCPE.h:57
dictionary corr
AlgebraicVector5 vector() const
#define N
Definition: blowfish.cc:9
float stripErrorSquared(const unsigned N, const float uProj, const SiStripDetId::SubDetector loc) const

◆ localParameters() [8/9]

StripClusterParameterEstimator::LocalValues StripCPEfromTrackAngle::localParameters ( const SiStripCluster cl,
AlgoParam const &  ap 
) const
overridevirtual

Reimplemented from StripCPE.

Definition at line 94 of file StripCPEfromTrackAngle.cc.

References StripCPE::AlgoParam::afullProjection, SiStripCluster::amplitudes(), SiStripCluster::barycenter(), chargeCK, siStripClusterTools::chargePerCM(), StripCPE::AlgoParam::corr, alignCSCRings::corr, SiStripCluster::isMerged(), legacy, legacyStripErrorSquared(), StripCPE::AlgoParam::loc, isotrackTrainRegressor::loc, StripCPE::AlgoParam::ltp, m_algo, maxChgOneMIP, mergeCK, N, StripCPE::AlgoParam::p, AlCaHLTBitMon_ParallelJobs::p, SiStripCluster::size(), nano_mu_digi_cff::strip, and stripErrorSquared().

95  {
96  auto const& p = par.p;
97  auto const& ltp = par.ltp;
98  auto loc = par.loc;
99  auto corr = par.corr;
100  auto afp = par.afullProjection;
101 
102  float uerr2 = 0;
103 
104  auto N = cluster.amplitudes().size();
105 
106  switch (m_algo) {
107  case Algo::chargeCK: {
108  auto dQdx = siStripClusterTools::chargePerCM(cluster, ltp, p.invThickness);
109  uerr2 = dQdx > maxChgOneMIP ? legacyStripErrorSquared(N, afp) : stripErrorSquared(N, afp, loc);
110  } break;
111  case Algo::legacy:
112  uerr2 = legacyStripErrorSquared(N, afp);
113  break;
114  case Algo::mergeCK:
115  uerr2 = cluster.isMerged() ? legacyStripErrorSquared(N, afp) : stripErrorSquared(N, afp, loc);
116  break;
117  }
118 
119  const float strip = cluster.barycenter() + corr;
120 
121  return std::make_pair(p.topology->localPosition(strip, ltp.vector()),
122  p.topology->localError(strip, uerr2, ltp.vector()));
123 }
float chargePerCM(DetId detid, Iter a, Iter b)
float legacyStripErrorSquared(const unsigned N, const float uProj) const
dictionary corr
#define N
Definition: blowfish.cc:9
float stripErrorSquared(const unsigned N, const float uProj, const SiStripDetId::SubDetector loc) const

◆ localParameters() [9/9]

StripClusterParameterEstimator::LocalValues StripCPEfromTrackAngle::localParameters ( const SiStripCluster cluster,
const GeomDetUnit det,
const LocalTrajectoryParameters ltp 
) const
overridevirtual

Reimplemented from StripClusterParameterEstimator.

Definition at line 125 of file StripCPEfromTrackAngle.cc.

References StripCPE::getAlgoParam(), and localParameters().

126  {
127  auto const& par = getAlgoParam(det, ltp);
128  return localParameters(cluster, par);
129 }
void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const override
AlgoParam getAlgoParam(const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
Definition: StripCPE.h:57

◆ stripErrorSquared()

float StripCPEfromTrackAngle::stripErrorSquared ( const unsigned  N,
const float  uProj,
const SiStripDetId::SubDetector  loc 
) const

Definition at line 32 of file StripCPEfromTrackAngle.cc.

References myMath::fast_expf(), ALCARECOEcalPhiSym_cff::float, isotrackTrainRegressor::loc, mHC_P, mLC_P, N, and x.

Referenced by localParameters().

34  {
35  auto fun = [&](float x) -> float { return mLC_P[0] * x * vdt::fast_expf(-x * mLC_P[1]) + mLC_P[2]; };
36  auto uerr = (N <= 4) ? fun(uProj) : mHC_P[loc - 3][0] + float(N) * mHC_P[loc - 3][1];
37  return uerr * uerr;
38 }
#define N
Definition: blowfish.cc:9
float fast_expf(float x)

Member Data Documentation

◆ m_algo

Algo StripCPEfromTrackAngle::m_algo
private

Definition at line 27 of file StripCPEfromTrackAngle.h.

Referenced by localParameters().

◆ maxChgOneMIP

const float StripCPEfromTrackAngle::maxChgOneMIP
private

Definition at line 23 of file StripCPEfromTrackAngle.h.

Referenced by localParameters().

◆ mHC_P

float StripCPEfromTrackAngle::mHC_P[4][2]
private

Definition at line 12 of file StripCPEfromTrackAngle.h.

Referenced by StripCPEfromTrackAngle(), and stripErrorSquared().

◆ mHC_P0

std::map<SiStripDetId::SubDetector, float> StripCPEfromTrackAngle::mHC_P0
private

Definition at line 15 of file StripCPEfromTrackAngle.h.

◆ mHC_P1

std::map<SiStripDetId::SubDetector, float> StripCPEfromTrackAngle::mHC_P1
private

Definition at line 16 of file StripCPEfromTrackAngle.h.

◆ mLC_P

float StripCPEfromTrackAngle::mLC_P[3]
private

Definition at line 11 of file StripCPEfromTrackAngle.h.

Referenced by StripCPEfromTrackAngle(), and stripErrorSquared().

◆ useLegacyError

const bool StripCPEfromTrackAngle::useLegacyError
private

Definition at line 19 of file StripCPEfromTrackAngle.h.