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 
8 
12 
17 
18 
20  :theConfig(conf)
21 {
22 
23  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
24  trackLabel = conf.getParameter<edm::InputTag>("TrackLabel");
25  beamSpotLabel = conf.getParameter<edm::InputTag>("beamSpotLabel");
26 
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 #ifdef __GXX_EXPERIMENTAL_CXX0X__
48  else if(clusteringAlgorithm == "DA_vect") {
49  theTrackClusterizer = new DAClusterizerInZ_vect(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
50  }
51 #endif
52 
53 
54  else{
55  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
56  }
57 
58 
59  // select and configure the vertex fitters
60  if (conf.exists("vertexCollections")){
61  std::vector<edm::ParameterSet> vertexCollections =conf.getParameter< std::vector<edm::ParameterSet> >("vertexCollections");
62 
63  for( std::vector< edm::ParameterSet >::const_iterator algoconf = vertexCollections.begin(); algoconf != vertexCollections.end(); algoconf++){
64 
66  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
67  if (fitterAlgorithm=="KalmanVertexFitter") {
68  algorithm.fitter= new KalmanVertexFitter();
69  } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
70  algorithm.fitter= new AdaptiveVertexFitter();
71  } else {
72  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
73  }
74  algorithm.label = algoconf->getParameter<std::string>("label");
75  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
76  algorithm.useBeamConstraint=algoconf->getParameter<bool>("useBeamConstraint");
77  algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
78  algorithms.push_back(algorithm);
79 
80  produces<reco::VertexCollection>(algorithm.label);
81  }
82  }else{
83  edm::LogWarning("MisConfiguration")<<"this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
84 
86  std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
87  if (fitterAlgorithm=="KalmanVertexFitter") {
88  algorithm.fitter= new KalmanVertexFitter();
89  } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
90  algorithm.fitter= new AdaptiveVertexFitter();
91  } else {
92  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
93  }
94  algorithm.label = "";
95  algorithm.minNdof = conf.getParameter<double>("minNdof");
96  algorithm.useBeamConstraint=conf.getParameter<bool>("useBeamConstraint");
97 
98  algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
99 
100  algorithms.push_back(algorithm);
101  produces<reco::VertexCollection>(algorithm.label);
102  }
103 
104 
105 }
106 
107 
109  if (theTrackFilter) delete theTrackFilter;
111  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
112  if (algorithm->fitter) delete algorithm->fitter;
113  if (algorithm->vertexSelector) delete algorithm->vertexSelector;
114  }
115 }
116 
117 
118 void
120 {
121 
122  // get the BeamSpot, it will alwys be needed, even when not used as a constraint
124  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
125  iEvent.getByLabel(beamSpotLabel,recoBeamSpotHandle);
126  if (recoBeamSpotHandle.isValid()){
127  beamSpot = *recoBeamSpotHandle;
128  }else{
129  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
130  }
131 
132  bool validBS = true;
133  VertexState beamVertexState(beamSpot);
134  if ( (beamVertexState.error().cxx() <= 0.) ||
135  (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 
142  // get RECO tracks from the event
143  // `tks` can be used as a ptr to a reco::TrackCollection
145  iEvent.getByLabel(trackLabel, tks);
146 
147 
148  // interface RECO tracks to vertex reconstruction
150  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
151  std::vector<reco::TransientTrack> t_tks = (*theB).build(tks, beamSpot);
152  if(fVerbose) {std::cout << "RecoVertex/PrimaryVertexProducer"
153  << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
154  }
155 
156 
157  // select tracks
158  std::vector<reco::TransientTrack> seltks = theTrackFilter->select( t_tks );
159 
160 
161  // clusterize tracks in Z
162  std::vector< std::vector<reco::TransientTrack> > clusters = theTrackClusterizer->clusterize(seltks);
163  if (fVerbose){std::cout << " clustering returned "<< clusters.size() << " clusters from " << seltks.size() << " selected tracks" <<std::endl;}
164 
165 
166  // vertex fits
167  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
168 
169 
170  std::auto_ptr<reco::VertexCollection> result(new reco::VertexCollection);
172 
173 
174  std::vector<TransientVertex> pvs;
175  for (std::vector< std::vector<reco::TransientTrack> >::const_iterator iclus
176  = clusters.begin(); iclus != clusters.end(); iclus++) {
177 
178 
180  if( algorithm->useBeamConstraint && validBS &&((*iclus).size()>1) ){
181 
182  v = algorithm->fitter->vertex(*iclus, beamSpot);
183 
184  }else if( !(algorithm->useBeamConstraint) && ((*iclus).size()>1) ) {
185 
186  v = algorithm->fitter->vertex(*iclus);
187 
188  }// else: no fit ==> v.isValid()=False
189 
190 
191  if (fVerbose){
192  if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " << v.position().z() << std::endl;
193  else std::cout <<"Invalid fitted vertex\n";
194  }
195 
196  if (v.isValid()
197  && (v.degreesOfFreedom()>=algorithm->minNdof)
198  && (!validBS || (*(algorithm->vertexSelector))(v,beamVertexState))
199  ) pvs.push_back(v);
200  }// end of cluster loop
201 
202  if(fVerbose){
203  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
204  }
205 
206 
207 
208 
209 
210  // sort vertices by pt**2 vertex (aka signal vertex tagging)
211  if(pvs.size()>1){
212  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
213  }
214 
215 
216 
217  // convert transient vertices returned by the theAlgo to (reco) vertices
218  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin();
219  iv != pvs.end(); iv++) {
220  reco::Vertex v = *iv;
221  vColl.push_back(v);
222  }
223 
224  if (vColl.empty()) {
225  GlobalError bse(beamSpot.rotatedCovariance3D());
226  if ( (bse.cxx() <= 0.) ||
227  (bse.cyy() <= 0.) ||
228  (bse.czz() <= 0.) ) {
230  we(0,0)=10000; we(1,1)=10000; we(2,2)=10000;
231  vColl.push_back(reco::Vertex(beamSpot.position(), we,0.,0.,0));
232  if(fVerbose){
233  std::cout <<"RecoVertex/PrimaryVertexProducer: "
234  << "Beamspot with invalid errors "<<bse.matrix()<<std::endl;
235  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
236  }
237  } else {
238  vColl.push_back(reco::Vertex(beamSpot.position(),
239  beamSpot.rotatedCovariance3D(),0.,0.,0));
240  if(fVerbose){
241  std::cout <<"RecoVertex/PrimaryVertexProducer: "
242  << " will put Vertex derived from BeamSpot into Event.\n";
243  }
244  }
245  }
246 
247  if(fVerbose){
248  int ivtx=0;
249  for(reco::VertexCollection::const_iterator v=vColl.begin();
250  v!=vColl.end(); ++v){
251  std::cout << "recvtx "<< ivtx++
252  << "#trk " << std::setw(3) << v->tracksSize()
253  << " chi2 " << std::setw(4) << v->chi2()
254  << " ndof " << std::setw(3) << v->ndof()
255  << " x " << std::setw(6) << v->position().x()
256  << " dx " << std::setw(6) << v->xError()
257  << " y " << std::setw(6) << v->position().y()
258  << " dy " << std::setw(6) << v->yError()
259  << " z " << std::setw(6) << v->position().z()
260  << " dz " << std::setw(6) << v->zError()
261  << std::endl;
262  }
263  }
264 
265 
266  *result = vColl;
267  iEvent.put(result, algorithm->label);
268  }
269 
270 }
271 
272 
273 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
Common base class.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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:243
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:94
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:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
tuple conf
Definition: dbtoconf.py:185
virtual void produce(edm::Event &, const edm::EventSetup &)
const T & get() const
Definition: EventSetup.h:55
std::vector< algo > algorithms
PrimaryVertexProducer(const edm::ParameterSet &)
GlobalError error() const
Definition: VertexState.h:34
tuple cout
Definition: gather_cfg.py:121
const Point & position() const
position
Definition: BeamSpot.h:63
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:79
T x() const
Definition: PV3DBase.h:62
bool isValid() const