#include <PixelVertexFinding/interface/PixelVertexProducer.h>
Public Member Functions | |
PixelVertexProducer (const edm::ParameterSet &) | |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
~PixelVertexProducer () | |
Private Attributes | |
edm::ParameterSet | conf_ |
DivisiveVertexFinder * | dvf_ |
double | ptMin_ |
int | verbose_ |
Implementation: This producer can use either the Divisive Primary Vertex Finder or the Histogramming Primary Vertex Finder (currently not implemented). It relies on the PixelTripletProducer and PixelTrackProducer having already been run upstream. This is code ported from ORCA originally written by S Cucciarelli, M Konecki, D Kotlinski.
Definition at line 34 of file PixelVertexProducer.h.
PixelVertexProducer::PixelVertexProducer | ( | const edm::ParameterSet & | conf | ) | [explicit] |
Definition at line 13 of file PixelVertexProducer.cc.
References conf_, dvf_, edm::ParameterSet::getParameter(), ptMin_, and verbose_.
00014 : conf_(conf), verbose_(0), dvf_(0), ptMin_(1.0) 00015 { 00016 // Register my product 00017 produces<reco::VertexCollection>(); 00018 00019 // Setup shop 00020 verbose_ = conf.getParameter<int>("Verbosity"); // 0 silent, 1 chatty, 2 loud 00021 std::string finder = conf.getParameter<std::string>("Finder"); // DivisiveVertexFinder 00022 bool useError = conf.getParameter<bool>("UseError"); // true 00023 bool wtAverage = conf.getParameter<bool>("WtAverage"); // true 00024 double zOffset = conf.getParameter<double>("ZOffset"); // 5.0 sigma 00025 double zSeparation = conf.getParameter<double>("ZSeparation"); // 0.05 cm 00026 int ntrkMin = conf.getParameter<int>("NTrkMin"); // 3 00027 // Tracking requirements before sending a track to be considered for vtx 00028 ptMin_ = conf_.getParameter<double>("PtMin"); // 1.0 GeV 00029 00030 if (finder == "DivisiveVertexFinder") { 00031 if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Using the DivisiveVertexFinder\n"; 00032 dvf_ = new DivisiveVertexFinder(zOffset, ntrkMin, useError, zSeparation, wtAverage, verbose_); 00033 } 00034 else { // Finder not supported, or you made a mistake in your request 00035 // throw an exception once I figure out how CMSSW does this 00036 } 00037 }
PixelVertexProducer::~PixelVertexProducer | ( | ) |
Definition at line 40 of file PixelVertexProducer.cc.
References dvf_.
00040 { 00041 delete dvf_; 00042 }
void PixelVertexProducer::produce | ( | edm::Event & | e, | |
const edm::EventSetup & | es | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 44 of file PixelVertexProducer.cc.
References reco::Vertex::add(), conf_, dvf_, reco::BeamSpot::dxdz(), reco::BeamSpot::dydz(), error, DivisiveVertexFinder::findVertexes(), DivisiveVertexFinder::findVertexesAlt(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), i, edm::Handle< T >::isValid(), it, edm::Handle< T >::product(), ptMin_, edm::RefVector< C, T, F >::push_back(), edm::Event::put(), edm::RefVector< C, T, F >::size(), funct::sqrt(), tracks, v, verbose_, x, reco::BeamSpot::x0(), y, reco::BeamSpot::y0(), z, and reco::BeamSpot::z0().
00044 { 00045 00046 // First fish the pixel tracks out of the event 00047 edm::Handle<reco::TrackCollection> trackCollection; 00048 edm::InputTag trackCollName = conf_.getParameter<edm::InputTag>("TrackCollection"); 00049 e.getByLabel(trackCollName,trackCollection); 00050 const reco::TrackCollection tracks = *(trackCollection.product()); 00051 if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Found " << tracks.size() << " tracks in TrackCollection called " << trackCollName << "\n"; 00052 00053 00054 // Second, make a collection of pointers to the tracks we want for the vertex finder 00055 reco::TrackRefVector trks; 00056 for (unsigned int i=0; i<tracks.size(); i++) { 00057 if (tracks[i].pt() > ptMin_) 00058 trks.push_back( reco::TrackRef(trackCollection, i) ); 00059 } 00060 if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Selected " << trks.size() << " of these tracks for vertexing\n"; 00061 00062 // Third, ship these tracks off to be vertexed 00063 std::auto_ptr<reco::VertexCollection> vertexes(new reco::VertexCollection); 00064 bool ok; 00065 if (conf_.getParameter<bool>("Method2")) { 00066 ok = dvf_->findVertexesAlt(trks, // input 00067 *vertexes); // output 00068 if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method2 returned status of " << ok; 00069 } 00070 else { 00071 ok = dvf_->findVertexes(trks, // input 00072 *vertexes); // output 00073 if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method1 returned status of " << ok; 00074 } 00075 00076 if (verbose_ > 0) { 00077 edm::LogInfo("PixelVertexProducer") << ": Found " << vertexes->size() << " vertexes\n"; 00078 for (unsigned int i=0; i<vertexes->size(); ++i) { 00079 edm::LogInfo("PixelVertexProducer") << "Vertex number " << i << " has " << (*vertexes)[i].tracksSize() << " tracks with a position of " << (*vertexes)[i].z() << " +- " << std::sqrt( (*vertexes)[i].error(2,2) ); 00080 } 00081 } 00082 00083 edm::Handle<reco::BeamSpot> bsHandle; 00084 edm::InputTag bsName = conf_.getParameter<edm::InputTag>("beamSpot"); 00085 e.getByLabel(bsName,bsHandle); 00086 if(bsHandle.isValid()) 00087 { 00088 const reco::BeamSpot & bs = *bsHandle; 00089 00090 for (unsigned int i=0; i<vertexes->size(); ++i) { 00091 double z=(*vertexes)[i].z(); 00092 double x=bs.x0()+bs.dxdz()*(z-bs.z0()); 00093 double y=bs.y0()+bs.dydz()*(z-bs.z0()); 00094 reco::Vertex v( reco::Vertex::Point(x,y,z), (*vertexes)[i].error(),(*vertexes)[i].chi2() , (*vertexes)[i].ndof() , (*vertexes)[i].tracksSize()); 00095 //Copy also the tracks 00096 for (std::vector<reco::TrackBaseRef >::const_iterator it = (*vertexes)[i].tracks_begin(); 00097 it !=(*vertexes)[i].tracks_end(); it++) { 00098 v.add( *it ); 00099 } 00100 (*vertexes)[i]=v; 00101 00102 } 00103 } 00104 else 00105 { 00106 edm::LogWarning("PixelVertexProducer") << "No beamspot found. Using returning vertexes with (0,0,Z) "; 00107 } 00108 00109 e.put(vertexes); 00110 }
edm::ParameterSet PixelVertexProducer::conf_ [private] |
Definition at line 42 of file PixelVertexProducer.h.
Referenced by PixelVertexProducer(), and produce().
DivisiveVertexFinder* PixelVertexProducer::dvf_ [private] |
Definition at line 45 of file PixelVertexProducer.h.
Referenced by PixelVertexProducer(), produce(), and ~PixelVertexProducer().
double PixelVertexProducer::ptMin_ [private] |
Definition at line 47 of file PixelVertexProducer.h.
Referenced by PixelVertexProducer(), and produce().
int PixelVertexProducer::verbose_ [private] |
Definition at line 44 of file PixelVertexProducer.h.
Referenced by PixelVertexProducer(), and produce().