CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Public Attributes | Private Attributes
PrimaryVertexProducerAlgorithm Class Reference

#include <RecoVertex/PrimaryVertexProducerAlgorithm/src/PrimaryVertexProducerAlgorithm.cc>

Inheritance diagram for PrimaryVertexProducerAlgorithm:
VertexReconstructor

Classes

struct  algo
 

Public Member Functions

PrimaryVertexProducerAlgorithmclone () const override
 
edm::ParameterSet config () const
 
 PrimaryVertexProducerAlgorithm (const edm::ParameterSet &)
 
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &tracks) const override
 
virtual std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &beamSpot, const std::string &label="") const
 
 ~PrimaryVertexProducerAlgorithm () override
 
- Public Member Functions inherited from VertexReconstructor
 VertexReconstructor ()
 
virtual std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &t, const reco::BeamSpot &) const
 
virtual std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &primaries, const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &spot) const
 
virtual ~VertexReconstructor ()
 

Public Attributes

edm::InputTag beamSpotLabel
 
edm::InputTag trackLabel
 

Private Attributes

std::vector< algoalgorithms
 
bool fVerbose
 
edm::ParameterSet theConfig
 
TrackClusterizerInZtheTrackClusterizer
 
TrackFilterForPVFindingBasetheTrackFilter
 

Detailed Description

Description: allow redoing the primary vertex reconstruction from a list of tracks, considered obsolete

Implementation: <Notes on="" implementation>="">

Definition at line 57 of file PrimaryVertexProducerAlgorithm.h.

Constructor & Destructor Documentation

PrimaryVertexProducerAlgorithm::PrimaryVertexProducerAlgorithm ( const edm::ParameterSet conf)
explicit

Definition at line 19 of file PrimaryVertexProducerAlgorithm.cc.

References qcdUeDQM_cfi::algorithm, algorithms, beamSpotLabel, edm::ParameterSet::exists(), PrimaryVertexProducerAlgorithm::algo::fitter, fVerbose, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), PrimaryVertexProducerAlgorithm::algo::label, PrimaryVertexProducerAlgorithm::algo::minNdof, AlCaHLTBitMon_QueryRunRegistry::string, theTrackClusterizer, theTrackFilter, trackLabel, PrimaryVertexProducerAlgorithm::algo::useBeamConstraint, HLT_2018_cff::vertexCollections, and PrimaryVertexProducerAlgorithm::algo::vertexSelector.

Referenced by clone().

