CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/RecoVertex/KalmanVertexFit/plugins/KVFTrackUpdate.cc

Go to the documentation of this file.
00001 #include "RecoVertex/KalmanVertexFit/plugins/KVFTrackUpdate.h"
00002 
00003 #include "DataFormats/VertexReco/interface/Vertex.h"
00004 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00005 #include "DataFormats/TrackReco/interface/Track.h"
00006 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00007 #include "DataFormats/Common/interface/Handle.h"
00008 #include "FWCore/Framework/interface/MakerMacros.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "FWCore/Framework/interface/ESHandle.h"
00011 
00012 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00013 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00014 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00015 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00016 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00017 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
00018 #include "RecoVertex/KalmanVertexFit/interface/SingleTrackVertexConstraint.h"
00019 
00020 #include <iostream>
00021 
00022 using namespace reco;
00023 using namespace edm;
00024 using namespace std;
00025 
00026 KVFTrackUpdate::KVFTrackUpdate(const edm::ParameterSet& iConfig)
00027 {
00028   trackLabel_ = iConfig.getParameter<edm::InputTag>("TrackLabel");
00029   beamSpotLabel = iConfig.getParameter<edm::InputTag>("beamSpotLabel");
00030 }
00031 
00032 
00033 KVFTrackUpdate::~KVFTrackUpdate() {
00034 }
00035 
00036 void KVFTrackUpdate::beginJob(){
00037 }
00038 
00039 
00040 void KVFTrackUpdate::endJob() {
00041 }
00042 
00043 //
00044 // member functions
00045 //
00046 
00047 void
00048 KVFTrackUpdate::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00049 {
00050 
00051 
00052 
00053   try {
00054     edm::LogInfo("RecoVertex/KVFTrackUpdate") 
00055       << "Reconstructing event number: " << iEvent.id() << "\n";
00056     
00057     // get RECO tracks from the event
00058     // `tks` can be used as a ptr to a reco::TrackCollection
00059     edm::Handle<reco::TrackCollection> tks;
00060     iEvent.getByLabel(trackLabel_, tks);
00061 
00062     edm::LogInfo("RecoVertex/KVFTrackUpdate") 
00063       << "Found: " << (*tks).size() << " reconstructed tracks" << "\n";
00064     std::cout << "got " << (*tks).size() << " tracks " << std::endl;
00065 
00066     // Transform Track to TransientTrack
00067 
00068     //get the builder:
00069     edm::ESHandle<TransientTrackBuilder> theB;
00070     iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
00071     //do the conversion:
00072     std::vector<TransientTrack> t_tks = (*theB).build(tks);
00073 
00074     edm::LogInfo("RecoVertex/KVFTrackUpdate") 
00075       << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
00076     
00077     GlobalPoint glbPos(0.,0.,0.);
00078 
00079     AlgebraicSymMatrix33 mat;
00080     mat[0][0] = (20.e-04)*(20.e-04);
00081     mat[1][1] = (20.e-04)*(20.e-04);
00082     mat[2][2] = (5.3)*(5.3);
00083     GlobalError glbErrPos(mat);
00084 
00085     reco::BeamSpot vertexBeamSpot;
00086     edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00087     iEvent.getByLabel(beamSpotLabel,recoBeamSpotHandle);
00088 
00089 
00090     SingleTrackVertexConstraint stvc;
00091     for (unsigned int i = 0; i<t_tks.size();i++) {
00092       SingleTrackVertexConstraint::BTFtuple a = 
00093         stvc.constrain(t_tks[i], glbPos, glbErrPos);
00094       std::cout << "Chi2: "<< a.get<2>()<<std::endl;
00095       if (recoBeamSpotHandle.isValid()){
00096         SingleTrackVertexConstraint::BTFtuple b =
00097           stvc.constrain(t_tks[i], *recoBeamSpotHandle);
00098         std::cout << "Chi2: "<< b.get<2>()<<std::endl;
00099       }
00100     }
00101   }
00102 
00103 
00104   catch (std::exception & err) {
00105     edm::LogInfo("RecoVertex/KVFTrackUpdate") 
00106       << "Exception during event number: " << iEvent.id() 
00107       << "\n" << err.what() << "\n";
00108   }
00109 
00110 }
00111 
00112 DEFINE_FWK_MODULE(KVFTrackUpdate);