CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoPixelVertexing/PixelVertexFinding/src/PixelVertexProducer.cc

Go to the documentation of this file.
00001 #include "RecoPixelVertexing/PixelVertexFinding/interface/PixelVertexProducer.h"
00002 #include "FWCore/Framework/interface/MakerMacros.h"
00003 #include "DataFormats/VertexReco/interface/Vertex.h"
00004 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00005 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00006 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
00007 #include "RecoPixelVertexing/PixelVertexFinding/interface/DivisiveVertexFinder.h"
00008 #include "FWCore/Utilities/interface/InputTag.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include <memory>
00011 #include <string>
00012 #include <cmath>
00013 
00014 PixelVertexProducer::PixelVertexProducer(const edm::ParameterSet& conf) 
00015   : conf_(conf), verbose_(0), dvf_(0), ptMin_(1.0)
00016 {
00017   // Register my product
00018   produces<reco::VertexCollection>();
00019 
00020   // Setup shop
00021   verbose_           = conf.getParameter<int>("Verbosity"); // 0 silent, 1 chatty, 2 loud
00022   std::string finder = conf.getParameter<std::string>("Finder"); // DivisiveVertexFinder
00023   bool useError      = conf.getParameter<bool>("UseError"); // true
00024   bool wtAverage     = conf.getParameter<bool>("WtAverage"); // true
00025   double zOffset     = conf.getParameter<double>("ZOffset"); // 5.0 sigma
00026   double zSeparation = conf.getParameter<double>("ZSeparation"); // 0.05 cm
00027   int ntrkMin        = conf.getParameter<int>("NTrkMin"); // 3
00028   // Tracking requirements before sending a track to be considered for vtx
00029   ptMin_ = conf_.getParameter<double>("PtMin"); // 1.0 GeV
00030 
00031   if (finder == "DivisiveVertexFinder") {
00032     if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Using the DivisiveVertexFinder\n";
00033     dvf_ = new DivisiveVertexFinder(zOffset, ntrkMin, useError, zSeparation, wtAverage, verbose_);
00034   }
00035   else { // Finder not supported, or you made a mistake in your request
00036     // throw an exception once I figure out how CMSSW does this
00037   }
00038 }
00039 
00040 
00041 PixelVertexProducer::~PixelVertexProducer() {
00042   delete dvf_;
00043 }
00044 
00045 void PixelVertexProducer::produce(edm::Event& e, const edm::EventSetup& es) {
00046 
00047   // First fish the pixel tracks out of the event
00048   edm::Handle<reco::TrackCollection> trackCollection;
00049   edm::InputTag trackCollName = conf_.getParameter<edm::InputTag>("TrackCollection");
00050   e.getByLabel(trackCollName,trackCollection);
00051   const reco::TrackCollection tracks = *(trackCollection.product());
00052   if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Found " << tracks.size() << " tracks in TrackCollection called " << trackCollName << "\n";
00053   
00054 
00055   // Second, make a collection of pointers to the tracks we want for the vertex finder
00056   reco::TrackRefVector trks;
00057   for (unsigned int i=0; i<tracks.size(); i++) {
00058     if (tracks[i].pt() > ptMin_)     
00059       trks.push_back( reco::TrackRef(trackCollection, i) );
00060   }
00061   if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << ": Selected " << trks.size() << " of these tracks for vertexing\n";
00062 
00063   edm::Handle<reco::BeamSpot> bsHandle;
00064   edm::InputTag bsName = conf_.getParameter<edm::InputTag>("beamSpot");
00065   e.getByLabel(bsName,bsHandle);
00066   math::XYZPoint myPoint(0.,0.,0.);
00067   if (bsHandle.isValid()) myPoint = math::XYZPoint(bsHandle->x0(),bsHandle->y0(), 0. ); //FIXME: fix last coordinate with vertex.z() at same time
00068 
00069   // Third, ship these tracks off to be vertexed
00070   std::auto_ptr<reco::VertexCollection> vertexes(new reco::VertexCollection);
00071   bool ok;
00072   if (conf_.getParameter<bool>("Method2")) {
00073     ok = dvf_->findVertexesAlt(trks,       // input
00074                              *vertexes,myPoint); // output
00075     if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method2 returned status of " << ok;
00076   }
00077   else {
00078     ok = dvf_->findVertexes(trks,       // input
00079                             *vertexes); // output
00080     if (verbose_ > 0) edm::LogInfo("PixelVertexProducer") << "Method1 returned status of " << ok;
00081   }
00082 
00083   if (verbose_ > 0) {
00084     edm::LogInfo("PixelVertexProducer") << ": Found " << vertexes->size() << " vertexes\n";
00085     for (unsigned int i=0; i<vertexes->size(); ++i) {
00086       edm::LogInfo("PixelVertexProducer") << "Vertex number " << i << " has " << (*vertexes)[i].tracksSize() << " tracks with a position of " << (*vertexes)[i].z() << " +- " << std::sqrt( (*vertexes)[i].covariance(2,2) );
00087     }
00088   }
00089 
00090 
00091   if(bsHandle.isValid())
00092     {
00093       const reco::BeamSpot & bs = *bsHandle;
00094       
00095       for (unsigned int i=0; i<vertexes->size(); ++i) {
00096         double z=(*vertexes)[i].z();
00097         double x=bs.x0()+bs.dxdz()*(z-bs.z0());
00098         double y=bs.y0()+bs.dydz()*(z-bs.z0()); 
00099         reco::Vertex v( reco::Vertex::Point(x,y,z), (*vertexes)[i].error(),(*vertexes)[i].chi2() , (*vertexes)[i].ndof() , (*vertexes)[i].tracksSize());
00100         //Copy also the tracks 
00101         for (std::vector<reco::TrackBaseRef >::const_iterator it = (*vertexes)[i].tracks_begin();
00102              it !=(*vertexes)[i].tracks_end(); it++) {
00103           v.add( *it );
00104         }
00105         (*vertexes)[i]=v;
00106         
00107       }
00108     }
00109   else
00110     {
00111       edm::LogWarning("PixelVertexProducer") << "No beamspot found. Using returning vertexes with (0,0,Z) ";
00112     } 
00113   
00114   if(!vertexes->size() && bsHandle.isValid()){
00115     
00116     const reco::BeamSpot & bs = *bsHandle;
00117       
00118       GlobalError bse(bs.rotatedCovariance3D());
00119       if ( (bse.cxx() <= 0.) ||
00120            (bse.cyy() <= 0.) ||
00121            (bse.czz() <= 0.) ) {
00122         AlgebraicSymMatrix33 we;
00123         we(0,0)=10000;
00124         we(1,1)=10000;
00125         we(2,2)=10000;
00126         vertexes->push_back(reco::Vertex(bs.position(), we,0.,0.,0));
00127         
00128         edm::LogInfo("PixelVertexProducer") <<"No vertices found. Beamspot with invalid errors " << bse.matrix() << std::endl
00129                                                << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n"
00130                                                << (*vertexes)[0].x() << "\n"
00131                                                << (*vertexes)[0].y() << "\n"
00132                                                << (*vertexes)[0].z() << "\n";
00133       } else {
00134         vertexes->push_back(reco::Vertex(bs.position(),
00135                                          bs.rotatedCovariance3D(),0.,0.,0));
00136         
00137         edm::LogInfo("PixelVertexProducer") << "No vertices found. Will put Vertex derived from BeamSpot into Event:\n"
00138                                                << (*vertexes)[0].x() << "\n"
00139                                                << (*vertexes)[0].y() << "\n"
00140                                                << (*vertexes)[0].z() << "\n";
00141       }
00142   }
00143       
00144   else if(!vertexes->size() && !bsHandle.isValid())
00145     {
00146       edm::LogWarning("PixelVertexProducer") << "No beamspot and no vertex found. No vertex returned.";
00147     }
00148   
00149   e.put(vertexes);
00150   
00151 }
00152 
00153