CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
StripCPEfromTrackAngle Class Reference

#include <StripCPEfromTrackAngle.h>

Inheritance diagram for StripCPEfromTrackAngle:
StripCPE StripClusterParameterEstimator

Public Member Functions

float legacyStripErrorSquared (const unsigned N, const float uProj) const
 
StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &, const GeomDetUnit &, const LocalTrajectoryParameters &) const
 
 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
 
StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &cl, const GeomDetUnit &) const
 
 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 Attributes

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

- Public Types inherited from StripClusterParameterEstimator
typedef std::pair< LocalPoint,
LocalError
LocalValues
 
typedef std::vector< LocalValuesVLocalValues
 
- 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.

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 
)
inline

Definition at line 34 of file StripCPEfromTrackAngle.h.

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

41  : StripCPE(conf, mag, geom, lorentz, backPlaneCorrection, confObj, latency )
42  , useLegacyError(conf.existsAs<bool>("useLegacyError") ? conf.getParameter<bool>("useLegacyError") : true)
43  , maxChgOneMIP(conf.existsAs<float>("maxChgOneMIP") ? conf.getParameter<double>("maxChgOneMIP") : -6000.)
44  {
45  mLC_P[0] = conf.existsAs<double>("mLC_P0") ? conf.getParameter<double>("mLC_P0") : -.326;
46  mLC_P[1] = conf.existsAs<double>("mLC_P1") ? conf.getParameter<double>("mLC_P1") : .618;
47  mLC_P[2] = conf.existsAs<double>("mLC_P2") ? conf.getParameter<double>("mLC_P2") : .300;
48 
49  mHC_P[SiStripDetId::TIB - 3][0] = conf.existsAs<double>("mTIB_P0") ? conf.getParameter<double>("mTIB_P0") : -.742 ;
50  mHC_P[SiStripDetId::TIB - 3][1] = conf.existsAs<double>("mTIB_P1") ? conf.getParameter<double>("mTIB_P1") : .202 ;
51  mHC_P[SiStripDetId::TID - 3][0] = conf.existsAs<double>("mTID_P0") ? conf.getParameter<double>("mTID_P0") : -1.026 ;
52  mHC_P[SiStripDetId::TID - 3][1] = conf.existsAs<double>("mTID_P1") ? conf.getParameter<double>("mTID_P1") : .253 ;
53  mHC_P[SiStripDetId::TOB - 3][0] = conf.existsAs<double>("mTOB_P0") ? conf.getParameter<double>("mTOB_P0") : -1.427 ;
54  mHC_P[SiStripDetId::TOB - 3][1] = conf.existsAs<double>("mTOB_P1") ? conf.getParameter<double>("mTOB_P1") : .433 ;
55  mHC_P[SiStripDetId::TEC - 3][0] = conf.existsAs<double>("mTEC_P0") ? conf.getParameter<double>("mTEC_P0") : -1.885 ;
56  mHC_P[SiStripDetId::TEC - 3][1] = conf.existsAs<double>("mTEC_P1") ? conf.getParameter<double>("mTEC_P1") : .471 ;
57  }
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:185
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 16 of file StripCPEfromTrackAngle.cc.

References constexpr, f, and myMath::fast_expf().

Referenced by localParameters().

16  {
17  if( (float(N)-uProj) > 3.5f )
18  return float(N*N)/12.f;
19  else {
20  static constexpr float P1=-0.339;
21  static constexpr float P2=0.90;
22  static constexpr float P3=0.279;
23  const float uerr = P1*uProj*vdt::fast_expf(-uProj*P2)+P3;
24  return uerr*uerr;
25  }
26 }
#define constexpr
double f[11][100]
#define N
Definition: blowfish.cc:9
float fast_expf(float x)
StripClusterParameterEstimator::LocalValues StripCPEfromTrackAngle::localParameters ( const SiStripCluster cluster,
const GeomDetUnit det,
const LocalTrajectoryParameters ltp 
) const
virtual

Reimplemented from StripClusterParameterEstimator.

Definition at line 29 of file StripCPEfromTrackAngle.cc.

