File indexing completed on 2024-04-06 12:29:06
0001 #include "RecoVertex/ConfigurableVertexReco/interface/ConfigurableTrimmedKalmanFinder.h"
0002 #include "RecoVertex/TrimmedKalmanVertexFinder/interface/KalmanTrimmedVertexFinder.h"
0003
0004 namespace {
0005 edm::ParameterSet mydefaults() {
0006 edm::ParameterSet ret;
0007 ret.addParameter<double>("ptcut", 0.);
0008 ret.addParameter<double>("trkcutpv", 0.05);
0009 ret.addParameter<double>("trkcutsv", 0.01);
0010 ret.addParameter<double>("vtxcut", 0.01);
0011 return ret;
0012 }
0013 }
0014
0015 ConfigurableTrimmedKalmanFinder::ConfigurableTrimmedKalmanFinder() : theRector(new KalmanTrimmedVertexFinder()) {}
0016
0017 void ConfigurableTrimmedKalmanFinder::configure(const edm::ParameterSet& n) {
0018 if (theRector)
0019 delete theRector;
0020 edm::ParameterSet m = n;
0021 m.augment(mydefaults());
0022 KalmanTrimmedVertexFinder* tmp = new KalmanTrimmedVertexFinder();
0023 tmp->setPtCut(m.getParameter<double>("ptcut"));
0024 tmp->setTrackCompatibilityCut(m.getParameter<double>("trkcutpv"));
0025 tmp->setTrackCompatibilityToSV(m.getParameter<double>("trkcutsv"));
0026 tmp->setVertexFitProbabilityCut(m.getParameter<double>("vtxcut"));
0027 theRector = tmp;
0028 }
0029
0030 ConfigurableTrimmedKalmanFinder::~ConfigurableTrimmedKalmanFinder() {
0031 if (theRector)
0032 delete theRector;
0033 }
0034
0035 ConfigurableTrimmedKalmanFinder::ConfigurableTrimmedKalmanFinder(const ConfigurableTrimmedKalmanFinder& o)
0036 : theRector(o.theRector->clone()) {}
0037
0038 ConfigurableTrimmedKalmanFinder* ConfigurableTrimmedKalmanFinder::clone() const {
0039 return new ConfigurableTrimmedKalmanFinder(*this);
0040 }
0041
0042 std::vector<TransientVertex> ConfigurableTrimmedKalmanFinder::vertices(const std::vector<reco::TransientTrack>& t,
0043 const reco::BeamSpot& s) const {
0044 return theRector->vertices(t, s);
0045 }
0046
0047 std::vector<TransientVertex> ConfigurableTrimmedKalmanFinder::vertices(const std::vector<reco::TransientTrack>& prims,
0048 const std::vector<reco::TransientTrack>& secs,
0049 const reco::BeamSpot& s) const {
0050 return theRector->vertices(prims, secs, s);
0051 }
0052
0053 std::vector<TransientVertex> ConfigurableTrimmedKalmanFinder::vertices(
0054 const std::vector<reco::TransientTrack>& t) const {
0055 return theRector->vertices(t);
0056 }
0057
0058 edm::ParameterSet ConfigurableTrimmedKalmanFinder::defaults() const { return mydefaults(); }
0059
0060 #include "RecoVertex/ConfigurableVertexReco/interface/ConfRecoBuilder.h"
0061
0062 namespace {
0063 const ConfRecoBuilder<ConfigurableTrimmedKalmanFinder> t("tkf", "Trimmed Kalman Vertex Finder");
0064 }