CMS 3D CMS Logo

PixelVertexProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PixelVertexProducer
4 // Class: PixelVertexProducer
5 //
18 //
19 // Original Author: Aaron Dominguez (UNL)
20 // Created: Thu May 25 10:17:32 CDT 2006
21 //
22 
37 #include <memory>
38 #include <string>
39 #include <cmath>
40 
42 public:
43  explicit PixelVertexProducer(const edm::ParameterSet&);
44  ~PixelVertexProducer() override;
45 
46  void produce(edm::Event&, const edm::EventSetup&) override;
47 
48 private:
49  // ----------member data ---------------------------
50  // Turn on debug printing if verbose_ > 0
51  const int verbose_;
52  // Tracking cuts before sending tracks to vertex algo
53  const double ptMin_;
54  const bool method2;
58 
60 };
61 
63  // 0 silent, 1 chatty, 2 loud
64  : verbose_(conf.getParameter<int>("Verbosity")),
65  // 1.0 GeV
66  ptMin_(conf.getParameter<double>("PtMin")),
67  method2(conf.getParameter<bool>("Method2")),
68  trackCollName(conf.getParameter<edm::InputTag>("TrackCollection")),
69  token_Tracks(consumes<reco::TrackCollection>(trackCollName)),
70  token_BeamSpot(consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpot"))) {
71  // Register my product
72  produces<reco::VertexCollection>();
73 
74  // Setup shop
75  std::string finder = conf.getParameter<std::string>("Finder"); // DivisiveVertexFinder
76  bool useError = conf.getParameter<bool>("UseError"); // true
77  bool wtAverage = conf.getParameter<bool>("WtAverage"); // true
78  double zOffset = conf.getParameter<double>("ZOffset"); // 5.0 sigma
79  double zSeparation = conf.getParameter<double>("ZSeparation"); // 0.05 cm
80  int ntrkMin = conf.getParameter<int>("NTrkMin"); // 3
81  // Tracking requirements before sending a track to be considered for vtx
82 
83  double track_pt_min = ptMin_;
84  double track_pt_max = 10.;
85  double track_chi2_max = 9999999.;
86  double track_prob_min = -1.;
87 
88  if (conf.exists("PVcomparer")) {
89  edm::ParameterSet PVcomparerPSet = conf.getParameter<edm::ParameterSet>("PVcomparer");
90  track_pt_min = PVcomparerPSet.getParameter<double>("track_pt_min");
91  if (track_pt_min != ptMin_) {
92  if (track_pt_min < ptMin_)
93  edm::LogInfo("PixelVertexProducer")
94  << "minimum track pT setting differs between PixelVertexProducer (" << ptMin_ << ") and PVcomparer ("
95  << track_pt_min << ") [PVcomparer considers tracks w/ lower threshold than PixelVertexProducer does] !!!";
96  else
97  edm::LogInfo("PixelVertexProducer") << "minimum track pT setting differs between PixelVertexProducer ("
98  << ptMin_ << ") and PVcomparer (" << track_pt_min << ") !!!";
99  }
100  track_pt_max = PVcomparerPSet.getParameter<double>("track_pt_max");
101  track_chi2_max = PVcomparerPSet.getParameter<double>("track_chi2_max");
102  track_prob_min = PVcomparerPSet.getParameter<double>("track_prob_min");
103  }
104 
105  if (finder == "DivisiveVertexFinder") {
106  if (verbose_ > 0)
107  edm::LogInfo("PixelVertexProducer") << ": Using the DivisiveVertexFinder\n";
109  track_pt_max,
112  zOffset,
113  ntrkMin,
114  useError,
115  zSeparation,
116  wtAverage,
117  verbose_);
118  } else { // Finder not supported, or you made a mistake in your request
119  // throw an exception once I figure out how CMSSW does this
120  }
121 }
122 
124 
126  // First fish the pixel tracks out of the event
128  e.getByToken(token_Tracks, trackCollection);
129  const reco::TrackCollection tracks = *(trackCollection.product());
130  if (verbose_ > 0)
131  edm::LogInfo("PixelVertexProducer") << ": Found " << tracks.size() << " tracks in TrackCollection called "
132  << trackCollName << "\n";
133 
134  // Second, make a collection of pointers to the tracks we want for the vertex finder
136  for (unsigned int i = 0; i < tracks.size(); i++) {
137  if (tracks[i].pt() > ptMin_)
139  }
140  if (verbose_ > 0)
141  edm::LogInfo("PixelVertexProducer") << ": Selected " << trks.size() << " of these tracks for vertexing\n";
142 
144  e.getByToken(token_BeamSpot, bsHandle);
145  math::XYZPoint myPoint(0., 0., 0.);
146  if (bsHandle.isValid())
147  //FIXME: fix last coordinate with vertex.z() at same time
148  myPoint = math::XYZPoint(bsHandle->x0(), bsHandle->y0(), 0.);
149 
150  // Third, ship these tracks off to be vertexed
151  auto vertexes = std::make_unique<reco::VertexCollection>();
152  bool ok;
153  if (method2) {
154  ok = dvf_->findVertexesAlt(trks, // input
155  *vertexes,
156  myPoint); // output
157  if (verbose_ > 0)
158  edm::LogInfo("PixelVertexProducer") << "Method2 returned status of " << ok;
159  } else {
160  ok = dvf_->findVertexes(trks, // input
161  *vertexes); // output
162  if (verbose_ > 0)
163  edm::LogInfo("PixelVertexProducer") << "Method1 returned status of " << ok;
164  }
165 
166  if (verbose_ > 0) {
167  edm::LogInfo("PixelVertexProducer") << ": Found " << vertexes->size() << " vertexes\n";
168  for (unsigned int i = 0; i < vertexes->size(); ++i) {
169  edm::LogInfo("PixelVertexProducer")
170  << "Vertex number " << i << " has " << (*vertexes)[i].tracksSize() << " tracks with a position of "
171  << (*vertexes)[i].z() << " +- " << std::sqrt((*vertexes)[i].covariance(2, 2));
172  }
173  }
174 
175  if (bsHandle.isValid()) {
176  const reco::BeamSpot& bs = *bsHandle;
177 
178  for (unsigned int i = 0; i < vertexes->size(); ++i) {
179  double z = (*vertexes)[i].z();
180  double x = bs.x0() + bs.dxdz() * (z - bs.z0());
181  double y = bs.y0() + bs.dydz() * (z - bs.z0());
183  (*vertexes)[i].error(),
184  (*vertexes)[i].chi2(),
185  (*vertexes)[i].ndof(),
186  (*vertexes)[i].tracksSize());
187  //Copy also the tracks
188  for (std::vector<reco::TrackBaseRef>::const_iterator it = (*vertexes)[i].tracks_begin();
189  it != (*vertexes)[i].tracks_end();
190  it++) {
191  v.add(*it);
192  }
193  (*vertexes)[i] = v;
194  }
195  } else {
196  edm::LogWarning("PixelVertexProducer") << "No beamspot found. Using returning vertexes with (0,0,Z) ";
197  }
198 
199  if (vertexes->empty() && bsHandle.isValid()) {
200  const reco::BeamSpot& bs = *bsHandle;
201 
202  GlobalError bse(bs.rotatedCovariance3D());
203  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
205  we(0, 0) = 10000;
206  we(1, 1) = 10000;
207  we(2, 2) = 10000;
208  vertexes->push_back(reco::Vertex(bs.position(), we, 0., 0., 0));
209 
210  edm::LogInfo("PixelVertexProducer")
211  << "No vertices found. Beamspot with invalid errors " << bse.matrix() << std::endl
212  << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n"
213  << (*vertexes)[0].x() << "\n"
214  << (*vertexes)[0].y() << "\n"
215  << (*vertexes)[0].z() << "\n";
216  } else {
217  vertexes->push_back(reco::Vertex(bs.position(), bs.rotatedCovariance3D(), 0., 0., 0));
218 
219  edm::LogInfo("PixelVertexProducer") << "No vertices found. Will put Vertex derived from BeamSpot into Event:\n"
220  << (*vertexes)[0].x() << "\n"
221  << (*vertexes)[0].y() << "\n"
222  << (*vertexes)[0].z() << "\n";
223  }
224  }
225 
226  else if (vertexes->empty() && !bsHandle.isValid()) {
227  edm::LogWarning("PixelVertexProducer") << "No beamspot and no vertex found. No vertex returned.";
228  }
229 
230  e.put(std::move(vertexes));
231 }
232 
DDAxes::y
HLT_FULL_cff.track_pt_min
track_pt_min
Definition: HLT_FULL_cff.py:256
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
AlgebraicSymMatrix33
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Definition: AlgebraicROOTObjects.h:21
electrons_cff.bool
bool
Definition: electrons_cff.py:393
mps_fire.i
i
Definition: mps_fire.py:428
HLT_FULL_cff.track_prob_min
track_prob_min
Definition: HLT_FULL_cff.py:255
MessageLogger.h
PixelVertexProducer::method2
const bool method2
Definition: PixelVertexProducer.cc:54
align::BeamSpot
Definition: StructureType.h:89
HLT_FULL_cff.finder
finder
Definition: HLT_FULL_cff.py:51924
HLT_FULL_cff.zSeparation
zSeparation
Definition: HLT_FULL_cff.py:104291
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
edm::EDGetTokenT< reco::TrackCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
DivisiveVertexFinder::findVertexesAlt
bool findVertexesAlt(const reco::TrackRefVector &trks, reco::VertexCollection &vertexes, const math::XYZPoint &bs)
Definition: DivisiveVertexFinder.cc:60
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
EDProducer.h
jets_cff.vertexes
vertexes
Definition: jets_cff.py:706
DDAxes::x
edm::RefVector< TrackCollection >
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
findQualityFiles.v
v
Definition: findQualityFiles.py:179
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
convertSQLiteXML.ok
bool ok
Definition: convertSQLiteXML.py:98
edm::Handle< reco::TrackCollection >
relativeConstraints.error
error
Definition: relativeConstraints.py:53
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
edm::Ref< TrackCollection >
ndof
Definition: HIMultiTrackSelector.h:49
DivisiveVertexFinder
Definition: DivisiveVertexFinder.h:34
MakerMacros.h
cms::cuda::bs
bs
Definition: HistoContainer.h:127
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
BeamSpot.h
PixelVertexProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PixelVertexProducer.cc:125
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
reco::BeamSpot
Definition: BeamSpot.h:21
DDAxes::z
DivisiveVertexFinder.h
PixelVertexProducer
Definition: PixelVertexProducer.cc:41
HLT_FULL_cff.track_chi2_max
track_chi2_max
Definition: HLT_FULL_cff.py:253
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
edm::ParameterSet::exists
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: ParameterSet.cc:681
HLT_FULL_cff.track_pt_max
track_pt_max
Definition: HLT_FULL_cff.py:254
edm::ParameterSet
Definition: ParameterSet.h:47
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
duplicaterechits_cfi.trackCollection
trackCollection
Definition: duplicaterechits_cfi.py:4
Event.h
PixelVertexProducer::token_BeamSpot
const edm::EDGetTokenT< reco::BeamSpot > token_BeamSpot
Definition: PixelVertexProducer.cc:57
PixelVertexProducer::token_Tracks
const edm::EDGetTokenT< reco::TrackCollection > token_Tracks
Definition: PixelVertexProducer.cc:56
createfilelist.int
int
Definition: createfilelist.py:10
Reconstruction_HI_cff.trackCollName
trackCollName
Definition: Reconstruction_HI_cff.py:34
GlobalError.h
edm::stream::EDProducer
Definition: EDProducer.h:38
GlobalErrorBase< double, ErrorMatrixTag >
edm::EventSetup
Definition: EventSetup.h:57
PixelVertexProducer::verbose_
const int verbose_
Definition: PixelVertexProducer.cc:51
edm::RefVector::push_back
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
InputTag.h
VertexFwd.h
reco::BeamSpot::x0
double x0() const
x coordinate
Definition: BeamSpot.h:61
PixelVertexProducer::trackCollName
const edm::InputTag trackCollName
Definition: PixelVertexProducer.cc:55
reco::Vertex::Point
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
eostools.move
def move(src, dest)
Definition: eostools.py:511
Frameworkfwd.h
DivisiveVertexFinder::findVertexes
bool findVertexes(const reco::TrackRefVector &trks, reco::VertexCollection &vertexes)
Run the divisive algorithm and return a vector of vertexes for the input track collection.
Definition: DivisiveVertexFinder.cc:35
EventSetup.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PixelVertexProducer::PixelVertexProducer
PixelVertexProducer(const edm::ParameterSet &)
Definition: PixelVertexProducer.cc:62
PixelVertexProducer::ptMin_
const double ptMin_
Definition: PixelVertexProducer.cc:53
PixelVertexProducer::~PixelVertexProducer
~PixelVertexProducer() override
Definition: PixelVertexProducer.cc:123
ParameterSet.h
PixelVertexProducer::dvf_
DivisiveVertexFinder * dvf_
Definition: PixelVertexProducer.cc:59
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
reco::BeamSpot::y0
double y0() const
y coordinate
Definition: BeamSpot.h:63
edm::Event
Definition: Event.h:73
edm::RefVector::size
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
edm::InputTag
Definition: InputTag.h:15
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
reco::Vertex
Definition: Vertex.h:35
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37