CMS 3D CMS Logo

AdaptiveVertexReconstructor Class Reference

#include <RecoVertex/AdaptiveVertexFinder/interface/AdaptiveVertexReconstructor.h>

Inheritance diagram for AdaptiveVertexReconstructor:

VertexReconstructor

List of all members.

Public Member Functions

 AdaptiveVertexReconstructor (const edm::ParameterSet &s)
 The ParameterSet should have the following defined: double primcut double seccut double minweight for descriptions see.
 AdaptiveVertexReconstructor (float primcut=2.0, float seccut=6.0, float minweight=0.5, bool smoothing=false)
virtual
AdaptiveVertexReconstructor
clone () const
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &primaries, const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &) const
 Reconstruct vertices, but exploit the fact that you know that some tracks cannot come from a secondary vertex.
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &, const reco::BeamSpot &) const
 Reconstruct vertices, exploiting the beamspot constraint for the primary vertex.
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &v) const
 Reconstruct vertices.
 ~AdaptiveVertexReconstructor ()

Private Member Functions

TransientVertex cleanUp (const TransientVertex &old) const
std::vector< TransientVertexcleanUpVertices (const std::vector< TransientVertex > &) const
 cleanup reconstructed vertices.
void erase (const TransientVertex &newvtx, std::set< reco::TransientTrack > &remainingtrks, float w) const
 contrary to what its name has you believe, erase removes all newvtx.originalTracks() above theMinWeight from remainingtrks.
void setupFitters (float primcut, float primT, float primr, float seccut, float secT, float secr, bool smoothing)
 setup the vertex fitters.
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &primaries, const std::vector< reco::TransientTrack > &trks, const reco::BeamSpot &, bool has_primaries, bool usespot) const
 the actual fit to avoid code duplication

Private Attributes

float theMinWeight
AdaptiveVertexFitterthePrimaryFitter
AdaptiveVertexFittertheSecondaryFitter
float theWeightThreshold


Detailed Description

Definition at line 9 of file AdaptiveVertexReconstructor.h.


Constructor & Destructor Documentation

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().

00090                                                                     :
00091   thePrimaryFitter ( 0 ), theSecondaryFitter ( 0 ),
00092        theMinWeight( min_weight ), theWeightThreshold ( 0.001 )
00093 {
00094   setupFitters ( primcut, 256., 0.25, seccut, 256., 0.25, smoothing );
00095 }

AdaptiveVertexReconstructor::~AdaptiveVertexReconstructor (  ) 

Definition at line 97 of file AdaptiveVertexReconstructor.cc.

References thePrimaryFitter, and theSecondaryFitter.

00098 {
00099   if ( thePrimaryFitter ) delete thePrimaryFitter;
00100   if ( theSecondaryFitter ) delete theSecondaryFitter;
00101 }

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 e, edm::ParameterSet::getParameter(), setupFitters(), theMinWeight, theWeightThreshold, and cms::Exception::what().

00137   : thePrimaryFitter(0), theSecondaryFitter(0), theMinWeight(0.5), theWeightThreshold ( 0.001 )
00138 {
00139   float primcut = 2.0;
00140   float seccut = 6.0;
00141   bool smoothing=false;
00142   // float primT = 4096.;
00143   // float primr = 0.125;
00144   float primT = 256.;
00145   float primr = 0.25;
00146   float secT = 256.;
00147   float secr = 0.25;
00148   
00149   try {
00150     primcut =  m.getParameter<double>("primcut");
00151     primT =  m.getParameter<double>("primT");
00152     primr =  m.getParameter<double>("primr");
00153     seccut =  m.getParameter<double>("seccut");
00154     secT =  m.getParameter<double>("secT");
00155     secr =  m.getParameter<double>("secr");
00156     theMinWeight = m.getParameter<double>("minweight");
00157     theWeightThreshold = m.getParameter<double>("weightthreshold");
00158     smoothing =  m.getParameter<bool>("smoothing");
00159   } catch ( edm::Exception & e ) {
00160     edm::LogError ("AdaptiveVertexReconstructor") << e.what();
00161   }
00162 
00163   setupFitters ( primcut, primT, primr, seccut, secT, secr, smoothing );
00164 }


Member Function Documentation

TransientVertex AdaptiveVertexReconstructor::cleanUp ( const TransientVertex old  )  const [private]

Definition at line 11 of file AdaptiveVertexReconstructor.cc.

References TransientVertex::degreesOfFreedom(), TransientVertex::hasPrior(), TransientVertex::hasRefittedTracks(), i, mp, TransientVertex::originalTracks(), TransientVertex::priorError(), TransientVertex::priorPosition(), TransientVertex::refittedTracks(), TransientVertex::totalChiSquared(), TransientVertex::trackWeight(), TransientVertex::vertexState(), and TransientVertex::weightMap().

