CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PrimaryVertexProducer.cc
Go to the documentation of this file.
2 
9 
13 
18 
19 
21  :theConfig(conf)
22 {
23 
24  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
25  trkToken = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackLabel"));
26  bsToken = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotLabel"));
27 
28  // select and configure the track selection
29  std::string trackSelectionAlgorithm=conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
30  if(trackSelectionAlgorithm=="filter"){
31  theTrackFilter= new TrackFilterForPVFinding( conf.getParameter<edm::ParameterSet>("TkFilterParameters") );
32  }else if (trackSelectionAlgorithm=="filterWithThreshold"){
34  }else{
35  throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " + trackSelectionAlgorithm);
36  }
37 
38 
39  // select and configure the track clusterizer
40  std::string clusteringAlgorithm=conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
41  if (clusteringAlgorithm=="gap"){
42  theTrackClusterizer = new GapClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
43  }else if(clusteringAlgorithm=="DA"){
44  theTrackClusterizer = new DAClusterizerInZ(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") {
48  theTrackClusterizer = new DAClusterizerInZ_vect(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
49  }
50 
51 
52  else{
53  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
54  }
55 
56 
57  // select and configure the vertex fitters
58  if (conf.exists("vertexCollections")){
59  std::vector<edm::ParameterSet> vertexCollections =conf.getParameter< std::vector<edm::ParameterSet> >("vertexCollections");
60 
61  for( std::vector< edm::ParameterSet >::const_iterator algoconf = vertexCollections.begin(); algoconf != vertexCollections.end(); algoconf++){
62 
64  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
65  if (fitterAlgorithm=="KalmanVertexFitter") {
66  algorithm.fitter= new KalmanVertexFitter();
67  } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
68  algorithm.fitter= new AdaptiveVertexFitter();
69  } else {
70  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
71  }
72  algorithm.label = algoconf->getParameter<std::string>("label");
73  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
74  algorithm.useBeamConstraint=algoconf->getParameter<bool>("useBeamConstraint");
75  algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
76  algorithms.push_back(algorithm);
77 
78  produces<reco::VertexCollection>(algorithm.label);
79  }
80  }else{
81  edm::LogWarning("MisConfiguration")<<"this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
82 
84  std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
85  if (fitterAlgorithm=="KalmanVertexFitter") {
86  algorithm.fitter= new KalmanVertexFitter();
87  } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
88  algorithm.fitter= new AdaptiveVertexFitter();
89  } else {
90  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
91  }
92  algorithm.label = "";
93  algorithm.minNdof = conf.getParameter<double>("minNdof");
94  algorithm.useBeamConstraint=conf.getParameter<bool>("useBeamConstraint");
95 
96  algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
97 
98  algorithms.push_back(algorithm);
99  produces<reco::VertexCollection>(algorithm.label);
100  }
101 
102 
103 }
104 
105 
107  if (theTrackFilter) delete theTrackFilter;
109  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
110  if (algorithm->fitter) delete algorithm->fitter;
111  if (algorithm->vertexSelector) delete algorithm->vertexSelector;
112  }
113 }
114 
115 
116 void
118 {
119 
120  // get the BeamSpot, it will alwys be needed, even when not used as a constraint
122  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
123  iEvent.getByToken(bsToken,recoBeamSpotHandle);
124  if (recoBeamSpotHandle.isValid()){
125  beamSpot = *recoBeamSpotHandle;
126  }else{
127  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
128  }
129 
130  bool validBS = true;
131  VertexState beamVertexState(beamSpot);
132  if ( (beamVertexState.error().cxx() <= 0.) ||
133  (beamVertexState.error().cyy() <= 0.) ||
134  (beamVertexState.error().czz() <= 0.) ) {
135  validBS = false;
136  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors "<<beamVertexState.error().matrix();
137  }
138 
139 
140  // get RECO tracks from the event
141  // `tks` can be used as a ptr to a reco::TrackCollection
143  iEvent.getByToken(trkToken, tks);
144 
145 
146  // interface RECO tracks to vertex reconstruction
148  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
149  std::vector<reco::TransientTrack> t_tks = (*theB).build(tks, beamSpot);
150  if(fVerbose) {std::cout << "RecoVertex/PrimaryVertexProducer"
151  << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
152  }
153 
154 
155  // select tracks
156  std::vector<reco::TransientTrack> && seltks = theTrackFilter->select( t_tks );
157 
158 
159  // clusterize tracks in Z
160  std::vector< std::vector<reco::TransientTrack> > && clusters = theTrackClusterizer->clusterize(seltks);
161  if (fVerbose){std::cout << " clustering returned "<< clusters.size() << " clusters from " << seltks.size() << " selected tracks" <<std::endl;}
162 
163 
164  // vertex fits
165  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
166 
167 
168  std::auto_ptr<reco::VertexCollection> result(new reco::VertexCollection);
169  reco::VertexCollection & vColl = (*result);
170 
171 
172  std::vector<TransientVertex> pvs;
173  for (std::vector< std::vector<reco::TransientTrack> >::const_iterator iclus
174  = clusters.begin(); iclus != clusters.end(); iclus++) {
175 
176 
178  if( algorithm->useBeamConstraint && validBS &&((*iclus).size()>1) ){
179 
180  v = algorithm->fitter->vertex(*iclus, beamSpot);
181 
182  }else if( !(algorithm->useBeamConstraint) && ((*iclus).size()>1) ) {
183 
184  v = algorithm->fitter->vertex(*iclus);
185 
186  }// else: no fit ==> v.isValid()=False
187 
188 
189  if (fVerbose){
190  if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " << v.position().z() << std::endl;
191  else std::cout <<"Invalid fitted vertex\n";
192  }
193 
194  if (v.isValid()
195  && (v.degreesOfFreedom()>=algorithm->minNdof)
196  && (!validBS || (*(algorithm->vertexSelector))(v,beamVertexState))
197  ) pvs.push_back(v);
198  }// end of cluster loop
199 
200  if(fVerbose){
201  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
202  }
203 
204 
205  if (clusters.size()>2 && clusters.size() > 2*pvs.size())
206  edm::LogWarning("PrimaryVertexProducer") << "more than half of candidate vertices lost " << pvs.size() << ' ' << clusters.size();
207 
208  if (pvs.empty() && seltks.size()>5)
209  edm::LogWarning("PrimaryVertexProducer") << "no vertex found with " << seltks.size() << " tracks and " << clusters.size() <<" vertex-candidates";
210 
211  // sort vertices by pt**2 vertex (aka signal vertex tagging)
212  if(pvs.size()>1){
213  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
214  }
215 
216 
217 
218  // convert transient vertices returned by the theAlgo to (reco) vertices
219  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin();
220  iv != pvs.end(); iv++) {
221  reco::Vertex v = *iv;
222  vColl.push_back(v);
223  }
224 
225  if (vColl.empty()) {
226  GlobalError bse(beamSpot.rotatedCovariance3D());
227  if ( (bse.cxx() <= 0.) ||
228  (bse.cyy() <= 0.) ||
229  (bse.czz() <= 0.) ) {
231  we(0,0)=10000; we(1,1)=10000; we(2,2)=10000;
232  vColl.push_back(reco::Vertex(beamSpot.position(), we,0.,0.,0));
233  if(fVerbose){
234  std::cout <<"RecoVertex/PrimaryVertexProducer: "
235  << "Beamspot with invalid errors "<<bse.matrix()<<std::endl;
236  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
237  }
238  } else {
239  vColl.push_back(reco::Vertex(beamSpot.position(),
240  beamSpot.rotatedCovariance3D(),0.,0.,0));
241  if(fVerbose){
242  std::cout <<"RecoVertex/PrimaryVertexProducer: "
243  << " will put Vertex derived from BeamSpot into Event.\n";
244  }
245  }
246  }
247 
248  if(fVerbose){
249  int ivtx=0;
250  for(reco::VertexCollection::const_iterator v=vColl.begin();
251  v!=vColl.end(); ++v){
252  std::cout << "recvtx "<< ivtx++
253  << "#trk " << std::setw(3) << v->tracksSize()
254  << " chi2 " << std::setw(4) << v->chi2()
255  << " ndof " << std::setw(3) << v->ndof()
256  << " x " << std::setw(6) << v->position().x()
257  << " dx " << std::setw(6) << v->xError()
258  << " y " << std::setw(6) << v->position().y()
259  << " dy " << std::setw(6) << v->yError()
260  << " z " << std::setw(6) << v->position().z()
261  << " dz " << std::setw(6) << v->zError()
262  << std::endl;
263  }
264  }
265 
266  iEvent.put(result, algorithm->label);
267  }
268 
269 }
270 
271 
272 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
tuple vertexCollections
edm::EDGetTokenT< reco::BeamSpot > bsToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
Common base class.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::TrackCollection > trkToken
const AlgebraicSymMatrix33 & matrix() const
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
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
int iEvent
Definition: GenABIO.cc:230
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
float degreesOfFreedom() const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
VertexCompatibleWithBeam * vertexSelector
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
TrackFilterForPVFindingBase * theTrackFilter
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
bool isValid() const
Definition: HandleBase.h:75
virtual void produce(edm::Event &, const edm::EventSetup &)
const T & get() const
Definition: EventSetup.h:56
std::vector< algo > algorithms
PrimaryVertexProducer(const edm::ParameterSet &)
GlobalError error() const
Definition: VertexState.h:55
tuple cout
Definition: gather_cfg.py:145
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