CMS 3D CMS Logo

PrimaryVertexProducer.cc
Go to the documentation of this file.
2 
9 
13 
18 
20 
22  :theConfig(conf)
23 {
24 
25  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
26 
27  trkToken = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackLabel"));
28  bsToken = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotLabel"));
29  f4D = false;
30 
31  // select and configure the track selection
32  std::string trackSelectionAlgorithm=conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
33  if(trackSelectionAlgorithm=="filter"){
34  theTrackFilter= new TrackFilterForPVFinding( conf.getParameter<edm::ParameterSet>("TkFilterParameters") );
35  }else if (trackSelectionAlgorithm=="filterWithThreshold"){
37  }else{
38  throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " + trackSelectionAlgorithm);
39  }
40 
41 
42  // select and configure the track clusterizer
43  std::string clusteringAlgorithm=conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
44  if (clusteringAlgorithm=="gap"){
45  theTrackClusterizer = new GapClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
46  }else if(clusteringAlgorithm=="DA"){
47  theTrackClusterizer = new DAClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
48  }
49  // provide the vectorized version of the clusterizer, if supported by the build
50  else if(clusteringAlgorithm == "DA_vect") {
51  theTrackClusterizer = new DAClusterizerInZ_vect(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
52  }
53  else if( clusteringAlgorithm=="DA2D" ) {
54  theTrackClusterizer = new DAClusterizerInZT(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
55  f4D = true;
56  }
57 
58  else{
59  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
60  }
61 
62  if( f4D ) {
63  trkTimesToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimesLabel") );
64  trkTimeResosToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimeResosLabel") );
65  }
66 
67 
68  // select and configure the vertex fitters
69  if (conf.exists("vertexCollections")){
70  std::vector<edm::ParameterSet> vertexCollections =conf.getParameter< std::vector<edm::ParameterSet> >("vertexCollections");
71 
72  for( std::vector< edm::ParameterSet >::const_iterator algoconf = vertexCollections.begin(); algoconf != vertexCollections.end(); algoconf++){
73 
75  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
76  if (fitterAlgorithm=="KalmanVertexFitter") {
77  algorithm.fitter= new KalmanVertexFitter();
78  } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
79  algorithm.fitter= new AdaptiveVertexFitter( GeometricAnnealing( algoconf->getParameter<double>("chi2cutoff")));
80  } else {
81  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
82  }
83  algorithm.label = algoconf->getParameter<std::string>("label");
84  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
85  algorithm.useBeamConstraint=algoconf->getParameter<bool>("useBeamConstraint");
86  algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
87  algorithms.push_back(algorithm);
88 
89  produces<reco::VertexCollection>(algorithm.label);
90  }
91  }else{
92  edm::LogWarning("MisConfiguration")<<"this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
93 
95  std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
96  if (fitterAlgorithm=="KalmanVertexFitter") {
97  algorithm.fitter= new KalmanVertexFitter();
98  } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
99  algorithm.fitter= new AdaptiveVertexFitter();
100  } else {
101  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
102  }
103  algorithm.label = "";
104  algorithm.minNdof = conf.getParameter<double>("minNdof");
105  algorithm.useBeamConstraint=conf.getParameter<bool>("useBeamConstraint");
106 
107  algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
108 
109  algorithms.push_back(algorithm);
110  produces<reco::VertexCollection>(algorithm.label);
111  }
112 
113 
114 }
115 
116 
118  if (theTrackFilter) delete theTrackFilter;
120  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
121  if (algorithm->fitter) delete algorithm->fitter;
122  if (algorithm->vertexSelector) delete algorithm->vertexSelector;
123  }
124 }
125 
126 
127 void
129 {
130 
131  // get the BeamSpot, it will alwys be needed, even when not used as a constraint
133  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
134  iEvent.getByToken(bsToken,recoBeamSpotHandle);
135  if (recoBeamSpotHandle.isValid()){
136  beamSpot = *recoBeamSpotHandle;
137  }else{
138  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
139  }
140 
141  bool validBS = true;
142  VertexState beamVertexState(beamSpot);
143  if ( (beamVertexState.error().cxx() <= 0.) ||
144  (beamVertexState.error().cyy() <= 0.) ||
145  (beamVertexState.error().czz() <= 0.) ) {
146  validBS = false;
147  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors "<<beamVertexState.error().matrix();
148  }
149 
150 
151  // get RECO tracks from the event
152  // `tks` can be used as a ptr to a reco::TrackCollection
154  iEvent.getByToken(trkToken, tks);
155 
156 
157  // interface RECO tracks to vertex reconstruction
159  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
160  std::vector<reco::TransientTrack> t_tks;
161 
162  if( f4D ) {
163  edm::Handle<edm::ValueMap<float> > trackTimesH;
164  edm::Handle<edm::ValueMap<float> > trackTimeResosH;
165  iEvent.getByToken(trkTimesToken, trackTimesH);
166  iEvent.getByToken(trkTimeResosToken, trackTimeResosH);
167  t_tks = (*theB).build(tks, beamSpot, *(trackTimesH.product()), *(trackTimeResosH.product()));
168  } else {
169  t_tks = (*theB).build(tks, beamSpot);
170  }
171  if(fVerbose) {std::cout << "RecoVertex/PrimaryVertexProducer"
172  << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
173  }
174 
175 
176  // select tracks
177  std::vector<reco::TransientTrack> && seltks = theTrackFilter->select( t_tks );
178 
179  // clusterize tracks in Z
180  std::vector< std::vector<reco::TransientTrack> > && clusters = theTrackClusterizer->clusterize(seltks);
181 
182  if (fVerbose){std::cout << " clustering returned "<< clusters.size() << " clusters from " << seltks.size() << " selected tracks" <<std::endl;}
183 
184 
185  // vertex fits
186  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
187 
188 
189  auto result = std::make_unique<reco::VertexCollection>();
190  reco::VertexCollection & vColl = (*result);
191 
192 
193  std::vector<TransientVertex> pvs;
194  for (std::vector< std::vector<reco::TransientTrack> >::const_iterator iclus
195  = clusters.begin(); iclus != clusters.end(); iclus++) {
196 
197  double meantime = 0.;
198  double expv_x2 = 0.;
199  double normw = 0.;
200  if( f4D ) {
201  for( const auto& tk : *iclus ) {
202  const double time = tk.timeExt();
203  const double inverr = 1.0/tk.dtErrorExt();
204  meantime += time*inverr;
205  expv_x2 += time*time*inverr;
206  normw += inverr;
207  }
208  meantime = meantime/normw;
209  expv_x2 = expv_x2/normw;
210  }
211  const double time_var = ( f4D ? expv_x2 - meantime*meantime : 0. );
212 
213 
215  if( algorithm->useBeamConstraint && validBS &&((*iclus).size()>1) ){
216 
217  v = algorithm->fitter->vertex(*iclus, beamSpot);
218 
219  if( f4D ) {
220  if( v.isValid() ) {
221  auto err = v.positionError().matrix4D();
222  err(3,3) = time_var/(double)iclus->size();
223  v = TransientVertex(v.position(),meantime,err,v.originalTracks(),v.totalChiSquared());
224  }
225  }
226 
227  }else if( !(algorithm->useBeamConstraint) && ((*iclus).size()>1) ) {
228 
229  v = algorithm->fitter->vertex(*iclus);
230 
231  if( f4D ) {
232  if( v.isValid() ) {
233  auto err = v.positionError().matrix4D();
234  err(3,3) = time_var/(double)iclus->size();
235  v = TransientVertex(v.position(),meantime,err,v.originalTracks(),v.totalChiSquared());
236  }
237  }
238 
239  }// else: no fit ==> v.isValid()=False
240 
241 
242  if (fVerbose){
243  if (v.isValid()) {
244  std::cout << "x,y,z";
245  if (f4D) std::cout << ",t";
246  std::cout << "=" << v.position().x() <<" " << v.position().y() << " " << v.position().z();
247  if (f4D) std::cout << " " << v.time();
248  std::cout << " cluster size = " << (*iclus).size() << std::endl;
249  }
250  else{
251  std::cout <<"Invalid fitted vertex, cluster size=" << (*iclus).size() << std::endl;
252  }
253  }
254 
255  if ( v.isValid()
256  && (v.degreesOfFreedom()>=algorithm->minNdof)
257  && (!validBS || (*(algorithm->vertexSelector))(v,beamVertexState))
258  ) pvs.push_back(v);
259  }// end of cluster loop
260 
261  if(fVerbose){
262  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
263  }
264 
265 
266  if (clusters.size()>2 && clusters.size() > 2*pvs.size())
267  edm::LogWarning("PrimaryVertexProducer") << "more than half of candidate vertices lost " << pvs.size() << ' ' << clusters.size();
268 
269  if (pvs.empty() && seltks.size()>5)
270  edm::LogWarning("PrimaryVertexProducer") << "no vertex found with " << seltks.size() << " tracks and " << clusters.size() <<" vertex-candidates";
271 
272  // sort vertices by pt**2 vertex (aka signal vertex tagging)
273  if(pvs.size()>1){
274  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
275  }
276 
277 
278 
279  // convert transient vertices returned by the theAlgo to (reco) vertices
280  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin();
281  iv != pvs.end(); iv++) {
282  reco::Vertex v = *iv;
283  vColl.push_back(v);
284  }
285 
286  if (vColl.empty()) {
287  GlobalError bse(beamSpot.rotatedCovariance3D());
288  if ( (bse.cxx() <= 0.) ||
289  (bse.cyy() <= 0.) ||
290  (bse.czz() <= 0.) ) {
292  we(0,0)=10000; we(1,1)=10000; we(2,2)=10000;
293  vColl.push_back(reco::Vertex(beamSpot.position(), we,0.,0.,0));
294  if(fVerbose){
295  std::cout <<"RecoVertex/PrimaryVertexProducer: "
296  << "Beamspot with invalid errors "<<bse.matrix()<<std::endl;
297  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
298  }
299  } else {
300  vColl.push_back(reco::Vertex(beamSpot.position(),
301  beamSpot.rotatedCovariance3D(),0.,0.,0));
302  if(fVerbose){
303  std::cout <<"RecoVertex/PrimaryVertexProducer: "
304  << " will put Vertex derived from BeamSpot into Event.\n";
305  }
306  }
307  }
308 
309  if(fVerbose){
310  int ivtx=0;
311  for(reco::VertexCollection::const_iterator v=vColl.begin();
312  v!=vColl.end(); ++v){
313  std::cout << "recvtx "<< ivtx++
314  << "#trk " << std::setw(3) << v->tracksSize()
315  << " chi2 " << std::setw(4) << v->chi2()
316  << " ndof " << std::setw(3) << v->ndof()
317  << " x " << std::setw(6) << v->position().x()
318  << " dx " << std::setw(6) << v->xError()
319  << " y " << std::setw(6) << v->position().y()
320  << " dy " << std::setw(6) << v->yError()
321  << " z " << std::setw(6) << v->position().z()
322  << " dz " << std::setw(6) << v->zError();
323  if( f4D ) {
324  std::cout << " t " << std::setw(6) << v->t()
325  << " dt " << std::setw(6) << v->tError();
326  }
327  std::cout << std::endl;
328  }
329  }
330 
331  iEvent.put(std::move(result), algorithm->label);
332  }
333 
334 }
335 
336 
337 //define this as a plug-in
GlobalError positionError() const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
vertexCollections
1 cm max separation between clusters
edm::EDGetTokenT< reco::BeamSpot > bsToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
Common base class.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const AlgebraicSymMatrix33 matrix() const
float totalChiSquared() const
edm::EDGetTokenT< reco::TrackCollection > trkToken
TrackClusterizerInZ * theTrackClusterizer
T y() const
Definition: PV3DBase.h:63
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
const AlgebraicSymMatrix44 & matrix4D() const
int iEvent
Definition: GenABIO.cc:230
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
std::vector< reco::TransientTrack > const & originalTracks() const
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
float degreesOfFreedom() const
VertexCompatibleWithBeam * vertexSelector
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:64
TrackFilterForPVFindingBase * theTrackFilter
bool isValid() const
Definition: HandleBase.h:74
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
virtual void produce(edm::Event &, const edm::EventSetup &)
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:56
std::vector< algo > algorithms
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
double time() const
PrimaryVertexProducer(const edm::ParameterSet &)
GlobalError error() const
Definition: VertexState.h:74
const Point & position() const
position
Definition: BeamSpot.h:62
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:78
T x() const
Definition: PV3DBase.h:62
bool isValid() const
def move(src, dest)
Definition: eostools.py:510