#include <RecoEgamma/EgammaPhotonAlgos/interface/ConversionTrackEcalImpactPoint.h>
Marinelli - Univ. of Notre Dame
Definition at line 39 of file ConversionTrackEcalImpactPoint.h.
ConversionTrackEcalImpactPoint::ConversionTrackEcalImpactPoint | ( | const MagneticField * | field | ) |
Definition at line 27 of file ConversionTrackEcalImpactPoint.cc.
References alongMomentum, dir_, forwardPropagator_, and theMF_.
00027 : 00028 00029 theMF_(field) 00030 { 00031 00032 forwardPropagator_ = new PropagatorWithMaterial ( dir_ = alongMomentum, 0.000511, theMF_ ); 00033 00034 }
ConversionTrackEcalImpactPoint::~ConversionTrackEcalImpactPoint | ( | ) |
Definition at line 36 of file ConversionTrackEcalImpactPoint.cc.
References forwardPropagator_.
00036 { 00037 00038 delete forwardPropagator_ ; 00039 00040 }
static const BoundCylinder& ConversionTrackEcalImpactPoint::barrel | ( | ) | [inline, static, private] |
Definition at line 86 of file ConversionTrackEcalImpactPoint.h.
References check(), and theBarrel_.
Referenced by find().
00086 { check(); return *theBarrel_;}
static float ConversionTrackEcalImpactPoint::barrelHalfLength | ( | ) | [inline, static, private] |
static float ConversionTrackEcalImpactPoint::barrelRadius | ( | ) | [inline, static, private] |
Hard-wired numbers defining the surfaces on which the crystal front faces lie.
Definition at line 71 of file ConversionTrackEcalImpactPoint.h.
Referenced by initialize().
Definition at line 77 of file ConversionTrackEcalImpactPoint.h.
References initialize(), and theInit_.
Referenced by barrel(), negativeEtaEndcap(), and positiveEtaEndcap().
00077 {if (!theInit_) initialize();}
static float ConversionTrackEcalImpactPoint::endcapRadius | ( | ) | [inline, static, private] |
static float ConversionTrackEcalImpactPoint::endcapZ | ( | ) | [inline, static, private] |
std::vector< math::XYZPoint > ConversionTrackEcalImpactPoint::find | ( | const std::vector< reco::TransientTrack > & | tracks, | |
const edm::Handle< edm::View< reco::CaloCluster > > & | bcHandle | |||
) |
Definition at line 42 of file ConversionTrackEcalImpactPoint.cc.
References barrel(), dPhi(), PV3DBase< T, PVType, FrameType >::eta(), forwardPropagator_, TrajectoryStateOnSurface::globalPosition(), i, TrajectoryStateOnSurface::isValid(), matchingBC(), matchingBC_, negativeEtaEndcap(), positiveEtaEndcap(), PropagatorWithMaterial::propagate(), HLT_VtxMuL3::result, funct::sqrt(), and stateAtECAL_.
Referenced by ConvertedPhotonProducer::buildCollections(), and SoftConversionProducer::produce().
00042 { 00043 00044 00045 std::vector<math::XYZPoint> result; 00046 // 00047 matchingBC_.clear(); 00048 00049 std::vector<reco::CaloClusterPtr> matchingBC(2); 00050 00051 00052 // 00053 00054 00055 00056 int iTrk=0; 00057 for ( std::vector<reco::TransientTrack>::const_iterator iTk=tracks.begin(); iTk!=tracks.end(); ++iTk) { 00058 00059 math::XYZPoint ecalImpactPosition(0.,0.,0.); 00060 const TrajectoryStateOnSurface myTSOS=(*iTk).innermostMeasurementState(); 00061 if ( !( myTSOS.isValid() ) ) continue; 00062 00063 stateAtECAL_= forwardPropagator_->propagate( myTSOS, barrel() ); 00064 00065 00066 if (!stateAtECAL_.isValid() || ( stateAtECAL_.isValid() && fabs(stateAtECAL_.globalPosition().eta() ) >1.479 ) ) { 00067 00068 00069 if ( (*iTk).innermostMeasurementState().globalPosition().eta() > 0.) { 00070 stateAtECAL_= forwardPropagator_->propagate( myTSOS, positiveEtaEndcap()); 00071 00072 } else { 00073 00074 stateAtECAL_= forwardPropagator_->propagate( myTSOS, negativeEtaEndcap()); 00075 00076 } 00077 } 00078 00079 00080 if ( stateAtECAL_.isValid() ) ecalImpactPosition = stateAtECAL_.globalPosition(); 00081 00082 00083 result.push_back(ecalImpactPosition ); 00084 00085 00086 if ( stateAtECAL_.isValid() ) { 00087 int goodBC=0; 00088 float bcDistanceToTrack=9999; 00089 reco::CaloClusterPtr matchedBCItr; 00090 int ibc=0; 00091 goodBC=0; 00092 00093 for (unsigned i = 0; i < bcHandle->size(); ++i ) { 00094 float dEta= bcHandle->ptrAt(i)->position().eta() - ecalImpactPosition.eta() ; 00095 float dPhi= bcHandle->ptrAt(i)->position().phi() - ecalImpactPosition.phi() ; 00096 if ( sqrt(dEta*dEta + dPhi*dPhi) < bcDistanceToTrack ) { 00097 goodBC=ibc; 00098 bcDistanceToTrack=sqrt(dEta*dEta + dPhi*dPhi); 00099 } 00100 ibc++; 00101 00102 } 00103 00104 matchingBC[iTrk]=bcHandle->ptrAt(goodBC); 00105 } 00106 00107 00108 iTrk++; 00109 } 00110 00111 matchingBC_=matchingBC; 00112 00113 return result; 00114 00115 }
Definition at line 120 of file ConversionTrackEcalImpactPoint.cc.
References barrelHalfLength(), barrelRadius(), endcapRadius(), endcapZ(), geometryDiff::epsilon, rot, theBarrel_, theInit_, theNegativeEtaEndcap_, and thePositiveEtaEndcap_.
Referenced by check().
00120 { 00121 00122 const float epsilon = 0.001; 00123 Surface::RotationType rot; // unit rotation matrix 00124 00125 00126 theBarrel_ = new BoundCylinder( Surface::PositionType(0,0,0), rot, 00127 SimpleCylinderBounds( barrelRadius()-epsilon, 00128 barrelRadius()+epsilon, 00129 -barrelHalfLength(), 00130 barrelHalfLength())); 00131 theNegativeEtaEndcap_ = 00132 new BoundDisk( Surface::PositionType( 0, 0, -endcapZ()), rot, 00133 SimpleDiskBounds( 0, endcapRadius(), -epsilon, epsilon)); 00134 00135 thePositiveEtaEndcap_ = 00136 new BoundDisk( Surface::PositionType( 0, 0, endcapZ()), rot, 00137 SimpleDiskBounds( 0, endcapRadius(), -epsilon, epsilon)); 00138 00139 theInit_ = true; 00140 00141 00142 }
std::vector<reco::CaloClusterPtr> ConversionTrackEcalImpactPoint::matchingBC | ( | ) | const [inline] |
Definition at line 51 of file ConversionTrackEcalImpactPoint.h.
References matchingBC_.
Referenced by ConvertedPhotonProducer::buildCollections(), and find().
00051 {return matchingBC_;}
static const BoundDisk& ConversionTrackEcalImpactPoint::negativeEtaEndcap | ( | ) | [inline, static, private] |
Definition at line 87 of file ConversionTrackEcalImpactPoint.h.
References check(), and theNegativeEtaEndcap_.
Referenced by find().
00087 { check(); return *theNegativeEtaEndcap_;}
static const BoundDisk& ConversionTrackEcalImpactPoint::positiveEtaEndcap | ( | ) | [inline, static, private] |
Definition at line 88 of file ConversionTrackEcalImpactPoint.h.
References check(), and thePositiveEtaEndcap_.
Referenced by find().
00088 { check(); return *thePositiveEtaEndcap_;}
void ConversionTrackEcalImpactPoint::setMagneticField | ( | const MagneticField * | mf | ) | [inline] |
Definition at line 52 of file ConversionTrackEcalImpactPoint.h.
References theMF_.
00052 { theMF_=mf;}
Definition at line 65 of file ConversionTrackEcalImpactPoint.h.
Referenced by ConversionTrackEcalImpactPoint().
PropagatorWithMaterial* ConversionTrackEcalImpactPoint::forwardPropagator_ [mutable, private] |
Definition at line 64 of file ConversionTrackEcalImpactPoint.h.
Referenced by ConversionTrackEcalImpactPoint(), find(), and ~ConversionTrackEcalImpactPoint().
std::vector<reco::CaloClusterPtr> ConversionTrackEcalImpactPoint::matchingBC_ [private] |
Definition at line 66 of file ConversionTrackEcalImpactPoint.h.
Referenced by find(), and matchingBC().
ReferenceCountingPointer< BoundCylinder > ConversionTrackEcalImpactPoint::theBarrel_ = 0 [static, private] |
Definition at line 82 of file ConversionTrackEcalImpactPoint.h.
Referenced by barrel(), and initialize().
bool ConversionTrackEcalImpactPoint::theInit_ = false [static, private] |
Definition at line 89 of file ConversionTrackEcalImpactPoint.h.
Referenced by check(), and initialize().
const MagneticField* ConversionTrackEcalImpactPoint::theMF_ [private] |
Definition at line 60 of file ConversionTrackEcalImpactPoint.h.
Referenced by ConversionTrackEcalImpactPoint(), and setMagneticField().
ReferenceCountingPointer< BoundDisk > ConversionTrackEcalImpactPoint::theNegativeEtaEndcap_ = 0 [static, private] |
Definition at line 83 of file ConversionTrackEcalImpactPoint.h.
Referenced by initialize(), and negativeEtaEndcap().
ReferenceCountingPointer< BoundDisk > ConversionTrackEcalImpactPoint::thePositiveEtaEndcap_ = 0 [static, private] |
Definition at line 84 of file ConversionTrackEcalImpactPoint.h.
Referenced by initialize(), and positiveEtaEndcap().