CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Types | 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
 
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 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

Definition at line 33 of file StripCPEfromTrackAngle.h.

Definition at line 32 of file StripCPEfromTrackAngle.h.

Definition at line 34 of file StripCPEfromTrackAngle.h.

Member Enumeration Documentation

enum StripCPEfromTrackAngle::Algo
strongprivate
Enumerator
legacy 
mergeCK 
chargeCK 

Definition at line 27 of file StripCPEfromTrackAngle.h.

27 { legacy, mergeCK, chargeCK };

Constructor & Destructor Documentation

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  {
19  mLC_P[0] = conf.existsAs<double>("mLC_P0") ? conf.getParameter<double>("mLC_P0") : -.326;
20  mLC_P[1] = conf.existsAs<double>("mLC_P1") ? conf.getParameter<double>("mLC_P1") : .618;
21  mLC_P[2] = conf.existsAs<double>("mLC_P2") ? conf.getParameter<double>("mLC_P2") : .300;
22 
23  mHC_P[SiStripDetId::TIB - 3][0] = conf.existsAs<double>("mTIB_P0") ? conf.getParameter<double>("mTIB_P0") : -.742 ;
24  mHC_P[SiStripDetId::TIB - 3][1] = conf.existsAs<double>("mTIB_P1") ? conf.getParameter<double>("mTIB_P1") : .202 ;
25  mHC_P[SiStripDetId::TID - 3][0] = conf.existsAs<double>("mTID_P0") ? conf.getParameter<double>("mTID_P0") : -1.026 ;
26  mHC_P[SiStripDetId::TID - 3][1] = conf.existsAs<double>("mTID_P1") ? conf.getParameter<double>("mTID_P1") : .253 ;
27  mHC_P[SiStripDetId::TOB - 3][0] = conf.existsAs<double>("mTOB_P0") ? conf.getParameter<double>("mTOB_P0") : -1.427 ;
28  mHC_P[SiStripDetId::TOB - 3][1] = conf.existsAs<double>("mTOB_P1") ? conf.getParameter<double>("mTOB_P1") : .433 ;
29  mHC_P[SiStripDetId::TEC - 3][0] = conf.existsAs<double>("mTEC_P0") ? conf.getParameter<double>("mTEC_P0") : -1.885 ;
30  mHC_P[SiStripDetId::TEC - 3][1] = conf.existsAs<double>("mTEC_P1") ? conf.getParameter<double>("mTEC_P1") : .471 ;
31  }
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
StripCPE(edm::ParameterSet &conf, const MagneticField &, const TrackerGeometry &, const SiStripLorentzAngle &, const SiStripBackPlaneCorrection &, const SiStripConfObject &, const SiStripLatency &)
Definition: StripCPE.cc:11

Member Function Documentation

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

Definition at line 40 of file StripCPEfromTrackAngle.cc.

References constexpr, f, myMath::fast_expf(), objects.autophobj::float, 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 }
#define constexpr
#define unlikely(x)
double f[11][100]
#define N
Definition: blowfish.cc:9
return(e1-e2)*(e1-e2)+dp *dp
float fast_expf(float x)
void StripCPEfromTrackAngle::localParameters ( AClusters const &  clusters,
ALocalValues retValues,
const GeomDetUnit gd,
const LocalTrajectoryParameters ltp 
) const
overridevirtual

Reimplemented from StripClusterParameterEstimator.

Definition at line 54 of file StripCPEfromTrackAngle.cc.

References chargeCK, siStripClusterTools::chargePerCM(), corr, lumiContext::fill, StripCPE::getAlgoParam(), mps_fire::i, legacy, legacyStripErrorSquared(), create_public_lumi_plots::loc, m_algo, maxChgOneMIP, mergeCK, N, AlCaHLTBitMon_ParallelJobs::p, DynArray< T >::size(), digi_MixPreMix_cfi::strip, stripErrorSquared(), mitigatedMETSequence_cff::U, and LocalTrajectoryParameters::vector().

Referenced by localParameters().

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

Reimplemented from StripCPE.

Definition at line 99 of file StripCPEfromTrackAngle.cc.

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

99  {
100  auto const & p = par.p;
101  auto const & ltp = par.ltp;
102  auto loc = par.loc;
103  auto corr = par.corr;
104  auto afp = par.afullProjection;
105 
106  float uerr2=0;
107 
108  auto N = cluster.amplitudes().size();
109 
110  switch (m_algo) {
111  case Algo::chargeCK :
112  {
113  auto dQdx = siStripClusterTools::chargePerCM(cluster, ltp, p.invThickness);
114  uerr2 = dQdx > maxChgOneMIP ? legacyStripErrorSquared(N,afp) : stripErrorSquared( N, afp,loc );
115  }
116  break;
117  case Algo::legacy :
118  uerr2 = legacyStripErrorSquared(N,afp);
119  break;
120  case Algo::mergeCK :
121  uerr2 = cluster.isMerged() ? legacyStripErrorSquared(N,afp) : stripErrorSquared( N, afp,loc );
122  break;
123  }
124 
125  const float strip = cluster.barycenter() + corr;
126 
127  return std::make_pair( p.topology->localPosition(strip, ltp.vector()),
128  p.topology->localError(strip, uerr2, ltp.vector()) );
129 }
float legacyStripErrorSquared(const unsigned N, const float uProj) const
float chargePerCM(DetId detid, Iter a, Iter b)
JetCorrectorParameters corr
Definition: classes.h:5
#define N
Definition: blowfish.cc:9
float stripErrorSquared(const unsigned N, const float uProj, const SiStripDetId::SubDetector loc) const
StripClusterParameterEstimator::LocalValues StripCPEfromTrackAngle::localParameters ( const SiStripCluster cluster,
const GeomDetUnit det,
const LocalTrajectoryParameters ltp 
) const
overridevirtual

Reimplemented from StripClusterParameterEstimator.

Definition at line 134 of file StripCPEfromTrackAngle.cc.

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

134  {
135 
136  auto const & par = getAlgoParam(det,ltp);
137  return localParameters(cluster,par);
138 }
AlgoParam getAlgoParam(const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
Definition: StripCPE.h:55
void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const override
float StripCPEfromTrackAngle::stripErrorSquared ( const unsigned  N,
const float  uProj,
const SiStripDetId::SubDetector  loc 
) const

Definition at line 33 of file StripCPEfromTrackAngle.cc.

References myMath::fast_expf(), objects.autophobj::float, mHC_P, mLC_P, and x.

Referenced by localParameters().

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

Member Data Documentation

Algo StripCPEfromTrackAngle::m_algo
private

Definition at line 29 of file StripCPEfromTrackAngle.h.

Referenced by localParameters().

const float StripCPEfromTrackAngle::maxChgOneMIP
private

Definition at line 25 of file StripCPEfromTrackAngle.h.

Referenced by localParameters().

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

Definition at line 14 of file StripCPEfromTrackAngle.h.

Referenced by StripCPEfromTrackAngle(), and stripErrorSquared().

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

Definition at line 17 of file StripCPEfromTrackAngle.h.

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

Definition at line 18 of file StripCPEfromTrackAngle.h.

float StripCPEfromTrackAngle::mLC_P[3]
private

Definition at line 13 of file StripCPEfromTrackAngle.h.

Referenced by StripCPEfromTrackAngle(), and stripErrorSquared().

const bool StripCPEfromTrackAngle::useLegacyError
private

Definition at line 21 of file StripCPEfromTrackAngle.h.