CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoHI/HiTracking/plugins/HIBestVertexProducer.cc

Go to the documentation of this file.
00001 #include "RecoHI/HiTracking/interface/HIBestVertexProducer.h"
00002 
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "FWCore/Framework/interface/Frameworkfwd.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 
00008 #include "DataFormats/VertexReco/interface/Vertex.h"
00009 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00010 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 
00013 #include <algorithm>
00014 #include <iostream>
00015 using namespace std;
00016 using namespace edm;
00017 
00018 /*****************************************************************************/
00019 HIBestVertexProducer::HIBestVertexProducer
00020 (const edm::ParameterSet& ps) : theConfig(ps),
00021   theMedianVertexCollection(ps.getParameter<edm::InputTag>("medianVertexCollection")),
00022   theAdaptiveVertexCollection(ps.getParameter<edm::InputTag>("adaptiveVertexCollection"))
00023 {
00024   produces<reco::VertexCollection>();
00025 }
00026 
00027 
00028 /*****************************************************************************/
00029 HIBestVertexProducer::~HIBestVertexProducer()
00030 { 
00031 }
00032 
00033 /*****************************************************************************/
00034 void HIBestVertexProducer::beginJob()
00035 {
00036 }
00037 
00038 /*****************************************************************************/
00039 void HIBestVertexProducer::produce
00040 (edm::Event& ev, const edm::EventSetup& es)
00041 {
00042   
00043   // 1. use best adaptive vertex preferentially
00044   // 2. use median vertex in case where adaptive algorithm failed
00045   // 3. use beamspot if netither vertexing method succeeds
00046 
00047   // New vertex collection
00048   std::auto_ptr<reco::VertexCollection> newVertexCollection(new reco::VertexCollection);
00049 
00050   //** Get precise adaptive vertex **/
00051   edm::Handle<reco::VertexCollection> vc1;
00052   ev.getByLabel(theAdaptiveVertexCollection, vc1);
00053   const reco::VertexCollection *vertices1 = vc1.product();
00054   
00055   if(vertices1->size()==0)
00056     LogError("HeavyIonVertexing") << "adaptive vertex collection is empty!" << endl;
00057 
00058   if(vertices1->begin()->zError()<3) { 
00059   
00060     reco::VertexCollection::const_iterator vertex1 = vertices1->begin();
00061     newVertexCollection->push_back(*vertex1);
00062 
00063     LogInfo("HeavyIonVertexing") << "adaptive vertex:\n vz = (" 
00064                                  << vertex1->x() << ", " << vertex1->y() << ", " << vertex1->z() << ")" 
00065                                  << "\n error = ("
00066                                  << vertex1->xError() << ", " << vertex1->yError() << ", " 
00067                                  << vertex1->zError() << ")" << endl;
00068   
00069   } else {
00070     
00071     //** Get fast median vertex **/
00072     edm::Handle<reco::VertexCollection> vc2;
00073     ev.getByLabel(theMedianVertexCollection, vc2);
00074     const reco::VertexCollection * vertices2 = vc2.product();
00075     
00076     //** Get beam spot position and error **/
00077     reco::BeamSpot beamSpot;
00078     edm::Handle<reco::BeamSpot> beamSpotHandle;
00079     ev.getByLabel("offlineBeamSpot", beamSpotHandle);
00080 
00081     if( beamSpotHandle.isValid() ) 
00082       beamSpot = *beamSpotHandle;
00083     else
00084       LogError("HeavyIonVertexing") << "no beamspot with name: 'offlineBeamSpot'" << endl;
00085 
00086     if(vertices2->size() > 0) { 
00087       
00088       reco::VertexCollection::const_iterator vertex2 = vertices2->begin();
00089       reco::Vertex::Error err;
00090       err(0,0)=pow(beamSpot.BeamWidthX(),2);
00091       err(1,1)=pow(beamSpot.BeamWidthY(),2);
00092       err(2,2)=pow(vertex2->zError(),2);
00093       reco::Vertex newVertex(reco::Vertex::Point(beamSpot.x0(),beamSpot.y0(),vertex2->z()),
00094                              err, 0, 1, 1);
00095       newVertexCollection->push_back(newVertex);  
00096 
00097       LogInfo("HeavyIonVertexing") << "median vertex + beamspot: \n position = (" 
00098                                    << newVertex.x() << ", " << newVertex.y() << ", " << newVertex.z() << ")" 
00099                                    << "\n error = ("
00100                                    << newVertex.xError() << ", " << newVertex.yError() << ", " 
00101                                    << newVertex.zError() << ")" << endl;
00102       
00103     } else { 
00104       
00105       reco::Vertex::Error err;
00106       err(0,0)=pow(beamSpot.BeamWidthX(),2);
00107       err(1,1)=pow(beamSpot.BeamWidthY(),2);
00108       err(2,2)=pow(beamSpot.sigmaZ(),2);
00109       reco::Vertex newVertex(beamSpot.position(), 
00110                              err, 0, 0, 1);
00111       newVertexCollection->push_back(newVertex);  
00112 
00113       LogInfo("HeavyIonVertexing") << "beam spot: \n position = (" 
00114                                    << newVertex.x() << ", " << newVertex.y() << ", " << newVertex.z() << ")" 
00115                                    << "\n error = ("
00116                                    << newVertex.xError() << ", " << newVertex.yError() << ", " 
00117                                    << newVertex.zError() << ")" << endl;
00118       
00119     }
00120     
00121   }
00122   
00123   // put new vertex collection into event
00124   ev.put(newVertexCollection);
00125   
00126 }
00127