19  : theConfig(conf) {
20  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
21  trackLabel = conf.getParameter<edm::InputTag>("TrackLabel");
22  beamSpotLabel = conf.getParameter<edm::InputTag>("beamSpotLabel");
23 
24  // select and configure the track selection
25  std::string trackSelectionAlgorithm =
26  conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
27  if (trackSelectionAlgorithm == "filter") {
28  theTrackFilter = new TrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
29  } else if (trackSelectionAlgorithm == "filterWithThreshold") {
31  } else {
32  throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " +
33  trackSelectionAlgorithm);
34  }
35 
36  // select and configure the track clusterizer
37  std::string clusteringAlgorithm =
38  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
39  if (clusteringAlgorithm == "gap") {
41  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
42  } else if (clusteringAlgorithm == "DA") {
44  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
45  }
46  // provide the vectorized version of the clusterizer, if supported by the build
47  else if (clusteringAlgorithm == "DA_vect") {
49  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
50  }
51 
52  else {
53  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
54  }
55 
56  // select and configure the vertex fitters
57  if (conf.exists("vertexCollections")) {
58  std::vector<edm::ParameterSet> vertexCollections =
59  conf.getParameter<std::vector<edm::ParameterSet> >("vertexCollections");
60 
61  for (std::vector<edm::ParameterSet>::const_iterator algoconf = vertexCollections.begin();
62  algoconf != vertexCollections.end();
63  algoconf++) {
65  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
66  if (fitterAlgorithm == "KalmanVertexFitter") {
67  algorithm.fitter = new KalmanVertexFitter();
68  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
69  algorithm.fitter = new AdaptiveVertexFitter();
70  } else {
71  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
72  }
73  algorithm.label = algoconf->getParameter<std::string>("label");
74  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
75  algorithm.useBeamConstraint = algoconf->getParameter<bool>("useBeamConstraint");
76  algorithm.vertexSelector =
77  new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
78  algorithms.push_back(algorithm);
79  }
80  } else {
81  edm::LogWarning("MisConfiguration")
82  << "this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
83 
85  std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
86  if (fitterAlgorithm == "KalmanVertexFitter") {
87  algorithm.fitter = new KalmanVertexFitter();
88  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
89  algorithm.fitter = new AdaptiveVertexFitter();
90  } else {
91  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
92  }
93  algorithm.label = "";
94  algorithm.minNdof = conf.getParameter<double>("minNdof");
95  algorithm.useBeamConstraint = conf.getParameter<bool>("useBeamConstraint");
96 
97  algorithm.vertexSelector = new VertexCompatibleWithBeam(
99  conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
100 
101  algorithms.push_back(algorithm);
102  }
103 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
Common base class.
bool exists(std::string const &parameterName) const
checks if a parameter exists
TrackFilterForPVFindingBase * theTrackFilter
PrimaryVertexProducerAlgorithm::~PrimaryVertexProducerAlgorithm ( )
override

Definition at line 105 of file PrimaryVertexProducerAlgorithm.cc.

References qcdUeDQM_cfi::algorithm, algorithms, theTrackClusterizer, and theTrackFilter.

105  {
106  if (theTrackFilter)
107  delete theTrackFilter;
109  delete theTrackClusterizer;
110  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
111  if (algorithm->fitter)
112  delete algorithm->fitter;
113  if (algorithm->vertexSelector)
114  delete algorithm->vertexSelector;
115  }
116 }
TrackFilterForPVFindingBase * theTrackFilter

Member Function Documentation

PrimaryVertexProducerAlgorithm* PrimaryVertexProducerAlgorithm::clone ( void  ) const
inlineoverridevirtual

Clone method

Implements VertexReconstructor.

Definition at line 70 of file PrimaryVertexProducerAlgorithm.h.

References PrimaryVertexProducerAlgorithm().

70 { return new PrimaryVertexProducerAlgorithm(*this); }
PrimaryVertexProducerAlgorithm(const edm::ParameterSet &)
edm::ParameterSet PrimaryVertexProducerAlgorithm::config ( void  ) const
inline

Definition at line 73 of file PrimaryVertexProducerAlgorithm.h.

References theConfig.

73 { return theConfig; }
std::vector< TransientVertex > PrimaryVertexProducerAlgorithm::vertices ( const std::vector< reco::TransientTrack > &  ) const
overridevirtual

Reconstruct vertices

Implements VertexReconstructor.

Definition at line 123 of file PrimaryVertexProducerAlgorithm.cc.

124  {
125  throw VertexException("PrimaryVertexProducerAlgorithm: cannot make a Primary Vertex without a beam spot");
126 
127  return std::vector<TransientVertex>();
128 }
Common base class.
std::vector< TransientVertex > PrimaryVertexProducerAlgorithm::vertices ( const std::vector< reco::TransientTrack > &  tracks,
const reco::BeamSpot beamSpot,
const std::string &  label = "" 
) const
virtual

Definition at line 130 of file PrimaryVertexProducerAlgorithm.cc.

References qcdUeDQM_cfi::algorithm, algorithms, TrackClusterizerInZ::clusterize(), bsc_activity_cfg::clusters, gather_cfg::cout, GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::cyy(), GlobalErrorBase< T, ErrorWeightType >::czz(), TransientVertex::degreesOfFreedom(), VertexState::error(), fVerbose, TransientVertex::isValid(), GlobalErrorBase< T, ErrorWeightType >::matrix(), TransientVertex::position(), FSQDQM_cfi::pvs, TrackFilterForPVFindingBase::select(), theTrackClusterizer, theTrackFilter, findQualityFiles::v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

