#include <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 12 of file VZeroFinder.cc.
References edm::EventSetup::get(), and edm::ParameterSet::getParameter().
{ // Get track-pair level cuts maxDca = pset.getParameter<double>("maxDca"); // Get mother cuts minCrossingRadius = pset.getParameter<double>("minCrossingRadius"); maxCrossingRadius = pset.getParameter<double>("maxCrossingRadius"); maxImpactMother = pset.getParameter<double>("maxImpactMother"); // Get magnetic field edm::ESHandle<MagneticField> magField; es.get<IdealMagneticFieldRecord>().get(magField); theMagField = magField.product(); }
VZeroFinder::~VZeroFinder | ( | ) |
Definition at line 30 of file VZeroFinder.cc.
{ }
bool VZeroFinder::checkTrackPair | ( | const reco::Track & | posTrack, |
const reco::Track & | negTrack, | ||
const reco::VertexCollection * | vertices, | ||
reco::VZeroData & | data | ||
) |
Definition at line 60 of file VZeroFinder.cc.
References alpha, 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(), mag(), PV3DBase< T, PVType, FrameType >::mag2(), maxCrossingRadius, maxDca, maxImpactMother, minCrossingRadius, reco::VZeroData::momenta, GlobalTrajectoryParameters::momentum(), AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::perp(), alignCSCRings::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().
{ /* LogTrace("MinBiasTracking") << " [VZeroFinder] tracks" << " +" << posTrack.algoName() << " " << posTrack.d0() << " -" << negTrack.algoName() << " " << negTrack.d0(); */ // Get trajectories GlobalTrajectoryParameters posGtp = getGlobalTrajectoryParameters(posTrack); GlobalTrajectoryParameters negGtp = getGlobalTrajectoryParameters(negTrack); // Two track minimum distance TwoTrackMinimumDistance theMinimum(TwoTrackMinimumDistance::SlowMode); // Closest approach theMinimum.calculate(posGtp,negGtp); // Closest points // std::pair<GlobalPoint, GlobalPoint> points = theMinimum.points(); // Momenta at closest point std::pair<GlobalVector,GlobalVector> momenta; momenta.first = rotate(posGtp.momentum(), theMinimum.firstAngle() ); momenta.second = rotate(negGtp.momentum(), theMinimum.secondAngle()); // dcaR float dca = theMinimum.distance(); GlobalPoint crossing = theMinimum.crossingPoint(); /* LogTrace("MinBiasTracking") << fixed << setprecision(2) << " [VZeroFinder] dca = "<<dca<<" (<"<<maxDca<<")"; LogTrace("MinBiasTracking") << fixed << setprecision(2) << " [VZeroFinder] crossR = "<<crossing.perp()<<" (" <<minCrossingRadius<<"< <"<<maxCrossingRadius<<")"; */ if(dca < maxDca && crossing.perp() > minCrossingRadius && crossing.perp() < maxCrossingRadius) { // Momentum of the mother GlobalVector momentum = momenta.first + momenta.second; float impact = -1.; if(vertices->size() > 0) { // Impact parameter of the mother wrt vertices, choose smallest for(reco::VertexCollection::const_iterator vertex = vertices->begin(); vertex!= vertices->end(); vertex++) { GlobalVector r(crossing.x(), crossing.y(), crossing.z() - vertex->position().z()); GlobalVector p(momentum.x(),momentum.y(),momentum.z()); GlobalVector b = r - (r*p)*p / p.mag2(); float im = b.mag(); if(im < impact || vertex == vertices->begin()) impact = im; } } else { // Impact parameter of the mother in the plane GlobalVector r_(crossing.x(),crossing.y(),0); GlobalVector p_(momentum.x(),momentum.y(),0); GlobalVector b_ = r_ - (r_*p_)*p_ / p_.mag2(); impact = b_.mag(); } /* LogTrace("MinBiasTracking") << fixed << setprecision(2) << " [VZeroFinder] impact = "<<impact<<" (<"<<maxImpactMother<<")"; */ if(impact < maxImpactMother) { // Armenteros variables float pt = (momenta.first.cross(momenta.second)).mag()/momentum.mag(); float alpha = (momenta.first.mag2() - momenta.second.mag2())/momentum.mag2(); // Fill data data.dca = dca; data.crossingPoint = crossing; data.impactMother = impact; data.armenterosPt = pt; data.armenterosAlpha = alpha; data.momenta = momenta; /* LogTrace("MinBiasTracking") << fixed << setprecision(2) << " [VZeroFinder] success {alpha = "<<alpha<<", pt = "<<pt<<"}"; */ return true; } } return false; }
GlobalTrajectoryParameters VZeroFinder::getGlobalTrajectoryParameters | ( | const reco::Track & | track | ) | [private] |
Definition at line 36 of file VZeroFinder.cc.
References reco::TrackBase::charge(), reco::TrackBase::momentum(), position, and reco::TrackBase::vertex().
Referenced by checkTrackPair().
{ GlobalPoint position(track.vertex().x(), track.vertex().y(), track.vertex().z()); GlobalVector momentum(track.momentum().x(), track.momentum().y(), track.momentum().z()); GlobalTrajectoryParameters gtp(position,momentum, track.charge(),theMagField); return gtp; }
GlobalVector VZeroFinder::rotate | ( | const GlobalVector & | p, |
double | a | ||
) |
Definition at line 53 of file VZeroFinder.cc.
References funct::cos(), PV3DBase< T, PVType, FrameType >::perp(), funct::sin(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by checkTrackPair().
float VZeroFinder::maxCrossingRadius [private] |
Definition at line 38 of file VZeroFinder.h.
Referenced by checkTrackPair().
float VZeroFinder::maxDca [private] |
Definition at line 38 of file VZeroFinder.h.
Referenced by checkTrackPair().
float VZeroFinder::maxImpactMother [private] |
Definition at line 38 of file VZeroFinder.h.
Referenced by checkTrackPair().
float VZeroFinder::minCrossingRadius [private] |
Definition at line 38 of file VZeroFinder.h.
Referenced by checkTrackPair().
const MagneticField* VZeroFinder::theMagField [private] |
Definition at line 43 of file VZeroFinder.h.