CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelVertexProducer.cc
Go to the documentation of this file.
9 #include <memory>
10 #include <string>
11 #include <cmath>
12 
14  : conf_(conf), verbose_(0), dvf_(0), ptMin_(1.0)
15 {
16  // Register my product
17  produces<reco::VertexCollection>();
18 
19  // Setup shop
20  verbose_ = conf.getParameter<int>("Verbosity"); // 0 silent, 1 chatty, 2 loud
21  std::string finder = conf.getParameter<std::string>("Finder"); // DivisiveVertexFinder
22  bool useError = conf.getParameter<bool>("UseError"); // true
23  bool wtAverage = conf.getParameter<bool>("WtAverage"); // true
24  double zOffset = conf.getParameter<double>("ZOffset"); // 5.0 sigma
25  double zSeparation = conf.getParameter<double>("ZSeparation"); // 0.05 cm
26  int ntrkMin = conf.getParameter<int>("NTrkMin"); // 3
27  // Tracking requirements before sending a track to be considered for vtx
28  ptMin_ = conf_.getParameter<double>("PtMin"); // 1.0 GeV
29 
30  if (finder == "DivisiveVertexFinder") {
31  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Using the DivisiveVertexFinder\n";
32  dvf_ = new DivisiveVertexFinder(zOffset, ntrkMin, useError, zSeparation, wtAverage, verbose_);
33  }
34  else { // Finder not supported, or you made a mistake in your request
35  // throw an exception once I figure out how CMSSW does this
36  }
37 }
38 
39 
41  delete dvf_;
42 }
43 
45 
46  // First fish the pixel tracks out of the event
47  edm::Handle<reco::TrackCollection> trackCollection;
48  edm::InputTag trackCollName = conf_.getParameter<edm::InputTag>("TrackCollection");
49  e.getByLabel(trackCollName,trackCollection);
50  const reco::TrackCollection tracks = *(trackCollection.product());
51  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Found " << tracks.size() << " tracks in TrackCollection called " << trackCollName << "\n";
52 
53 
54  // Second, make a collection of pointers to the tracks we want for the vertex finder
56  for (unsigned int i=0; i<tracks.size(); i++) {
57  if (tracks[i].pt() > ptMin_)
58  trks.push_back( reco::TrackRef(trackCollection, i) );
59  }
60  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Selected " << trks.size() << " of these tracks for vertexing\n";
61 
63  edm::InputTag bsName = conf_.getParameter<edm::InputTag>("beamSpot");
64  e.getByLabel(bsName,bsHandle);
65  math::XYZPoint myPoint(0.,0.,0.);
66  if (bsHandle.isValid()) myPoint = math::XYZPoint(bsHandle->x0(),bsHandle->y0(), 0. ); //FIXME: fix last coordinate with vertex.z() at same time
67 
68  // Third, ship these tracks off to be vertexed
69  std::auto_ptr<reco::VertexCollection> vertexes(new reco::VertexCollection);
70  bool ok;
71  if (conf_.getParameter<bool>("Method2")) {
72  ok = dvf_->findVertexesAlt(trks, // input
73  *vertexes,myPoint); // output
74  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method2 returned status of " << ok;
75  }
76  else {
77  ok = dvf_->findVertexes(trks, // input
78  *vertexes); // output
79  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method1 returned status of " << ok;
80  }
81 
82  if (verbose_ > 0) {
83  edm::LogInfo("PixelVertexProducer") << ": Found " << vertexes->size() << " vertexes\n";
84  for (unsigned int i=0; i<vertexes->size(); ++i) {
85  edm::LogInfo("PixelVertexProducer") << "Vertex number " << i << " has " << (*vertexes)[i].tracksSize() << " tracks with a position of " << (*vertexes)[i].z() << " +- " << std::sqrt( (*vertexes)[i].covariance(2,2) );
86  }
87  }
88 
89 
90  if(bsHandle.isValid())
91  {
92  const reco::BeamSpot & bs = *bsHandle;
93 
94  for (unsigned int i=0; i<vertexes->size(); ++i) {
95  double z=(*vertexes)[i].z();
96  double x=bs.x0()+bs.dxdz()*(z-bs.z0());
97  double y=bs.y0()+bs.dydz()*(z-bs.z0());
98  reco::Vertex v( reco::Vertex::Point(x,y,z), (*vertexes)[i].error(),(*vertexes)[i].chi2() , (*vertexes)[i].ndof() , (*vertexes)[i].tracksSize());
99  //Copy also the tracks
100  for (std::vector<reco::TrackBaseRef >::const_iterator it = (*vertexes)[i].tracks_begin();
101  it !=(*vertexes)[i].tracks_end(); it++) {
102  v.add( *it );
103  }
104  (*vertexes)[i]=v;
105 
106  }
107  }
108  else
109  {
110  edm::LogWarning("PixelVertexProducer") << "No beamspot found. Using returning vertexes with (0,0,Z) ";
111  }
112 
113  e.put(vertexes);
114 }
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:69
int i
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &, const edm::EventSetup &)
PixelVertexProducer(const edm::ParameterSet &)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool findVertexesAlt(const reco::TrackRefVector &trks, reco::VertexCollection &vertexes, const math::XYZPoint &bs)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::ParameterSet conf_
double double double z
double dydz() const
dydz slope
Definition: BeamSpot.h:85
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
T sqrt(T t)
Definition: SSEVec.h:28
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
bool findVertexes(const reco::TrackRefVector &trks, reco::VertexCollection &vertexes)
Run the divisive algorithm and return a vector of vertexes for the input track collection.
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
double dxdz() const
dxdz slope
Definition: BeamSpot.h:83
tuple conf
Definition: dbtoconf.py:185
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
tuple tracks
Definition: testEve_cfg.py:39
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
T const * product() const
Definition: Handle.h:74
DivisiveVertexFinder * dvf_
double y0() const
y coordinate
Definition: BeamSpot.h:67
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:60
size_type size() const
Size of the RefVector.
Definition: RefVector.h:84
Definition: DDAxes.h:10
mathSSE::Vec4< T > v
double x0() const
x coordinate
Definition: BeamSpot.h:65