CMS 3D CMS Logo

List of all members | Public Member Functions | Static Private Member Functions | Private Attributes | Static Private Attributes
ConversionTrackEcalImpactPoint Class Reference

#include <ConversionTrackEcalImpactPoint.h>

Public Member Functions

 ConversionTrackEcalImpactPoint (const MagneticField *field)
 
std::vector< math::XYZPointFfind (const std::vector< reco::TransientTrack > &tracks, const edm::Handle< edm::View< reco::CaloCluster > > &bcHandle)
 
std::vector< reco::CaloClusterPtrmatchingBC () const
 
void setMagneticField (const MagneticField *mf)
 
 ~ConversionTrackEcalImpactPoint ()
 

Static Private Member Functions

static const BoundCylinderbarrel ()
 
static const BoundDisknegativeEtaEndcap ()
 
static const BoundDiskpositiveEtaEndcap ()
 

Private Attributes

PropagationDirection dir_
 
PropagatorWithMaterialforwardPropagator_
 
std::vector< reco::CaloClusterPtrmatchingBC_
 
TrajectoryStateOnSurface stateAtECAL_
 
const MagneticFieldtheMF_
 

Static Private Attributes

static const ReferenceCountingPointer< BoundCylindertheBarrel_ = initBarrel()
 
static const ReferenceCountingPointer< BoundDisktheNegativeEtaEndcap_ = initNegative()
 
static const ReferenceCountingPointer< BoundDiskthePositiveEtaEndcap_ = initPositive()
 

Detailed Description

Author
N. Marinelli - Univ. of Notre Dame
Version

Definition at line 37 of file ConversionTrackEcalImpactPoint.h.

Constructor & Destructor Documentation

◆ ConversionTrackEcalImpactPoint()

ConversionTrackEcalImpactPoint::ConversionTrackEcalImpactPoint ( const MagneticField field)

◆ ~ConversionTrackEcalImpactPoint()

ConversionTrackEcalImpactPoint::~ConversionTrackEcalImpactPoint ( )

Definition at line 54 of file ConversionTrackEcalImpactPoint.cc.

References forwardPropagator_.

54 { delete forwardPropagator_; }

Member Function Documentation

◆ barrel()

static const BoundCylinder& ConversionTrackEcalImpactPoint::barrel ( )
inlinestaticprivate

Definition at line 62 of file ConversionTrackEcalImpactPoint.h.

References theBarrel_.

Referenced by find().

62 { return *theBarrel_; }
static const ReferenceCountingPointer< BoundCylinder > theBarrel_

◆ find()

std::vector< math::XYZPointF > ConversionTrackEcalImpactPoint::find ( const std::vector< reco::TransientTrack > &  tracks,
const edm::Handle< edm::View< reco::CaloCluster > > &  bcHandle 
)

Definition at line 56 of file ConversionTrackEcalImpactPoint.cc.

References barrel(), HLT_2024v13_cff::dEta, HLT_2024v13_cff::dPhi, PV3DBase< T, PVType, FrameType >::eta(), forwardPropagator_, TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), matchingBC(), matchingBC_, negativeEtaEndcap(), positiveEtaEndcap(), PropagatorWithMaterial::propagate(), mps_fire::result, mathSSE::sqrt(), stateAtECAL_, HLT_2024v13_cff::track, and DiMuonV_cfg::tracks.

Referenced by ConvertedPhotonProducer::buildCollections().

57  {
58  std::vector<math::XYZPointF> result;
59  //
60  matchingBC_.clear();
61 
62  std::vector<reco::CaloClusterPtr> matchingBC(2);
63 
64  //
65 
66  int iTrk = 0;
67  for (auto const& track : tracks) {
68  math::XYZPointF ecalImpactPosition(0., 0., 0.);
69  const TrajectoryStateOnSurface myTSOS = track.innermostMeasurementState();
70  if (!(myTSOS.isValid()))
71  continue;
72 
74 
75  if (!stateAtECAL_.isValid() || (stateAtECAL_.isValid() && fabs(stateAtECAL_.globalPosition().eta()) > 1.479)) {
76  if (track.innermostMeasurementState().globalPosition().eta() > 0.) {
78 
79  } else {
81  }
82  }
83 
84  if (stateAtECAL_.isValid())
85  ecalImpactPosition = stateAtECAL_.globalPosition();
86 
87  result.push_back(ecalImpactPosition);
88 
89  if (stateAtECAL_.isValid()) {
90  int goodBC = 0;
91  float bcDistanceToTrack = 9999;
92  reco::CaloClusterPtr matchedBCItr;
93  int ibc = 0;
94  goodBC = 0;
95 
96  for (auto const& bc : *bcHandle) {
97  float dEta = bc.position().eta() - ecalImpactPosition.eta();
98  float dPhi = bc.position().phi() - ecalImpactPosition.phi();
99  if (sqrt(dEta * dEta + dPhi * dPhi) < bcDistanceToTrack) {
100  goodBC = ibc;
101  bcDistanceToTrack = sqrt(dEta * dEta + dPhi * dPhi);
102  }
103  ibc++;
104  }
105 
106  matchingBC[iTrk] = bcHandle->ptrAt(goodBC);
107  }
108 
109  iTrk++;
110  }
111 
113 
114  return result;
115 }
T eta() const
Definition: PV3DBase.h:73
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
static const BoundDisk & positiveEtaEndcap()
std::vector< reco::CaloClusterPtr > matchingBC() const
std::vector< reco::CaloClusterPtr > matchingBC_
GlobalPoint globalPosition() const
static const BoundDisk & negativeEtaEndcap()
T sqrt(T t)
Definition: SSEVec.h:23
static const BoundCylinder & barrel()
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:50

