Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:09

0001 #include "RecoVertex/KalmanVertexFit/plugins/KVFTrackUpdate.h"
0002 
0003 #include "DataFormats/VertexReco/interface/Vertex.h"
0004 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0005 #include "DataFormats/TrackReco/interface/Track.h"
0006 #include "DataFormats/Common/interface/Handle.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0010 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0011 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0012 #include "RecoVertex/KalmanVertexFit/interface/SingleTrackVertexConstraint.h"
0013 #include "FWCore/Utilities/interface/InputTag.h"
0014 
0015 #include <iostream>
0016 
0017 using namespace reco;
0018 using namespace edm;
0019 using namespace std;
0020 
0021 KVFTrackUpdate::KVFTrackUpdate(const edm::ParameterSet& iConfig)
0022     : estoken_TTB(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))) {
0023   token_tracks = consumes<TrackCollection>(iConfig.getParameter<InputTag>("TrackLabel"));
0024   token_beamSpot = consumes<BeamSpot>(iConfig.getParameter<InputTag>("beamSpotLabel"));
0025 }
0026 
0027 KVFTrackUpdate::~KVFTrackUpdate() {}
0028 
0029 void KVFTrackUpdate::beginJob() {}
0030 
0031 void KVFTrackUpdate::endJob() {}
0032 
0033 //
0034 // member functions
0035 //
0036 
0037 void KVFTrackUpdate::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0038   try {
0039     edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Reconstructing event number: " << iEvent.id() << "\n";
0040 
0041     // get RECO tracks from the event
0042     // `tks` can be used as a ptr to a reco::TrackCollection
0043     edm::Handle<reco::TrackCollection> tks;
0044     iEvent.getByToken(token_tracks, tks);
0045 
0046     edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Found: " << (*tks).size() << " reconstructed tracks"
0047                                               << "\n";
0048     edm::LogPrint("RecoVertex/KVFTrackUpdate") << "got " << (*tks).size() << " tracks " << std::endl;
0049 
0050     // Transform Track to TransientTrack
0051 
0052     //get the builder:
0053     const auto& theB = &iSetup.getData(estoken_TTB);
0054     //do the conversion:
0055     std::vector<TransientTrack> t_tks = theB->build(tks);
0056 
0057     edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Found: " << t_tks.size() << " reconstructed tracks"
0058                                               << "\n";
0059 
0060     GlobalPoint glbPos(0., 0., 0.);
0061 
0062     AlgebraicSymMatrix33 mat;
0063     mat[0][0] = (20.e-04) * (20.e-04);
0064     mat[1][1] = (20.e-04) * (20.e-04);
0065     mat[2][2] = (5.3) * (5.3);
0066     GlobalError glbErrPos(mat);
0067 
0068     reco::BeamSpot vertexBeamSpot;
0069     edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
0070     iEvent.getByToken(token_beamSpot, recoBeamSpotHandle);
0071 
0072     SingleTrackVertexConstraint stvc;
0073     for (unsigned int i = 0; i < t_tks.size(); i++) {
0074       SingleTrackVertexConstraint::BTFtuple a = stvc.constrain(t_tks[i], glbPos, glbErrPos);
0075       edm::LogPrint("RecoVertex/KVFTrackUpdate") << "Chi2: " << std::get<2>(a) << std::endl;
0076       if (recoBeamSpotHandle.isValid()) {
0077         SingleTrackVertexConstraint::BTFtuple b = stvc.constrain(t_tks[i], *recoBeamSpotHandle);
0078         edm::LogPrint("RecoVertex/KVFTrackUpdate") << "Chi2: " << std::get<2>(b) << std::endl;
0079       }
0080     }
0081   }
0082 
0083   catch (std::exception& err) {
0084     edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Exception during event number: " << iEvent.id() << "\n"
0085                                               << err.what() << "\n";
0086   }
0087 }
0088 
0089 DEFINE_FWK_MODULE(KVFTrackUpdate);