Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /**\class TrackMix TrackMix.cc RecoVertex/TrackMix/src/TrackMix.cc
0002 
0003  Description: Simple test to see that reco::Vertex can store tracks from different
0004  Collections and of different types in the same vertex.
0005 */
0006 
0007 // system include files
0008 #include <memory>
0009 #include <iostream>
0010 
0011 // user include files
0012 #include "DataFormats/Common/interface/Handle.h"
0013 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0014 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
0015 #include "DataFormats/TrackReco/interface/Track.h"
0016 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0017 #include "DataFormats/VertexReco/interface/Vertex.h"
0018 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0019 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0020 #include "FWCore/Framework/interface/ESHandle.h"
0021 #include "FWCore/Framework/interface/Event.h"
0022 #include "FWCore/Framework/interface/EventSetup.h"
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/MakerMacros.h"
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0027 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0028 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0029 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0030 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0031 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0032 
0033 // ROOT includes
0034 #include <TFile.h>
0035 
0036 /**
0037    * This is a very simple test analyzer mean to test the KalmanVertexFitter
0038    */
0039 
0040 class TrackMix : public edm::one::EDAnalyzer<> {
0041 public:
0042   explicit TrackMix(const edm::ParameterSet&);
0043   ~TrackMix();
0044 
0045   virtual void analyze(const edm::Event&, const edm::EventSetup&);
0046 
0047 private:
0048   const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> estoken_ttk;
0049   edm::EDGetTokenT<edm::View<reco::Track> > token_gsf, token_ckf;
0050 };
0051 
0052 using namespace reco;
0053 using namespace edm;
0054 using namespace std;
0055 
0056 TrackMix::TrackMix(const edm::ParameterSet& iConfig)
0057     : estoken_ttk(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))) {
0058   token_gsf = consumes<edm::View<reco::Track> >(iConfig.getParameter<string>("gsfTrackLabel"));
0059   token_ckf = consumes<edm::View<reco::Track> >(iConfig.getParameter<string>("ckfTrackLabel"));
0060 }
0061 
0062 TrackMix::~TrackMix() = default;
0063 
0064 void TrackMix::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0065   try {
0066     edm::LogInfo("RecoVertex/TrackMix") << "Reconstructing event number: " << iEvent.id() << "\n";
0067 
0068     // get RECO tracks from the event
0069     // `tks` can be used as a ptr to a reco::TrackCollection
0070     edm::Handle<edm::View<reco::Track> > tks;
0071     iEvent.getByToken(token_gsf, tks);
0072     edm::Handle<edm::View<reco::Track> > tks2;
0073     iEvent.getByToken(token_ckf, tks2);
0074 
0075     edm::LogPrint("TrackMix") << "got " << (*tks).size() << " gsf tracks " << endl;
0076     edm::LogPrint("TrackMix") << "got " << (*tks2).size() << " ckf tracks " << endl;
0077 
0078     // Transform Track to TransientTrack
0079 
0080     //get the builder:
0081     edm::ESHandle<TransientTrackBuilder> theB = iSetup.getHandle(estoken_ttk);
0082     //do the conversion:
0083     vector<TransientTrack> t_tks = (*theB).build(tks);
0084     vector<TransientTrack> t_tks2 = (*theB).build(tks2);
0085     t_tks.insert(t_tks.end(), t_tks2.begin(), t_tks2.end());
0086 
0087     edm::LogPrint("TrackMix") << "Total: " << t_tks.size() << " reconstructed tracks"
0088                               << "\n";
0089 
0090     // Call the KalmanVertexFitter if more than 1 track
0091     if (t_tks.size() > 1) {
0092       KalmanVertexFitter kvf(true);
0093       TransientVertex tv = kvf.vertex(t_tks);
0094 
0095       edm::LogPrint("TrackMix") << "Position: " << Vertex::Point(tv.position()) << "\n";
0096 
0097       reco::Vertex v1 = tv;
0098       reco::Vertex::trackRef_iterator v1TrackIter;
0099       reco::Vertex::trackRef_iterator v1TrackBegin = v1.tracks_begin();
0100       reco::Vertex::trackRef_iterator v1TrackEnd = v1.tracks_end();
0101       edm::LogPrint("TrackMix") << v1.position() << v1.tracksSize() << endl;
0102       for (v1TrackIter = v1TrackBegin; v1TrackIter != v1TrackEnd; v1TrackIter++) {
0103         edm::LogPrint("TrackMix") << "pt" << (**v1TrackIter).pt() << endl;
0104         edm::LogPrint("TrackMix") << " weight " << v1.trackWeight(*v1TrackIter) << endl;
0105         edm::LogPrint("TrackMix") << " ref " << v1.refittedTrack(*v1TrackIter).pt() << endl;
0106       }
0107     }
0108   } catch (cms::Exception& err) {
0109     edm::LogError("TrackMix") << "Exception during event number: " << iEvent.id() << "\n" << err.what() << "\n";
0110   }
0111 }
0112 
0113 DEFINE_FWK_MODULE(TrackMix);