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,
00027 reco::VertexCollection &vertexes){
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 v.add(reco::TrackBaseRef(trks[i]));
00042 }
00043
00044 vertexes.push_back(v);
00045
00046 return true;
00047 }
00048
00049 bool DivisiveVertexFinder::findVertexesAlt(const reco::TrackRefVector &trks,
00050 reco::VertexCollection &vertexes){
00051 std::vector< PVCluster > in;
00052 std::pair< std::vector< PVCluster >, std::vector< const reco::Track* > > out;
00053
00054
00055
00056 std::map< const reco::Track*, reco::TrackRef > mapa;
00057
00058 for (unsigned int i=0; i<trks.size(); ++i) {
00059 std::vector< const reco::Track* > temp;
00060 temp.clear();
00061 temp.push_back( &(*trks[i]) );
00062
00063 in.push_back( PVCluster( Measurement1D(trks[i]->dz(), trks[i]->dzError() ), temp ) );
00064 mapa[temp[0]] = trks[i];
00065 }
00066
00067 if (verbose_ > 0 ) {
00068 edm::LogInfo("DivisiveVertexFinder") << "size of input vector of clusters " << in.size();
00069 for (unsigned int i=0; i<in.size(); ++i) {
00070 edm::LogInfo("DivisiveVertexFinder") << "Track " << i << " addr " << in[i].tracks()[0]
00071 << " dz " << in[i].tracks()[0]->dz()
00072 << " +- " << in[i].tracks()[0]->dzError()
00073 << " prodID " << mapa[in[i].tracks()[0]].id()
00074 << " dz from RefTrack " << mapa[in[i].tracks()[0]]->dz()
00075 << " +- " << mapa[in[i].tracks()[0]]->dzError();
00076 }
00077 }
00078
00079
00080 out = divmeth_(in);
00081
00082 if (verbose_ > 0) edm::LogInfo("DivisiveVertexFinder") << " DivisiveClusterizer1D found "
00083 << out.first.size() << " vertexes";
00084
00085
00086 for (unsigned int iv=0; iv<out.first.size(); ++iv) {
00087 reco::Vertex::Error err;
00088 err(2,2) = out.first[iv].position().error()*out.first[iv].position().error();
00089
00090 reco::Vertex v( reco::Vertex::Point(0,0,out.first[iv].position().value()), err, 0, 1, out.second.size() );
00091 if (verbose_ > 0 ) edm::LogInfo("DivisiveVertexFinder") << " DivisiveClusterizer1D vertex " << iv
00092 << " has " << out.first[iv].tracks().size()
00093 << " tracks and a position of " << v.z()
00094 << " +- " << std::sqrt(v.error(2,2));
00095 for (unsigned int itrk=0; itrk<out.first[iv].tracks().size(); ++itrk) {
00096 v.add( reco::TrackBaseRef(mapa[out.first[iv].tracks()[itrk]] ) );
00097 }
00098 vertexes.push_back(v);
00099 }
00100
00101
00102 std::sort(vertexes.begin(), vertexes.end(), PVClusterComparer());
00103
00104 return true;
00105 }