References funct::abs(), SiStripCluster::amplitudes(), StripCPE::Param::backplanecorrection, SiStripCluster::barycenter(), siStripClusterTools::chargePerCM(), StripCPE::Param::coveredStrips(), StripCPE::Param::drift, f, GeomDet::geographicalId(), SiStripCluster::isMerged(), legacyStripErrorSquared(), StripTopology::localError(), StripTopology::localPosition(), maxChgOneMIP, StripCPE::Param::maxLength, LocalTrajectoryParameters::momentum(), N, AlCaHLTBitMon_ParallelJobs::p, StripCPE::param(), LocalTrajectoryParameters::position(), stripErrorSquared(), SiStripDetId::subDetector(), StripCPE::Param::thickness, StripCPE::Param::topology, useLegacyError, and LocalTrajectoryParameters::vector().

29  {
30 
31  StripCPE::Param const & p = param(det);
32  SiStripDetId ssdid = SiStripDetId( det.geographicalId() );
33 
34  LocalVector track = ltp.momentum();
35  track *=
36  (track.z()<0) ? std::abs(p.thickness/track.z()) :
37  (track.z()>0) ? -std::abs(p.thickness/track.z()) :
38  p.maxLength/track.mag() ;
39 
40  const unsigned N = cluster.amplitudes().size();
41  const float fullProjection = p.coveredStrips( track+p.drift, ltp.position());
42  float uerr2;
43  if (useLegacyError) {
44  uerr2 = legacyStripErrorSquared(N,std::abs(fullProjection));
45  } else if (maxChgOneMIP < 0.0) {
46  uerr2 = cluster.isMerged() ? legacyStripErrorSquared(N,std::abs(fullProjection)) : stripErrorSquared( N, std::abs(fullProjection),ssdid.subDetector() );
47  } else {
48  float dQdx = siStripClusterTools::chargePerCM(ssdid, cluster, ltp);
49  uerr2 = dQdx > maxChgOneMIP ? legacyStripErrorSquared(N,std::abs(fullProjection)) : stripErrorSquared( N, std::abs(fullProjection),ssdid.subDetector() );
50  }
51  const float strip = cluster.barycenter() - 0.5f*(1.f-p.backplanecorrection) * fullProjection
52  + 0.5f*p.coveredStrips(track, ltp.position());
53 
54  return std::make_pair( p.topology->localPosition(strip, ltp.vector()),
55  p.topology->localError(strip, uerr2, ltp.vector()) );
56 }
float legacyStripErrorSquared(const unsigned N, const float uProj) const
bool isMerged() const
LocalPoint position() const
Local x and y position coordinates.
float thickness
Definition: StripCPE.h:45
StripTopology const * topology
Definition: StripCPE.h:43
float chargePerCM(DetId detid, Iter a, Iter b)
float backplanecorrection
Definition: StripCPE.h:47
AlgebraicVector5 vector() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:77
double f[11][100]
float barycenter() const
LocalVector momentum() const
Momentum vector in the local frame.
float coveredStrips(const LocalVector &, const LocalPoint &) const
Definition: StripCPE.cc:76
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
#define N
Definition: blowfish.cc:9
virtual LocalError localError(float strip, float stripErr2) const =0
float maxLength
Definition: StripCPE.h:45
Param const & param(const GeomDetUnit &det) const
Definition: StripCPE.h:51
LocalVector drift
Definition: StripCPE.h:44
virtual LocalPoint localPosition(float strip) const =0
const std::vector< uint8_t > & amplitudes() const
float stripErrorSquared(const unsigned N, const float uProj, const SiStripDetId::SubDetector loc) const
float StripCPEfromTrackAngle::stripErrorSquared ( const unsigned  N,
const float  uProj,
const SiStripDetId::SubDetector  loc 
) const

Definition at line 7 of file StripCPEfromTrackAngle.cc.

References edm::hlt::Exception, myMath::fast_expf(), mHC_P, mLC_P, SiStripDetId::UNKNOWN, and x.

Referenced by localParameters().

7  {
8  if( loc == SiStripDetId::UNKNOWN)
9  throw cms::Exception("StripCPEfromTrackAngle::stripErrorSquared", "Incompatible sub-detector.");
10 
11  auto fun = [&] (float x) -> float { return mLC_P[0]*x*vdt::fast_expf(-x*mLC_P[1])+mLC_P[2];};
12  auto uerr = (N <= 4) ? fun(uProj) : mHC_P[loc-3][0]+float(N)*mHC_P[loc-3][1];
13  return uerr*uerr;
14 }
#define N
Definition: blowfish.cc:9
float fast_expf(float x)
Definition: DDAxes.h:10

Member Data Documentation

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.

Referenced by localParameters().