#include <TracksClusteringFromDisplacedSeed.h>
Classes | |
struct | Cluster |
Public Member Functions | |
std::vector< Cluster > | clusters (const reco::Vertex &pv, const std::vector< reco::TransientTrack > &selectedTracks) |
TracksClusteringFromDisplacedSeed (const edm::ParameterSet ¶ms) | |
Private Member Functions | |
std::pair< std::vector < reco::TransientTrack > , GlobalPoint > | nearTracks (const reco::TransientTrack &seed, const std::vector< reco::TransientTrack > &tracks, const reco::Vertex &primaryVertex) const |
bool | trackFilter (const reco::TrackRef &track) const |
Private Attributes | |
double | clusterMaxDistance |
double | clusterMaxSignificance |
double | clusterMinAngleCosine |
double | clusterScale |
double | min3DIPSignificance |
double | min3DIPValue |
Definition at line 25 of file TracksClusteringFromDisplacedSeed.h.
TracksClusteringFromDisplacedSeed::TracksClusteringFromDisplacedSeed | ( | const edm::ParameterSet & | params | ) |
Definition at line 5 of file TracksClusteringFromDisplacedSeed.cc.
: // maxNTracks(params.getParameter<unsigned int>("maxNTracks")), min3DIPSignificance(params.getParameter<double>("seedMin3DIPSignificance")), min3DIPValue(params.getParameter<double>("seedMin3DIPValue")), clusterMaxDistance(params.getParameter<double>("clusterMaxDistance")), clusterMaxSignificance(params.getParameter<double>("clusterMaxSignificance")), //3 clusterScale(params.getParameter<double>("clusterScale")),//10. clusterMinAngleCosine(params.getParameter<double>("clusterMinAngleCosine")) //0.0 { }
std::vector< TracksClusteringFromDisplacedSeed::Cluster > TracksClusteringFromDisplacedSeed::clusters | ( | const reco::Vertex & | pv, |
const std::vector< reco::TransientTrack > & | selectedTracks | ||
) |
Definition at line 88 of file TracksClusteringFromDisplacedSeed.cc.
References IPTools::absoluteImpactParameter3D(), gather_cfg::cout, i, min3DIPSignificance, min3DIPValue, nearTracks(), dt_dqm_sourceclient_common_cff::reco, alignCSCRings::s, TracksClusteringFromDisplacedSeed::Cluster::seedingTrack, TracksClusteringFromDisplacedSeed::Cluster::seedPoint, and TracksClusteringFromDisplacedSeed::Cluster::tracks.
{ using namespace reco; std::vector<TransientTrack> seeds; for(std::vector<TransientTrack>::const_iterator it = selectedTracks.begin(); it != selectedTracks.end(); it++){ std::pair<bool,Measurement1D> ip = IPTools::absoluteImpactParameter3D(*it,pv); if(ip.first && ip.second.value() >= min3DIPValue && ip.second.significance() >= min3DIPSignificance) { #ifdef VTXDEBUG std::cout << "new seed " << it-selectedTracks.begin() << " ref " << it->trackBaseRef().key() << " " << ip.second.value() << " " << ip.second.significance() << " " << it->track().hitPattern().trackerLayersWithMeasurement() << " " << it->track().pt() << " " << it->track().eta() << std::endl; #endif seeds.push_back(*it); } } std::vector< Cluster > clusters; int i = 0; for(std::vector<TransientTrack>::const_iterator s = seeds.begin(); s != seeds.end(); ++s, ++i) { #ifdef VTXDEBUG std::cout << "Seed N. "<<i << std::endl; #endif // VTXDEBUG std::pair<std::vector<reco::TransientTrack>,GlobalPoint> ntracks = nearTracks(*s,selectedTracks,pv); // std::cout << ntracks.first.size() << " " << ntracks.first.size() << std::endl; // if(ntracks.first.size() == 0 || ntracks.first.size() > maxNTracks ) continue; ntracks.first.push_back(*s); Cluster aCl; aCl.seedingTrack = *s; aCl.seedPoint = ntracks.second; aCl.tracks = ntracks.first; clusters.push_back(aCl); } return clusters; }
std::pair< std::vector< reco::TransientTrack >, GlobalPoint > TracksClusteringFromDisplacedSeed::nearTracks | ( | const reco::TransientTrack & | seed, |
const std::vector< reco::TransientTrack > & | tracks, | ||
const reco::Vertex & | primaryVertex | ||
) | const [private] |
Definition at line 18 of file TracksClusteringFromDisplacedSeed.cc.
References IPTools::absoluteImpactParameter3D(), TwoTrackMinimumDistance::calculate(), TrajectoryStateOnSurface::cartesianError(), clusterMaxDistance, clusterMaxSignificance, clusterMinAngleCosine, clusterScale, gather_cfg::cout, CommonMethods::cp(), TwoTrackMinimumDistance::crossingPoint(), VertexDistance3D::distance(), TwoTrackMinimumDistance::distance(), TrajectoryStateOnSurface::globalDirection(), reco::TransientTrack::impactPointState(), m, mag(), TwoTrackMinimumDistance::points(), CartesianTrajectoryError::position(), reco::Vertex::position(), query::result, Measurement1D::significance(), mathSSE::sqrt(), groupFilesInBlocks::tt, csvLumiCalc::unit, Vector3DBase< T, FrameTag >::unit(), w(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by clusters().
{ VertexDistance3D distanceComputer; GlobalPoint pv(primaryVertex.position().x(),primaryVertex.position().y(),primaryVertex.position().z()); std::vector<reco::TransientTrack> result; TwoTrackMinimumDistance dist; GlobalPoint seedingPoint; float sumWeights=0; std::pair<bool,Measurement1D> ipSeed = IPTools::absoluteImpactParameter3D(seed,primaryVertex); float pvDistance = ipSeed.second.value(); // float densityFactor = 2./sqrt(20.*tracks.size()); // assuming all tracks being in 2 narrow jets of cone 0.3 float densityFactor = 2./sqrt(20.*80); // assuming 80 tracks being in 2 narrow jets of cone 0.3 for(std::vector<reco::TransientTrack>::const_iterator tt = tracks.begin();tt!=tracks.end(); ++tt ) { if(*tt==seed) continue; std::pair<bool,Measurement1D> ip = IPTools::absoluteImpactParameter3D(*tt,primaryVertex); if(dist.calculate(tt->impactPointState(),seed.impactPointState())) { GlobalPoint ttPoint = dist.points().first; GlobalError ttPointErr = tt->impactPointState().cartesianError().position(); GlobalPoint seedPosition = dist.points().second; GlobalError seedPositionErr = seed.impactPointState().cartesianError().position(); Measurement1D m = distanceComputer.distance(VertexState(seedPosition,seedPositionErr), VertexState(ttPoint, ttPointErr)); GlobalPoint cp(dist.crossingPoint()); float distanceFromPV = (dist.points().second-pv).mag(); float distance = dist.distance(); GlobalVector trackDir2D(tt->impactPointState().globalDirection().x(),tt->impactPointState().globalDirection().y(),0.); GlobalVector seedDir2D(seed.impactPointState().globalDirection().x(),seed.impactPointState().globalDirection().y(),0.); //SK:UNUSED// float dotprodTrackSeed2D = trackDir2D.unit().dot(seedDir2D.unit()); float dotprodTrack = (dist.points().first-pv).unit().dot(tt->impactPointState().globalDirection().unit()); float dotprodSeed = (dist.points().second-pv).unit().dot(seed.impactPointState().globalDirection().unit()); float w = distanceFromPV*distanceFromPV/(pvDistance*distance); bool selected = (m.significance() < clusterMaxSignificance && dotprodSeed > clusterMinAngleCosine && //Angles between PV-PCAonSeed vectors and seed directions dotprodTrack > clusterMinAngleCosine && //Angles between PV-PCAonTrack vectors and track directions // dotprodTrackSeed2D > clusterMinAngleCosine && //Angle between track and seed // distance*clusterScale*tracks.size() < (distanceFromPV+pvDistance)*(distanceFromPV+pvDistance)/pvDistance && // cut scaling with track density distance*clusterScale < densityFactor*distanceFromPV && // cut scaling with track density distance < clusterMaxDistance); // absolute distance cut #ifdef VTXDEBUG std::cout << tt->trackBaseRef().key() << " : " << (selected?"+":" ")<< " " << m.significance() << " < " << clusterMaxSignificance << " && " << dotprodSeed << " > " << clusterMinAngleCosine << " && " << dotprodTrack << " > " << clusterMinAngleCosine << " && " << dotprodTrackSeed2D << " > " << clusterMinAngleCosine << " && " << distance*clusterScale << " < " << densityFactor*distanceFromPV << " crossingtoPV: " << distanceFromPV << " dis*scal " << distance*clusterScale << " < " << densityFactor*distanceFromPV << " dist: " << distance << " < " << clusterMaxDistance << std::endl; // cut scaling with track density #endif if(selected) { result.push_back(*tt); seedingPoint = GlobalPoint(cp.x()*w+seedingPoint.x(),cp.y()*w+seedingPoint.y(),cp.z()*w+seedingPoint.z()); sumWeights+=w; } } } seedingPoint = GlobalPoint(seedingPoint.x()/sumWeights,seedingPoint.y()/sumWeights,seedingPoint.z()/sumWeights); return std::pair<std::vector<reco::TransientTrack>,GlobalPoint>(result,seedingPoint); }
bool TracksClusteringFromDisplacedSeed::trackFilter | ( | const reco::TrackRef & | track | ) | const [private] |
double TracksClusteringFromDisplacedSeed::clusterMaxDistance [private] |
Definition at line 49 of file TracksClusteringFromDisplacedSeed.h.
Referenced by nearTracks().
double TracksClusteringFromDisplacedSeed::clusterMaxSignificance [private] |
Definition at line 50 of file TracksClusteringFromDisplacedSeed.h.
Referenced by nearTracks().
double TracksClusteringFromDisplacedSeed::clusterMinAngleCosine [private] |
Definition at line 52 of file TracksClusteringFromDisplacedSeed.h.
Referenced by nearTracks().
double TracksClusteringFromDisplacedSeed::clusterScale [private] |
Definition at line 51 of file TracksClusteringFromDisplacedSeed.h.
Referenced by nearTracks().
double TracksClusteringFromDisplacedSeed::min3DIPSignificance [private] |
Definition at line 47 of file TracksClusteringFromDisplacedSeed.h.
Referenced by clusters().
double TracksClusteringFromDisplacedSeed::min3DIPValue [private] |
Definition at line 48 of file TracksClusteringFromDisplacedSeed.h.
Referenced by clusters().