CMS 3D CMS Logo

PixelVertexProducer.cc
Go to the documentation of this file.
9 #include <memory>
10 #include <string>
11 #include <cmath>
12 
14  : verbose_(conf.getParameter<int>("Verbosity") ) // 0 silent, 1 chatty, 2 loud
15  , ptMin_ (conf.getParameter<double>("PtMin") ) // 1.0 GeV
16  , method2 ( conf.getParameter<bool>("Method2") )
17  , trackCollName ( conf.getParameter<edm::InputTag>("TrackCollection") )
18  , token_Tracks ( consumes<reco::TrackCollection>(trackCollName) )
19  , token_BeamSpot ( consumes<reco::BeamSpot> (conf.getParameter<edm::InputTag>("beamSpot") ) )
20 {
21  // Register my product
22  produces<reco::VertexCollection>();
23 
24  // Setup shop
25  std::string finder = conf.getParameter<std::string>("Finder"); // DivisiveVertexFinder
26  bool useError = conf.getParameter<bool>("UseError"); // true
27  bool wtAverage = conf.getParameter<bool>("WtAverage"); // true
28  double zOffset = conf.getParameter<double>("ZOffset"); // 5.0 sigma
29  double zSeparation = conf.getParameter<double>("ZSeparation"); // 0.05 cm
30  int ntrkMin = conf.getParameter<int>("NTrkMin"); // 3
31  // Tracking requirements before sending a track to be considered for vtx
32 
33 
34  double track_pt_min = ptMin_;
35  double track_pt_max = 10.;
36  double track_chi2_max = 9999999.;
37  double track_prob_min = -1.;
38 
39  if ( conf.exists("PVcomparer") ) {
40  edm::ParameterSet PVcomparerPSet = conf.getParameter<edm::ParameterSet>("PVcomparer");
41  track_pt_min = PVcomparerPSet.getParameter<double>("track_pt_min");
42  if (track_pt_min != ptMin_) {
43  if (track_pt_min < ptMin_)
44  edm::LogInfo("PixelVertexProducer") << "minimum track pT setting differs between PixelVertexProducer (" << ptMin_ << ") and PVcomparer (" << track_pt_min << ") [PVcomparer considers tracks w/ lower threshold than PixelVertexProducer does] !!!";
45  else
46  edm::LogInfo("PixelVertexProducer") << "minimum track pT setting differs between PixelVertexProducer (" << ptMin_ << ") and PVcomparer (" << track_pt_min << ") !!!";
47  }
48  track_pt_max = PVcomparerPSet.getParameter<double>("track_pt_max");
49  track_chi2_max = PVcomparerPSet.getParameter<double>("track_chi2_max");
50  track_prob_min = PVcomparerPSet.getParameter<double>("track_prob_min");
51  }
52 
53  if (finder == "DivisiveVertexFinder") {
54  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Using the DivisiveVertexFinder\n";
55  dvf_ = new DivisiveVertexFinder(track_pt_min,track_pt_max,track_chi2_max,track_prob_min,zOffset, ntrkMin, useError, zSeparation, wtAverage, verbose_);
56  }
57  else { // Finder not supported, or you made a mistake in your request
58  // throw an exception once I figure out how CMSSW does this
59  }
60 }
61 
62 
64  delete dvf_;
65 }
66 
68 
69  // First fish the pixel tracks out of the event
71  e.getByToken(token_Tracks,trackCollection);
72  const reco::TrackCollection tracks = *(trackCollection.product());
73  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Found " << tracks.size() << " tracks in TrackCollection called " << trackCollName << "\n";
74 
75 
76  // Second, make a collection of pointers to the tracks we want for the vertex finder
78  for (unsigned int i=0; i<tracks.size(); i++) {
79  if (tracks[i].pt() > ptMin_)
80  trks.push_back( reco::TrackRef(trackCollection, i) );
81  }
82  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Selected " << trks.size() << " of these tracks for vertexing\n";
83 
85  e.getByToken(token_BeamSpot,bsHandle);
86  math::XYZPoint myPoint(0.,0.,0.);
87  if (bsHandle.isValid()) myPoint = math::XYZPoint(bsHandle->x0(),bsHandle->y0(), 0. ); //FIXME: fix last coordinate with vertex.z() at same time
88 
89  // Third, ship these tracks off to be vertexed
90  auto vertexes = std::make_unique<reco::VertexCollection>();
91  bool ok;
92  if (method2) {
93  ok = dvf_->findVertexesAlt(trks, // input
94  *vertexes,myPoint); // output
95  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method2 returned status of " << ok;
96  }
97  else {
98  ok = dvf_->findVertexes(trks, // input
99  *vertexes); // output
100  if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method1 returned status of " << ok;
101  }
102 
103  if (verbose_ > 0) {
104  edm::LogInfo("PixelVertexProducer") << ": Found " << vertexes->size() << " vertexes\n";
105  for (unsigned int i=0; i<vertexes->size(); ++i) {
106  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) );
107  }
108  }
109 
110 
111  if(bsHandle.isValid())
112  {
113  const reco::BeamSpot & bs = *bsHandle;
114 
115  for (unsigned int i=0; i<vertexes->size(); ++i) {
116  double z=(*vertexes)[i].z();
117  double x=bs.x0()+bs.dxdz()*(z-bs.z0());
118  double y=bs.y0()+bs.dydz()*(z-bs.z0());
119  reco::Vertex v( reco::Vertex::Point(x,y,z), (*vertexes)[i].error(),(*vertexes)[i].chi2() , (*vertexes)[i].ndof() , (*vertexes)[i].tracksSize());
120  //Copy also the tracks
121  for (std::vector<reco::TrackBaseRef >::const_iterator it = (*vertexes)[i].tracks_begin();
122  it !=(*vertexes)[i].tracks_end(); it++) {
123  v.add( *it );
124  }
125  (*vertexes)[i]=v;
126 
127  }
128  }
129  else
130  {
131  edm::LogWarning("PixelVertexProducer") << "No beamspot found. Using returning vertexes with (0,0,Z) ";
132  }
133 
134  if(!vertexes->size() && bsHandle.isValid()){
135 
136  const reco::BeamSpot & bs = *bsHandle;
137 
139  if ( (bse.cxx() <= 0.) ||
140  (bse.cyy() <= 0.) ||
141  (bse.czz() <= 0.) ) {
143  we(0,0)=10000;
144  we(1,1)=10000;
145  we(2,2)=10000;
146  vertexes->push_back(reco::Vertex(bs.position(), we,0.,0.,0));
147 
148  edm::LogInfo("PixelVertexProducer") <<"No vertices found. Beamspot with invalid errors " << bse.matrix() << std::endl
149  << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n"
150  << (*vertexes)[0].x() << "\n"
151  << (*vertexes)[0].y() << "\n"
152  << (*vertexes)[0].z() << "\n";
153  } else {
154  vertexes->push_back(reco::Vertex(bs.position(),
155  bs.rotatedCovariance3D(),0.,0.,0));
156 
157  edm::LogInfo("PixelVertexProducer") << "No vertices found. Will put Vertex derived from BeamSpot into Event:\n"
158  << (*vertexes)[0].x() << "\n"
159  << (*vertexes)[0].y() << "\n"
160  << (*vertexes)[0].z() << "\n";
161  }
162  }
163 
164  else if(!vertexes->size() && !bsHandle.isValid())
165  {
166  edm::LogWarning("PixelVertexProducer") << "No beamspot and no vertex found. No vertex returned.";
167  }
168 
169  e.put(std::move(vertexes));
170 
171 }
172 
173 
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
const edm::EDGetTokenT< reco::TrackCollection > token_Tracks
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
PixelVertexProducer(const edm::ParameterSet &)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool findVertexesAlt(const reco::TrackRefVector &trks, reco::VertexCollection &vertexes, const math::XYZPoint &bs)
bool exists(std::string const &parameterName) const
checks if a parameter exists
const edm::InputTag trackCollName
const edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot
double dydz() const
dydz slope
Definition: BeamSpot.h:84
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
T sqrt(T t)
Definition: SSEVec.h:18
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
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:74
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
T const * product() const
Definition: Handle.h:81
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
DivisiveVertexFinder * dvf_
fixed size matrix
HLT enums.
double y0() const
y coordinate
Definition: BeamSpot.h:66
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
const Point & position() const
position
Definition: BeamSpot.h:62
virtual void produce(edm::Event &, const edm::EventSetup &) override
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:78
def move(src, dest)
Definition: eostools.py:510
double x0() const
x coordinate
Definition: BeamSpot.h:64