Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:06

0001 #include "RecoVertex/ConfigurableVertexReco/interface/ConfigurableAdaptiveReconstructor.h"
0002 #include "RecoVertex/AdaptiveVertexFinder/interface/AdaptiveVertexReconstructor.h"
0003 
0004 using namespace std;
0005 
0006 namespace {
0007   edm::ParameterSet mydefaults() {
0008     edm::ParameterSet ret;
0009     ret.addParameter<double>("primcut", 2.0);
0010     ret.addParameter<double>("primT", 256.0);
0011     ret.addParameter<double>("primr", 0.25);
0012     ret.addParameter<double>("seccut", 6.0);
0013     ret.addParameter<double>("secT", 256.0);
0014     ret.addParameter<double>("secr", 0.25);
0015     ret.addParameter<double>("minweight", 0.5);
0016     ret.addParameter<double>("weightthreshold", 0.001);
0017     ret.addParameter<bool>("smoothing", false);
0018     return ret;
0019   }
0020 }  // namespace
0021 
0022 ConfigurableAdaptiveReconstructor::ConfigurableAdaptiveReconstructor() : theRector(new AdaptiveVertexReconstructor()) {}
0023 
0024 void ConfigurableAdaptiveReconstructor::configure(const edm::ParameterSet& n) {
0025   edm::ParameterSet m = n;
0026   m.augment(mydefaults());
0027   if (theRector)
0028     delete theRector;
0029   theRector = new AdaptiveVertexReconstructor(m);
0030 }
0031 
0032 ConfigurableAdaptiveReconstructor::~ConfigurableAdaptiveReconstructor() {
0033   if (theRector)
0034     delete theRector;
0035 }
0036 
0037 ConfigurableAdaptiveReconstructor::ConfigurableAdaptiveReconstructor(const ConfigurableAdaptiveReconstructor& o)
0038     : theRector(o.theRector->clone()) {}
0039 
0040 ConfigurableAdaptiveReconstructor* ConfigurableAdaptiveReconstructor::clone() const {
0041   return new ConfigurableAdaptiveReconstructor(*this);
0042 }
0043 
0044 vector<TransientVertex> ConfigurableAdaptiveReconstructor::vertices(const std::vector<reco::TransientTrack>& t) const {
0045   return theRector->vertices(t);
0046 }
0047 
0048 vector<TransientVertex> ConfigurableAdaptiveReconstructor::vertices(const std::vector<reco::TransientTrack>& t,
0049                                                                     const reco::BeamSpot& s) const {
0050   return theRector->vertices(t, s);
0051 }
0052 
0053 vector<TransientVertex> ConfigurableAdaptiveReconstructor::vertices(const std::vector<reco::TransientTrack>& prims,
0054                                                                     const std::vector<reco::TransientTrack>& secs,
0055                                                                     const reco::BeamSpot& s) const {
0056   return theRector->vertices(prims, secs, s);
0057 }
0058 
0059 edm::ParameterSet ConfigurableAdaptiveReconstructor::defaults() const { return mydefaults(); }
0060 
0061 #include "RecoVertex/ConfigurableVertexReco/interface/ConfRecoBuilder.h"
0062 
0063 namespace {
0064   const ConfRecoBuilder<ConfigurableAdaptiveReconstructor> t("avr", "Adaptive Vertex Reconstructor [ = Iterative avf]");
0065   // ConfRecoBuilder < ConfigurableAdaptiveReconstructor > s ( "default", "Adaptive Vertex Reconstructor [ = Iterative avf]" );
0066 }  // namespace