CMS 3D CMS Logo

PrimaryVertexProducerAlgorithm.cc
Go to the documentation of this file.
3 
9 
13 
18 
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 }
104 
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 }
117 
118 //
119 // member functions
120 //
121 
122 // obsolete method, unfortunately required through inheritance from VertexReconstructor
123 std::vector<TransientVertex> PrimaryVertexProducerAlgorithm::vertices(
124  const std::vector<reco::TransientTrack>& tracks) const {
125  throw VertexException("PrimaryVertexProducerAlgorithm: cannot make a Primary Vertex without a beam spot");
126 
127  return std::vector<TransientVertex>();
128 }
129 
130 std::vector<TransientVertex> PrimaryVertexProducerAlgorithm::vertices(const std::vector<reco::TransientTrack>& t_tks,
131  const reco::BeamSpot& beamSpot,
132  const std::string& label) const {
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 }
AdaptiveVertexFitter
Definition: AdaptiveVertexFitter.h:29
TrackFilterForPVFindingBase::select
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
Handle.h
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
pwdgSkimBPark_cfi.beamSpot
beamSpot
Definition: pwdgSkimBPark_cfi.py:5
MessageLogger.h
ESHandle.h
VertexException
Common base class.
Definition: VertexException.h:12
gather_cfg.cout
cout
Definition: gather_cfg.py:144
GlobalErrorBase::matrix
const AlgebraicSymMatrix33 matrix() const
Definition: GlobalErrorBase.h:121
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
TransientTrack.h
findQualityFiles.v
v
Definition: findQualityFiles.py:179
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
PrimaryVertexProducerAlgorithm::~PrimaryVertexProducerAlgorithm
~PrimaryVertexProducerAlgorithm() override
Definition: PrimaryVertexProducerAlgorithm.cc:105
PrimaryVertexProducerAlgorithm::trackLabel
edm::InputTag trackLabel
Definition: PrimaryVertexProducerAlgorithm.h:74
MakerMacros.h
GlobalErrorBase::cyy
T cyy() const
Definition: GlobalErrorBase.h:101
TrackFwd.h
VertexDistanceXY.h
BeamSpot.h
reco::BeamSpot
Definition: BeamSpot.h:21
VertexState::error
GlobalError error() const
Definition: VertexState.h:64
GlobalErrorBase::cxx
T cxx() const
Definition: GlobalErrorBase.h:97
PrimaryVertexProducerAlgorithm::theTrackFilter
TrackFilterForPVFindingBase * theTrackFilter
Definition: PrimaryVertexProducerAlgorithm.h:80
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
edm::ParameterSet::exists
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: ParameterSet.cc:681
TransientTrackBuilder.h
edm::ParameterSet
Definition: ParameterSet.h:47
TrackClusterizerInZ::clusterize
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
qcdUeDQM_cfi.algorithm
algorithm
Definition: qcdUeDQM_cfi.py:32
PrimaryVertexProducerAlgorithm::algo
Definition: PrimaryVertexProducerAlgorithm.h:84
VertexCompatibleWithBeam
Definition: VertexCompatibleWithBeam.h:15
DAClusterizerInZ
Definition: DAClusterizerInZ.h:18
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
TransientVertex
Definition: TransientVertex.h:18
VertexDistanceXY
Definition: VertexDistanceXY.h:11
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
TransientTrackRecord.h
FSQDQM_cfi.pvs
pvs
Definition: FSQDQM_cfi.py:12
VertexFwd.h
TransientVertex.h
DAClusterizerInZ_vect
Definition: DAClusterizerInZ_vect.h:20
PrimaryVertexProducerAlgorithm::algorithms
std::vector< algo > algorithms
Definition: PrimaryVertexProducerAlgorithm.h:92
HITrackFilterForPVFinding
Definition: HITrackFilterForPVFinding.h:13
PrimaryVertexProducerAlgorithm::PrimaryVertexProducerAlgorithm
PrimaryVertexProducerAlgorithm(const edm::ParameterSet &)
Definition: PrimaryVertexProducerAlgorithm.cc:19
VertexState
Definition: VertexState.h:13
HLT_FULL_cff.vertexCollections
vertexCollections
Definition: HLT_FULL_cff.py:34798
PrimaryVertexProducerAlgorithm::theTrackClusterizer
TrackClusterizerInZ * theTrackClusterizer
Definition: PrimaryVertexProducerAlgorithm.h:81
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
GapClusterizerInZ
Definition: GapClusterizerInZ.h:14
GlobalErrorBase::czz
T czz() const
Definition: GlobalErrorBase.h:107
TrackFilterForPVFinding
Definition: TrackFilterForPVFinding.h:14
VertexHigherPtSquared
Definition: VertexHigherPtSquared.h:13
PrimaryVertexProducerAlgorithm::fVerbose
bool fVerbose
Definition: PrimaryVertexProducerAlgorithm.h:95
dummy
Definition: DummySelector.h:38
PrimaryVertexProducerAlgorithm.h
PrimaryVertexProducerAlgorithm::vertices
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
Definition: PrimaryVertexProducerAlgorithm.cc:123
PrimaryVertexProducerAlgorithm::beamSpotLabel
edm::InputTag beamSpotLabel
Definition: PrimaryVertexProducerAlgorithm.h:75
edm::InputTag
Definition: InputTag.h:15
label
const char * label
Definition: PFTauDecayModeTools.cc:11
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22