00012 {
00013   vector < reco::TransientTrack > trks = old.originalTracks();
00014   vector < reco::TransientTrack > newtrks;
00015   TransientVertex::TransientTrackToFloatMap mp;
00016   static const float minweight = 1.e-8; // discard all tracks with lower weight
00017   for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
00018         i!=trks.end() ; ++i )
00019   {
00020     if ( old.trackWeight ( *i ) > minweight )
00021     {
00022       newtrks.push_back ( *i );
00023       mp[*i]=old.trackWeight ( *i );
00024     }
00025   }
00026 
00027   TransientVertex ret;
00028 
00029   if ( old.hasPrior() )
00030   {
00031     VertexState priorstate ( old.priorPosition(), old.priorError() );
00032     ret=TransientVertex ( priorstate, old.vertexState(), newtrks,
00033         old.totalChiSquared(), old.degreesOfFreedom() );
00034   } else {
00035     ret=TransientVertex ( old.vertexState(), newtrks,
00036                           old.totalChiSquared(), old.degreesOfFreedom() );
00037   }
00038   ret.weightMap ( mp ); // set weight map
00039 
00040   if ( old.hasRefittedTracks() )
00041   {
00042     // we have refitted tracks -- copy them!
00043     vector < reco::TransientTrack > newrfs;
00044     vector < reco::TransientTrack > oldrfs=old.refittedTracks();
00045     vector < reco::TransientTrack >::const_iterator origtrkiter=trks.begin();
00046     for ( vector< reco::TransientTrack >::const_iterator i=oldrfs.begin(); i!=oldrfs.end() ; ++i )
00047     {
00048       if ( old.trackWeight ( *origtrkiter ) > minweight )
00049       {
00050         newrfs.push_back ( *i );
00051       }
00052       origtrkiter++;
00053     }
00054     if ( newrfs.size() ) ret.refittedTracks ( newrfs ); // copy refitted tracks
00055   }
00056 
00057   if ( ret.refittedTracks().size() > ret.originalTracks().size() )
00058   {
00059     edm::LogError("AdaptiveVertexReconstructor" )
00060       << "More refitted tracks than original tracks!";
00061   }
00062 
00063   return ret;
00064 }

std::vector<TransientVertex> AdaptiveVertexReconstructor::cleanUpVertices ( const std::vector< TransientVertex > &   )  const [private]

cleanup reconstructed vertices.

discard all with too few significant tracks.

virtual AdaptiveVertexReconstructor* AdaptiveVertexReconstructor::clone ( void   )  const [inline, virtual]

Implements VertexReconstructor.

Definition at line 45 of file AdaptiveVertexReconstructor.h.

References AdaptiveVertexReconstructor().

00045                                                       {
00046     return new AdaptiveVertexReconstructor( * this );
00047   }

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().

00106 {
00107   VertexSmoother<5> * smoother ;
00108   if ( smoothing )
00109   {
00110     smoother = new KalmanVertexSmoother();
00111   } else {
00112     smoother = new DummyVertexSmoother<5>();
00113   }
00114 
00115   if ( thePrimaryFitter ) delete thePrimaryFitter;
00116   if ( theSecondaryFitter ) delete theSecondaryFitter;
00117 
00118   /*
00119   edm::LogError ("AdaptiveVertexReconstructor" )
00120     << "Tini and r are hardcoded now!";
00121     */
00122   thePrimaryFitter = new AdaptiveVertexFitter ( GeometricAnnealing ( primcut, primT, primr ), DefaultLinearizationPointFinder(),
00123       KalmanVertexUpdator<5>(), KalmanVertexTrackCompatibilityEstimator<5>(), *smoother );
00124   thePrimaryFitter->setWeightThreshold ( theWeightThreshold );
00125   // if the primary fails, sth is wrong, so here we set a threshold on the weight.
00126   theSecondaryFitter = new AdaptiveVertexFitter ( GeometricAnnealing ( seccut, secT, secr ), DefaultLinearizationPointFinder(),
00127       KalmanVertexUpdator<5>(), 
00128       KalmanVertexTrackCompatibilityEstimator<5>(), *smoother );
00129   theSecondaryFitter->setWeightThreshold ( 0. );
00130   // need to set it or else we have 
00131   // unwanted exceptions to deal with.
00132   // cleanup can come later!
00133   delete smoother;
00134 }

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 > &  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.

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 > &   )  const [virtual]

Reconstruct vertices.

Implements VertexReconstructor.


Member Data Documentation

float AdaptiveVertexReconstructor::theMinWeight [private]

Definition at line 84 of file AdaptiveVertexReconstructor.h.

Referenced by AdaptiveVertexReconstructor().

AdaptiveVertexFitter* AdaptiveVertexReconstructor::thePrimaryFitter [private]

Definition at line 80 of file AdaptiveVertexReconstructor.h.

Referenced by setupFitters(), and ~AdaptiveVertexReconstructor().

AdaptiveVertexFitter* AdaptiveVertexReconstructor::theSecondaryFitter [private]

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(), and setupFitters().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:13:52 2009 for CMSSW by  doxygen 1.5.4