CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoVertex/PrimaryVertexProducer/plugins/PrimaryVertexProducer.cc

Go to the documentation of this file.
00001 #include "RecoVertex/PrimaryVertexProducer/interface/PrimaryVertexProducer.h"
00002 
00003 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00005 #include "DataFormats/Common/interface/EDProduct.h"
00006 #include "DataFormats/Common/interface/Handle.h"
00007 #include "FWCore/Framework/interface/MakerMacros.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 
00010 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00011 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00012 
00013 #include "FWCore/Framework/interface/ESHandle.h"
00014 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00015 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00016 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00017 
00018 //using namespace reco;
00019 
00020 //
00021 // constants, enums and typedefs
00022 //
00023 
00024 //
00025 // static data member definitions
00026 //
00027 
00028 //
00029 // constructors and destructor
00030 //
00031 PrimaryVertexProducer::PrimaryVertexProducer(const edm::ParameterSet& conf)
00032   : theAlgo(conf), theConfig(conf)
00033 {
00034   edm::LogInfo("PVDebugInfo") 
00035     << "Initializing PV producer " << "\n";
00036   fVerbose=conf.getUntrackedParameter<bool>("verbose", false);
00037   trackLabel = conf.getParameter<edm::InputTag>("TrackLabel");
00038   beamSpotLabel = conf.getParameter<edm::InputTag>("beamSpotLabel");
00039  
00040   produces<reco::VertexCollection>();
00041 
00042 }
00043 
00044 
00045 PrimaryVertexProducer::~PrimaryVertexProducer() {}
00046 
00047 //
00048 // member functions
00049 //
00050 
00051 // ------------ method called to produce the data  ------------
00052 void
00053 PrimaryVertexProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00054 {
00055   using namespace edm;
00056 
00057   std::auto_ptr<reco::VertexCollection> result(new reco::VertexCollection);
00058   reco::VertexCollection vColl;
00059 
00060   // get the BeamSpot, it will alwys be needed, even when not used as a constraint
00061   reco::BeamSpot vertexBeamSpot;
00062   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00063   iEvent.getByLabel(beamSpotLabel,recoBeamSpotHandle);
00064   if (recoBeamSpotHandle.isValid()){
00065     vertexBeamSpot = *recoBeamSpotHandle;
00066   }else{
00067     edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
00068   }
00069 
00070 
00071 
00072   // get RECO tracks from the event
00073   // `tks` can be used as a ptr to a reco::TrackCollection
00074   edm::Handle<reco::TrackCollection> tks;
00075   iEvent.getByLabel(trackLabel, tks);
00076 
00077 
00078   // interface RECO tracks to vertex reconstruction
00079   edm::ESHandle<TransientTrackBuilder> theB;
00080   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
00081   std::vector<reco::TransientTrack> t_tks = (*theB).build(tks, vertexBeamSpot);
00082   if(fVerbose) {std::cout << "RecoVertex/PrimaryVertexProducer"
00083                      << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
00084   }
00085 
00086 
00087   // call vertex reconstruction
00088   std::vector<TransientVertex> t_vts = theAlgo.vertices(t_tks, vertexBeamSpot);
00089   if(fVerbose){
00090     std::cout <<"RecoVertex/PrimaryVertexProducer: "
00091               << " found " << t_vts.size() << " reconstructed vertices" << "\n";
00092   }
00093    
00094   // convert transient vertices returned by the theAlgo to (reco) vertices
00095   for (std::vector<TransientVertex>::const_iterator iv = t_vts.begin();
00096        iv != t_vts.end(); iv++) {
00097     reco::Vertex v = *iv;
00098     vColl.push_back(v);
00099   }
00100 
00101   if (vColl.empty()) {
00102     GlobalError bse(vertexBeamSpot.rotatedCovariance3D());
00103     if ( (bse.cxx() <= 0.) || 
00104         (bse.cyy() <= 0.) ||
00105         (bse.czz() <= 0.) ) {
00106       AlgebraicSymMatrix33 we;
00107       we(0,0)=10000; we(1,1)=10000; we(2,2)=10000;
00108       vColl.push_back(reco::Vertex(vertexBeamSpot.position(), we,0.,0.,0));
00109       if(fVerbose){
00110         std::cout <<"RecoVertex/PrimaryVertexProducer: "
00111                   << "Beamspot with invalid errors "<<bse.matrix()<<std::endl;
00112         std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
00113       }
00114     } else {
00115       vColl.push_back(reco::Vertex(vertexBeamSpot.position(), 
00116                                  vertexBeamSpot.rotatedCovariance3D(),0.,0.,0));
00117       if(fVerbose){
00118         std::cout <<"RecoVertex/PrimaryVertexProducer: "
00119                   << " will put Vertex derived from BeamSpot into Event.\n";
00120       }
00121     }
00122   }
00123 
00124   if(fVerbose){
00125     int ivtx=0;
00126     for(reco::VertexCollection::const_iterator v=vColl.begin(); 
00127         v!=vColl.end(); ++v){
00128       std::cout << "recvtx "<< ivtx++ 
00129                 << "#trk " << std::setw(3) << v->tracksSize()
00130                 << " chi2 " << std::setw(4) << v->chi2() 
00131                 << " ndof " << std::setw(3) << v->ndof() 
00132                 << " x "  << std::setw(6) << v->position().x() 
00133                 << " dx " << std::setw(6) << v->xError()
00134                 << " y "  << std::setw(6) << v->position().y() 
00135                 << " dy " << std::setw(6) << v->yError()
00136                 << " z "  << std::setw(6) << v->position().z() 
00137                 << " dz " << std::setw(6) << v->zError()
00138                 << std::endl;
00139     }
00140   }
00141 
00142   
00143   *result = vColl;
00144   //  iEvent.put(result, "PrimaryVertex");
00145   iEvent.put(result);
00146   
00147 }
00148 
00149 
00150 //define this as a plug-in
00151 DEFINE_FWK_MODULE(PrimaryVertexProducer);