◆ matchingBC()

std::vector<reco::CaloClusterPtr> ConversionTrackEcalImpactPoint::matchingBC ( ) const
inline

Definition at line 46 of file ConversionTrackEcalImpactPoint.h.

References matchingBC_.

Referenced by ConvertedPhotonProducer::buildCollections(), and find().

46 { return matchingBC_; }
std::vector< reco::CaloClusterPtr > matchingBC_

◆ negativeEtaEndcap()

static const BoundDisk& ConversionTrackEcalImpactPoint::negativeEtaEndcap ( )
inlinestaticprivate

Definition at line 63 of file ConversionTrackEcalImpactPoint.h.

References theNegativeEtaEndcap_.

Referenced by find().

63 { return *theNegativeEtaEndcap_; }
static const ReferenceCountingPointer< BoundDisk > theNegativeEtaEndcap_

◆ positiveEtaEndcap()

static const BoundDisk& ConversionTrackEcalImpactPoint::positiveEtaEndcap ( )
inlinestaticprivate

Definition at line 64 of file ConversionTrackEcalImpactPoint.h.

References thePositiveEtaEndcap_.

Referenced by find().

64 { return *thePositiveEtaEndcap_; }
static const ReferenceCountingPointer< BoundDisk > thePositiveEtaEndcap_

◆ setMagneticField()

void ConversionTrackEcalImpactPoint::setMagneticField ( const MagneticField mf)
inline

Definition at line 47 of file ConversionTrackEcalImpactPoint.h.

References theMF_.

47 { theMF_ = mf; }

Member Data Documentation

◆ dir_

PropagationDirection ConversionTrackEcalImpactPoint::dir_
private

Definition at line 55 of file ConversionTrackEcalImpactPoint.h.

Referenced by ConversionTrackEcalImpactPoint().

◆ forwardPropagator_

PropagatorWithMaterial* ConversionTrackEcalImpactPoint::forwardPropagator_
private

◆ matchingBC_

std::vector<reco::CaloClusterPtr> ConversionTrackEcalImpactPoint::matchingBC_
private

Definition at line 56 of file ConversionTrackEcalImpactPoint.h.

Referenced by find(), and matchingBC().

◆ stateAtECAL_

TrajectoryStateOnSurface ConversionTrackEcalImpactPoint::stateAtECAL_
private

Definition at line 52 of file ConversionTrackEcalImpactPoint.h.

Referenced by find().

◆ theBarrel_

const ReferenceCountingPointer< BoundCylinder > ConversionTrackEcalImpactPoint::theBarrel_ = initBarrel()
staticprivate

Definition at line 58 of file ConversionTrackEcalImpactPoint.h.

Referenced by barrel().

◆ theMF_

const MagneticField* ConversionTrackEcalImpactPoint::theMF_
private

◆ theNegativeEtaEndcap_

const ReferenceCountingPointer< BoundDisk > ConversionTrackEcalImpactPoint::theNegativeEtaEndcap_ = initNegative()
staticprivate

Definition at line 59 of file ConversionTrackEcalImpactPoint.h.

Referenced by negativeEtaEndcap().

◆ thePositiveEtaEndcap_

const ReferenceCountingPointer< BoundDisk > ConversionTrackEcalImpactPoint::thePositiveEtaEndcap_ = initPositive()
staticprivate

Definition at line 60 of file ConversionTrackEcalImpactPoint.h.

Referenced by positiveEtaEndcap().