132  {
133  bool validBS = true;
134  VertexState beamVertexState(beamSpot);
135  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
136  (beamVertexState.error().czz() <= 0.)) {
137  validBS = false;
138  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
139  }
140 
141  // // get RECO tracks from the event
142  // // `tks` can be used as a ptr to a reco::TrackCollection
143  // edm::Handle<reco::TrackCollection> tks;
144  // iEvent.getByLabel(trackLabel, tks);
145 
146  // select tracks
147  std::vector<reco::TransientTrack> seltks = theTrackFilter->select(t_tks);
148 
149  // clusterize tracks in Z
150  std::vector<std::vector<reco::TransientTrack> > clusters = theTrackClusterizer->clusterize(seltks);
151  if (fVerbose) {
152  std::cout << " clustering returned " << clusters.size() << " clusters from " << seltks.size()
153  << " selected tracks" << std::endl;
154  }
155 
156  // vertex fits
157  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
158  if (!(algorithm->label == label))
159  continue;
160 
161  //std::auto_ptr<reco::VertexCollection> result(new reco::VertexCollection);
162  // reco::VertexCollection vColl;
163 
164  std::vector<TransientVertex> pvs;
165  for (std::vector<std::vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin();
166  iclus != clusters.end();
167  iclus++) {
169  if (algorithm->useBeamConstraint && validBS && ((*iclus).size() > 1)) {
170  v = algorithm->fitter->vertex(*iclus, beamSpot);
171 
172  } else if (!(algorithm->useBeamConstraint) && ((*iclus).size() > 1)) {
173  v = algorithm->fitter->vertex(*iclus);
174 
175  } // else: no fit ==> v.isValid()=False
176 
177  if (fVerbose) {
178  if (v.isValid())
179  std::cout << "x,y,z=" << v.position().x() << " " << v.position().y() << " " << v.position().z() << std::endl;
180  else
181  std::cout << "Invalid fitted vertex\n";
182  }
183 
184  if (v.isValid() && (v.degreesOfFreedom() >= algorithm->minNdof) &&
185  (!validBS || (*(algorithm->vertexSelector))(v, beamVertexState)))
186  pvs.push_back(v);
187  } // end of cluster loop
188 
189  if (fVerbose) {
190  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
191  }
192 
193  // sort vertices by pt**2 vertex (aka signal vertex tagging)
194  if (pvs.size() > 1) {
195  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
196  }
197 
198  return pvs;
199  }
200 
201  std::vector<TransientVertex> dummy;
202  return dummy; //avoid compiler warning, should never be here
203 }
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
T y() const
Definition: PV3DBase.h:60
char const * label
float degreesOfFreedom() const
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:61
TrackFilterForPVFindingBase * theTrackFilter
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
T x() const
Definition: PV3DBase.h:59
bool isValid() const

Member Data Documentation

std::vector<algo> PrimaryVertexProducerAlgorithm::algorithms
private
edm::InputTag PrimaryVertexProducerAlgorithm::beamSpotLabel

Definition at line 75 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm().

bool PrimaryVertexProducerAlgorithm::fVerbose
private

Definition at line 95 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm(), and vertices().

edm::ParameterSet PrimaryVertexProducerAlgorithm::theConfig
private

Definition at line 94 of file PrimaryVertexProducerAlgorithm.h.

Referenced by config().

TrackClusterizerInZ* PrimaryVertexProducerAlgorithm::theTrackClusterizer
private
TrackFilterForPVFindingBase* PrimaryVertexProducerAlgorithm::theTrackFilter
private
edm::InputTag PrimaryVertexProducerAlgorithm::trackLabel

Definition at line 74 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm().