CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoPixelVertexing/PixelVertexFinding/src/DivisiveVertexFinder.cc

Go to the documentation of this file.
00001 #include "RecoPixelVertexing/PixelVertexFinding/interface/DivisiveVertexFinder.h"
00002 #include "RecoPixelVertexing/PixelVertexFinding/interface/PVPositionBuilder.h"
00003 #include "RecoPixelVertexing/PixelVertexFinding/interface/PVCluster.h"
00004 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "RecoPixelVertexing/PixelVertexFinding/interface/PVClusterComparer.h"
00007 #include <utility>
00008 #include <vector>
00009 #include <map>
00010 #include <algorithm>
00011 #include <cmath>
00012 
00013 DivisiveVertexFinder::DivisiveVertexFinder(double zOffset, int ntrkMin, 
00014                                            bool useError, double zSeparation, bool wtAverage,
00015                                            int verbosity)
00016   : zOffset_(zOffset), zSeparation_(zSeparation), ntrkMin_(ntrkMin), useError_(useError),
00017     wtAverage_(wtAverage),
00018     divmeth_(zOffset, ntrkMin, useError, zSeparation, wtAverage),
00019     verbose_(verbosity)
00020 {
00021 
00022 }
00023 
00024 DivisiveVertexFinder::~DivisiveVertexFinder(){}
00025 
00026 bool DivisiveVertexFinder::findVertexes(const reco::TrackRefVector &trks,  // input
00027                                         reco::VertexCollection &vertexes){ // output
00028   PVPositionBuilder pos;
00029   Measurement1D vz;
00030   if (wtAverage_) {
00031     vz = pos.wtAverage(trks);
00032   }
00033   else {
00034     vz = pos.average(trks);
00035   }
00036   reco::Vertex::Error err;
00037   err(2,2) = vz.error()*vz.error();
00038 
00039   reco::Vertex v( reco::Vertex::Point(0,0,vz.value()), err, 0, 1, trks.size() );
00040     for (unsigned int i=0; i<trks.size(); i++) {
00041       double vz = trks[i]->vz();
00042       if(std::isnan(vz)) continue;
00043       v.add(reco::TrackBaseRef(trks[i]));
00044     }
00045     
00046   vertexes.push_back(v);
00047 
00048   return true;
00049 }
00050 
00051 bool DivisiveVertexFinder::findVertexesAlt(const reco::TrackRefVector &trks,  // input
00052                                            reco::VertexCollection &vertexes, const math::XYZPoint & bs){ // output
00053   std::vector< PVCluster > in;
00054   std::pair< std::vector< PVCluster >, std::vector< const reco::Track* > > out;
00055   
00056   // Convert input track collection into container needed by Wolfgang's templated code
00057   // Need to save a map to reconvert from bare pointers, oy vey
00058   std::map< const reco::Track*, reco::TrackRef > mapa; 
00059   //  std::vector< std::vector< const reco::Track* > > trkps;
00060   for (unsigned int i=0; i<trks.size(); ++i) {
00061     double vz = trks[i]->vz();
00062     if(std::isnan(vz)) continue;
00063     std::vector< const reco::Track* > temp;
00064     temp.clear();
00065     temp.push_back( &(*trks[i]) );
00066 
00067     in.push_back( PVCluster( Measurement1D(trks[i]->dz(bs), trks[i]->dzError() ), temp ) );
00068     mapa[temp[0]] = trks[i];
00069   }
00070 
00071   if (verbose_ > 0 ) {
00072     edm::LogInfo("DivisiveVertexFinder") << "size of input vector of clusters " << in.size();
00073     for (unsigned int i=0; i<in.size(); ++i) {
00074       edm::LogInfo("DivisiveVertexFinder") << "Track " << i << " addr " << in[i].tracks()[0] 
00075                                            << " dz " << in[i].tracks()[0]->dz(bs)
00076                                            << " +- " << in[i].tracks()[0]->dzError()
00077                                            << " prodID " << mapa[in[i].tracks()[0]].id()
00078                                            << " dz from RefTrack " << mapa[in[i].tracks()[0]]->dz(bs)
00079                                            << " +- " << mapa[in[i].tracks()[0]]->dzError();
00080     }
00081   }
00082 
00083   // Run the darn thing
00084   divmeth_.setBeamSpot(bs);
00085   out = divmeth_(in);
00086 
00087   if (verbose_ > 0) edm::LogInfo("DivisiveVertexFinder") << " DivisiveClusterizer1D found " 
00088                                                          << out.first.size() << " vertexes";
00089 
00090   // Now convert the output yet again into something we can safely store in the event
00091   for (unsigned int iv=0; iv<out.first.size(); ++iv) { // loop over output vertexes
00092     reco::Vertex::Error err;
00093     err(2,2) = out.first[iv].position().error()*out.first[iv].position().error();
00094     
00095     reco::Vertex v( reco::Vertex::Point(0,0,out.first[iv].position().value()), err, 0, 1, out.second.size() );
00096     if (verbose_ > 0 ) edm::LogInfo("DivisiveVertexFinder") << " DivisiveClusterizer1D vertex " << iv 
00097                                                             << " has " << out.first[iv].tracks().size()
00098                                                             << " tracks and a position of " << v.z() 
00099                                                             << " +- " << std::sqrt(v.covariance(2,2));
00100     for (unsigned int itrk=0; itrk<out.first[iv].tracks().size(); ++itrk) {
00101       v.add( reco::TrackBaseRef(mapa[out.first[iv].tracks()[itrk]] ) );
00102     }
00103     vertexes.push_back(v); // Done with horrible conversion, save it
00104   }
00105 
00106   // Finally, sort the vertexes in decreasing sumPtSquared
00107   std::sort(vertexes.begin(), vertexes.end(), PVClusterComparer());
00108 
00109   return true;
00110 }