#include <RecoVZero/VZeroFinding/interface/VZeroFinder.h>
Public Member Functions | |
bool | checkTrackPair (const reco::Track &posTrack, const reco::Track &negTrack, const reco::VertexCollection *vertices, reco::VZeroData &data) |
GlobalVector | rotate (const GlobalVector &p, double a) |
VZeroFinder (const edm::EventSetup &es, const edm::ParameterSet &pset) | |
~VZeroFinder () | |
Private Member Functions | |
GlobalTrajectoryParameters | getGlobalTrajectoryParameters (const reco::Track &track) |
Private Attributes | |
float | maxCrossingRadius |
float | maxDca |
float | maxImpactMother |
float | minCrossingRadius |
const MagneticField * | theMagField |
Definition at line 22 of file VZeroFinder.h.
VZeroFinder::VZeroFinder | ( | const edm::EventSetup & | es, | |
const edm::ParameterSet & | pset | |||
) |
Definition at line 14 of file VZeroFinder.cc.
References edm::EventSetup::get(), and edm::ParameterSet::getParameter().
00016 { 00017 // Get track-pair level cuts 00018 maxDca = pset.getParameter<double>("maxDca"); 00019 00020 // Get mother cuts 00021 minCrossingRadius = pset.getParameter<double>("minCrossingRadius"); 00022 maxCrossingRadius = pset.getParameter<double>("maxCrossingRadius"); 00023 maxImpactMother = pset.getParameter<double>("maxImpactMother"); 00024 00025 // Get magnetic field 00026 edm::ESHandle<MagneticField> magField; 00027 es.get<IdealMagneticFieldRecord>().get(magField); 00028 theMagField = magField.product(); 00029 }
VZeroFinder::~VZeroFinder | ( | ) |
bool VZeroFinder::checkTrackPair | ( | const reco::Track & | posTrack, | |
const reco::Track & | negTrack, | |||
const reco::VertexCollection * | vertices, | |||
reco::VZeroData & | data | |||
) |
Definition at line 62 of file VZeroFinder.cc.
References reco::VZeroData::armenterosAlpha, reco::VZeroData::armenterosPt, b, begin, TwoTrackMinimumDistance::calculate(), TwoTrackMinimumDistance::crossingPoint(), reco::VZeroData::crossingPoint, reco::VZeroData::dca, TwoTrackMinimumDistance::distance(), TwoTrackMinimumDistance::firstAngle(), getGlobalTrajectoryParameters(), reco::VZeroData::impactMother, PV3DBase< T, PVType, FrameType >::mag(), muonGeometry::mag(), PV3DBase< T, PVType, FrameType >::mag2(), maxCrossingRadius, maxDca, maxImpactMother, minCrossingRadius, reco::VZeroData::momenta, GlobalTrajectoryParameters::momentum(), p, PV3DBase< T, PVType, FrameType >::perp(), TwoTrackMinimumDistance::points(), r, rotate(), TwoTrackMinimumDistance::secondAngle(), TwoTrackMinimumDistance::SlowMode, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by VZeroProducer::produce().
00066 { 00067 /* 00068 LogTrace("MinBiasTracking") << " [VZeroFinder] tracks" 00069 << " +" << posTrack.algoName() << " " << posTrack.d0() 00070 << " -" << negTrack.algoName() << " " << negTrack.d0(); 00071 */ 00072 00073 // Get trajectories 00074 GlobalTrajectoryParameters posGtp = getGlobalTrajectoryParameters(posTrack); 00075 GlobalTrajectoryParameters negGtp = getGlobalTrajectoryParameters(negTrack); 00076 00077 // Two track minimum distance 00078 TwoTrackMinimumDistance theMinimum(TwoTrackMinimumDistance::SlowMode); 00079 00080 // Closest approach 00081 theMinimum.calculate(posGtp,negGtp); 00082 00083 // Closest points 00084 pair<GlobalPoint, GlobalPoint> points = theMinimum.points(); 00085 00086 // Momenta at closest point 00087 pair<GlobalVector,GlobalVector> momenta; 00088 momenta.first = rotate(posGtp.momentum(), theMinimum.firstAngle() ); 00089 momenta.second = rotate(negGtp.momentum(), theMinimum.secondAngle()); 00090 00091 // dcaR 00092 float dca = theMinimum.distance(); 00093 GlobalPoint crossing = theMinimum.crossingPoint(); 00094 00095 /* 00096 LogTrace("MinBiasTracking") << fixed << setprecision(2) 00097 << " [VZeroFinder] dca = "<<dca<<" (<"<<maxDca<<")"; 00098 LogTrace("MinBiasTracking") << fixed << setprecision(2) 00099 << " [VZeroFinder] crossR = "<<crossing.perp()<<" (" 00100 <<minCrossingRadius<<"< <"<<maxCrossingRadius<<")"; 00101 */ 00102 00103 if(dca < maxDca && 00104 crossing.perp() > minCrossingRadius && 00105 crossing.perp() < maxCrossingRadius) 00106 { 00107 // Momentum of the mother 00108 GlobalVector momentum = momenta.first + momenta.second; 00109 float impact = -1.; 00110 00111 if(vertices->size() > 0) 00112 { 00113 // Impact parameter of the mother wrt vertices, choose smallest 00114 for(reco::VertexCollection::const_iterator 00115 vertex = vertices->begin(); vertex!= vertices->end(); vertex++) 00116 { 00117 GlobalVector r(crossing.x(), 00118 crossing.y(), 00119 crossing.z() - vertex->position().z()); 00120 GlobalVector p(momentum.x(),momentum.y(),momentum.z()); 00121 00122 GlobalVector b = r - (r*p)*p / p.mag2(); 00123 00124 float im = b.mag(); 00125 if(im < impact || vertex == vertices->begin()) 00126 impact = im; 00127 } 00128 } 00129 else 00130 { 00131 // Impact parameter of the mother in the plane 00132 GlobalVector r_(crossing.x(),crossing.y(),0); 00133 GlobalVector p_(momentum.x(),momentum.y(),0); 00134 00135 GlobalVector b_ = r_ - (r_*p_)*p_ / p_.mag2(); 00136 impact = b_.mag(); 00137 } 00138 00139 /* 00140 LogTrace("MinBiasTracking") << fixed << setprecision(2) 00141 << " [VZeroFinder] impact = "<<impact<<" (<"<<maxImpactMother<<")"; 00142 */ 00143 00144 if(impact < maxImpactMother) 00145 { 00146 // Armenteros variables 00147 float pt = 00148 (momenta.first.cross(momenta.second)).mag()/momentum.mag(); 00149 float alpha = 00150 (momenta.first.mag2() - momenta.second.mag2())/momentum.mag2(); 00151 00152 // Fill data 00153 data.dca = dca; 00154 data.crossingPoint = crossing; 00155 data.impactMother = impact; 00156 data.armenterosPt = pt; 00157 data.armenterosAlpha = alpha; 00158 data.momenta = momenta; 00159 00160 /* 00161 LogTrace("MinBiasTracking") << fixed << setprecision(2) 00162 << " [VZeroFinder] success {alpha = "<<alpha<<", pt = "<<pt<<"}"; 00163 */ 00164 00165 return true; 00166 } 00167 } 00168 00169 return false; 00170 }
GlobalTrajectoryParameters VZeroFinder::getGlobalTrajectoryParameters | ( | const reco::Track & | track | ) | [private] |
Definition at line 38 of file VZeroFinder.cc.
References reco::TrackBase::charge(), reco::TrackBase::momentum(), theMagField, and reco::TrackBase::vertex().
Referenced by checkTrackPair().
00039 { 00040 GlobalPoint position(track.vertex().x(), 00041 track.vertex().y(), 00042 track.vertex().z()); 00043 00044 GlobalVector momentum(track.momentum().x(), 00045 track.momentum().y(), 00046 track.momentum().z()); 00047 00048 GlobalTrajectoryParameters gtp(position,momentum, 00049 track.charge(),theMagField); 00050 00051 return gtp; 00052 }
GlobalVector VZeroFinder::rotate | ( | const GlobalVector & | p, | |
double | a | |||
) |
Definition at line 55 of file VZeroFinder.cc.
References funct::cos(), PV3DBase< T, PVType, FrameType >::perp(), funct::sin(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by checkTrackPair().
00056 { 00057 double pt = p.perp(); 00058 return GlobalVector( -pt*cos(a), -pt*sin(a), p.z()); 00059 }
float VZeroFinder::maxCrossingRadius [private] |
float VZeroFinder::maxDca [private] |
float VZeroFinder::maxImpactMother [private] |
float VZeroFinder::minCrossingRadius [private] |
const MagneticField* VZeroFinder::theMagField [private] |