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