#include <AdaptiveVertexReconstructor.h>
Public Member Functions | |
AdaptiveVertexReconstructor (float primcut=2.0, float seccut=6.0, float minweight=0.5, bool smoothing=false) | |
AdaptiveVertexReconstructor (const edm::ParameterSet &s) | |
virtual AdaptiveVertexReconstructor * | clone () const |
std::vector< TransientVertex > | vertices (const std::vector< reco::TransientTrack > &, const reco::BeamSpot &) const |
std::vector< TransientVertex > | vertices (const std::vector< reco::TransientTrack > &v) const |
std::vector< TransientVertex > | vertices (const std::vector< reco::TransientTrack > &primaries, const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &) const |
~AdaptiveVertexReconstructor () | |
Private Member Functions | |
TransientVertex | cleanUp (const TransientVertex &old) const |
std::vector< TransientVertex > | cleanUpVertices (const std::vector< TransientVertex > &) const |
void | erase (const TransientVertex &newvtx, std::set< reco::TransientTrack > &remainingtrks, float w) const |
void | setupFitters (float primcut, float primT, float primr, float seccut, float secT, float secr, bool smoothing) |
std::vector< TransientVertex > | vertices (const std::vector< reco::TransientTrack > &primaries, const std::vector< reco::TransientTrack > &trks, const reco::BeamSpot &, bool has_primaries, bool usespot) const |
Private Attributes | |
float | theMinWeight |
AdaptiveVertexFitter * | thePrimaryFitter |
AdaptiveVertexFitter * | theSecondaryFitter |
float | theWeightThreshold |
Definition at line 9 of file AdaptiveVertexReconstructor.h.
AdaptiveVertexReconstructor::AdaptiveVertexReconstructor | ( | float | primcut = 2.0 , |
float | seccut = 6.0 , |
||
float | minweight = 0.5 , |
||
bool | smoothing = false |
||
) |
Definition at line 89 of file AdaptiveVertexReconstructor.cc.
References setupFitters().
Referenced by clone().
: thePrimaryFitter ( 0 ), theSecondaryFitter ( 0 ), theMinWeight( min_weight ), theWeightThreshold ( 0.001 ) { setupFitters ( primcut, 256., 0.25, seccut, 256., 0.25, smoothing ); }
AdaptiveVertexReconstructor::~AdaptiveVertexReconstructor | ( | ) |
Definition at line 97 of file AdaptiveVertexReconstructor.cc.
References thePrimaryFitter, and theSecondaryFitter.
{ if ( thePrimaryFitter ) delete thePrimaryFitter; if ( theSecondaryFitter ) delete theSecondaryFitter; }
AdaptiveVertexReconstructor::AdaptiveVertexReconstructor | ( | const edm::ParameterSet & | s | ) |
The ParameterSet should have the following defined: double primcut double seccut double minweight for descriptions see
Definition at line 136 of file AdaptiveVertexReconstructor.cc.
References edm::ParameterSet::getParameter(), setupFitters(), theMinWeight, theWeightThreshold, and cms::Exception::what().
: thePrimaryFitter(0), theSecondaryFitter(0), theMinWeight(0.5), theWeightThreshold ( 0.001 ) { float primcut = 2.0; float seccut = 6.0; bool smoothing=false; // float primT = 4096.; // float primr = 0.125; float primT = 256.; float primr = 0.25; float secT = 256.; float secr = 0.25; try { primcut = m.getParameter<double>("primcut"); primT = m.getParameter<double>("primT"); primr = m.getParameter<double>("primr"); seccut = m.getParameter<double>("seccut"); secT = m.getParameter<double>("secT"); secr = m.getParameter<double>("secr"); theMinWeight = m.getParameter<double>("minweight"); theWeightThreshold = m.getParameter<double>("weightthreshold"); smoothing = m.getParameter<bool>("smoothing"); } catch ( edm::Exception & e ) { edm::LogError ("AdaptiveVertexReconstructor") << e.what(); } setupFitters ( primcut, primT, primr, seccut, secT, secr, smoothing ); }
TransientVertex AdaptiveVertexReconstructor::cleanUp | ( | const TransientVertex & | old | ) | const [private] |
Definition at line 11 of file AdaptiveVertexReconstructor.cc.
References TransientVertex::degreesOfFreedom(), TransientVertex::hasPrior(), TransientVertex::hasRefittedTracks(), i, TransientVertex::originalTracks(), TransientVertex::priorError(), TransientVertex::priorPosition(), TransientVertex::refittedTracks(), runTheMatrix::ret, TransientVertex::totalChiSquared(), TransientVertex::trackWeight(), TransientVertex::vertexState(), and TransientVertex::weightMap().
{ vector < reco::TransientTrack > trks = old.originalTracks(); vector < reco::TransientTrack > newtrks; TransientVertex::TransientTrackToFloatMap mp; static const float minweight = 1.e-8; // discard all tracks with lower weight for ( vector< reco::TransientTrack >::const_iterator i=trks.begin(); i!=trks.end() ; ++i ) { if ( old.trackWeight ( *i ) > minweight ) { newtrks.push_back ( *i ); mp[*i]=old.trackWeight ( *i ); } } TransientVertex ret; if ( old.hasPrior() ) { VertexState priorstate ( old.priorPosition(), old.priorError() ); ret=TransientVertex ( priorstate, old.vertexState(), newtrks, old.totalChiSquared(), old.degreesOfFreedom() ); } else { ret=TransientVertex ( old.vertexState(), newtrks, old.totalChiSquared(), old.degreesOfFreedom() ); } ret.weightMap ( mp ); // set weight map if ( old.hasRefittedTracks() ) { // we have refitted tracks -- copy them! vector < reco::TransientTrack > newrfs; vector < reco::TransientTrack > oldrfs=old.refittedTracks(); vector < reco::TransientTrack >::const_iterator origtrkiter=trks.begin(); for ( vector< reco::TransientTrack >::const_iterator i=oldrfs.begin(); i!=oldrfs.end() ; ++i ) { if ( old.trackWeight ( *origtrkiter ) > minweight ) { newrfs.push_back ( *i ); } origtrkiter++; } if ( newrfs.size() ) ret.refittedTracks ( newrfs ); // copy refitted tracks } if ( ret.refittedTracks().size() > ret.originalTracks().size() ) { edm::LogError("AdaptiveVertexReconstructor" ) << "More refitted tracks than original tracks!"; } return ret; }
vector< TransientVertex > AdaptiveVertexReconstructor::cleanUpVertices | ( | const std::vector< TransientVertex > & | ) | const [private] |
cleanup reconstructed vertices. discard all with too few significant tracks.
Definition at line 268 of file AdaptiveVertexReconstructor.cc.
References i, runTheMatrix::ret, and theWeightThreshold.
{ vector < TransientVertex > ret; for ( vector< TransientVertex >::const_iterator i=old.begin(); i!=old.end() ; ++i ) { if (!(i->hasTrackWeight())) { // if we dont have track weights, we take the vtx ret.push_back ( *i ); continue; } // maybe this should be replaced with asking for the ndf ... // e.g. if ( ndf > - 1. ) int nsig=0; // number of significant tracks. TransientVertex::TransientTrackToFloatMap wm = i->weightMap(); for ( TransientVertex::TransientTrackToFloatMap::const_iterator w=wm.begin(); w!=wm.end() ; ++w ) { if (w->second > theWeightThreshold) nsig++; } if ( nsig > 1 ) ret.push_back ( *i ); } return ret; }
virtual AdaptiveVertexReconstructor* AdaptiveVertexReconstructor::clone | ( | void | ) | const [inline, virtual] |
Implements VertexReconstructor.
Definition at line 45 of file AdaptiveVertexReconstructor.h.
References AdaptiveVertexReconstructor().
{ return new AdaptiveVertexReconstructor( * this ); }
void AdaptiveVertexReconstructor::erase | ( | const TransientVertex & | newvtx, |
std::set< reco::TransientTrack > & | remainingtrks, | ||
float | w | ||
) | const [private] |
contrary to what its name has you believe, ::erase removes all newvtx.originalTracks() above theMinWeight from remainingtrks.
void AdaptiveVertexReconstructor::setupFitters | ( | float | primcut, |
float | primT, | ||
float | primr, | ||
float | seccut, | ||
float | secT, | ||
float | secr, | ||
bool | smoothing | ||
) | [private] |
setup the vertex fitters.
Definition at line 103 of file AdaptiveVertexReconstructor.cc.
References AdaptiveVertexFitter::setWeightThreshold(), thePrimaryFitter, theSecondaryFitter, and theWeightThreshold.
Referenced by AdaptiveVertexReconstructor().
{ VertexSmoother<5> * smoother ; if ( smoothing ) { smoother = new KalmanVertexSmoother(); } else { smoother = new DummyVertexSmoother<5>(); } if ( thePrimaryFitter ) delete thePrimaryFitter; if ( theSecondaryFitter ) delete theSecondaryFitter; /* edm::LogError ("AdaptiveVertexReconstructor" ) << "Tini and r are hardcoded now!"; */ thePrimaryFitter = new AdaptiveVertexFitter ( GeometricAnnealing ( primcut, primT, primr ), DefaultLinearizationPointFinder(), KalmanVertexUpdator<5>(), KalmanVertexTrackCompatibilityEstimator<5>(), *smoother ); thePrimaryFitter->setWeightThreshold ( theWeightThreshold ); // if the primary fails, sth is wrong, so here we set a threshold on the weight. theSecondaryFitter = new AdaptiveVertexFitter ( GeometricAnnealing ( seccut, secT, secr ), DefaultLinearizationPointFinder(), KalmanVertexUpdator<5>(), KalmanVertexTrackCompatibilityEstimator<5>(), *smoother ); theSecondaryFitter->setWeightThreshold ( 0. ); // need to set it or else we have // unwanted exceptions to deal with. // cleanup can come later! delete smoother; }
std::vector<TransientVertex> AdaptiveVertexReconstructor::vertices | ( | const std::vector< reco::TransientTrack > & | primaries, |
const std::vector< reco::TransientTrack > & | trks, | ||
const reco::BeamSpot & | , | ||
bool | has_primaries, | ||
bool | usespot | ||
) | const [private] |
the actual fit to avoid code duplication
std::vector<TransientVertex> AdaptiveVertexReconstructor::vertices | ( | const std::vector< reco::TransientTrack > & | ) | const [virtual] |
Reconstruct vertices
Implements VertexReconstructor.
std::vector<TransientVertex> AdaptiveVertexReconstructor::vertices | ( | const std::vector< reco::TransientTrack > & | t, |
const reco::BeamSpot & | |||
) | const [virtual] |
Reconstruct vertices, exploiting the beamspot constraint for the primary vertex
Reimplemented from VertexReconstructor.
std::vector<TransientVertex> AdaptiveVertexReconstructor::vertices | ( | const std::vector< reco::TransientTrack > & | primaries, |
const std::vector< reco::TransientTrack > & | tracks, | ||
const reco::BeamSpot & | spot | ||
) | const [virtual] |
Reconstruct vertices, but exploit the fact that you know that some tracks cannot come from a secondary vertex. primaries Tracks that _cannot_ come from a secondary vertex (but can, in principle, be non-primaries, also). tracks These are the tracks that are of unknown origin. These tracks are subjected to pattern recognition. spot A beamspot constraint is mandatory in this method.
Reimplemented from VertexReconstructor.
float AdaptiveVertexReconstructor::theMinWeight [private] |
Definition at line 84 of file AdaptiveVertexReconstructor.h.
Referenced by AdaptiveVertexReconstructor().
Definition at line 80 of file AdaptiveVertexReconstructor.h.
Referenced by setupFitters(), and ~AdaptiveVertexReconstructor().
Definition at line 81 of file AdaptiveVertexReconstructor.h.
Referenced by setupFitters(), and ~AdaptiveVertexReconstructor().
float AdaptiveVertexReconstructor::theWeightThreshold [private] |
Definition at line 86 of file AdaptiveVertexReconstructor.h.
Referenced by AdaptiveVertexReconstructor(), cleanUpVertices(), and setupFitters().