Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // user includes
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0004 #include "RecoVertex/ConfigurableVertexReco/interface/ConfigurableVertexReconstructor.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/ESHandle.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "MagneticField/Engine/interface/MagneticField.h"
0010 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0011 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0012 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0013 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0014 #include "DataFormats/Provenance/interface/RunLumiEventNumber.h"
0015 
0016 // system includes
0017 #include <vector>
0018 #include <iostream>
0019 
0020 class CVRTest : public edm::one::EDAnalyzer<> {
0021   /**
0022    *  Class that glues the combined btagging algorithm to the framework
0023    */
0024 public:
0025   explicit CVRTest(const edm::ParameterSet&);
0026   ~CVRTest();
0027 
0028   virtual void analyze(const edm::Event&, const edm::EventSetup&);
0029 
0030 private:
0031   void discussPrimary(const edm::Event&) const;
0032 
0033 private:
0034   const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> estoken_ttk;
0035   const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> estoken_mf;
0036 
0037   ConfigurableVertexReconstructor* vrec_;
0038   std::string trackcoll_;
0039   std::string vertexcoll_;
0040   std::string beamspot_;
0041 };
0042 
0043 using namespace std;
0044 using namespace reco;
0045 using namespace edm;
0046 
0047 namespace {
0048   void printTSOS(const TrajectoryStateOnSurface& tsos) {
0049     cout << tsos.globalPosition() << " , " << tsos.globalMomentum() << endl;
0050   }
0051 
0052   void printVertex(const TransientVertex& vtx) {
0053     cout << " `- pos=(" << vtx.position().x() << ", " << vtx.position().y() << ", " << vtx.position().z()
0054          << ") chi2=" << vtx.totalChiSquared() << " ndf=" << vtx.degreesOfFreedom() << " hr=" << vtx.hasRefittedTracks()
0055          << endl;
0056     if (vtx.originalTracks().size() && vtx.hasRefittedTracks()) {
0057       cout << "    `- 1st trk: ";
0058       reco::TransientTrack t = vtx.originalTracks()[0];
0059       TrajectoryStateOnSurface tsos = t.impactPointState();
0060       printTSOS(tsos);
0061       if (vtx.refittedTracks().size()) {
0062         cout << "     `- 1st refttd: ";
0063         reco::TransientTrack t2 = vtx.refittedTracks()[0];
0064         printTSOS(t2.impactPointState());
0065       }
0066     }
0067   }
0068 
0069   void printVertices(const vector<TransientVertex>& vtces) {
0070     cout << "[CVRTest] " << vtces.size() << " vertices." << endl;
0071     for (vector<TransientVertex>::const_iterator i = vtces.begin(); i != vtces.end(); ++i) {
0072       printVertex(*i);
0073       cout << endl;
0074     }
0075   }
0076 
0077   void discussBeamSpot(const reco::BeamSpot& bs) {
0078     cout << "[CVRTest] beamspot at " << bs.position() << endl;
0079     reco::BeamSpot::Covariance3DMatrix cov = bs.rotatedCovariance3D();
0080     cout << "[CVRTest] cov=" << cov << endl;
0081   }
0082 }  // namespace
0083 
0084 CVRTest::CVRTest(const edm::ParameterSet& iconfig)
0085     : estoken_ttk(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
0086       estoken_mf(esConsumes()),
0087       trackcoll_(iconfig.getParameter<string>("trackcoll")),
0088       vertexcoll_(iconfig.getParameter<string>("vertexcoll")),
0089       beamspot_(iconfig.getParameter<string>("beamspot")) {
0090   edm::ParameterSet vtxconfig = iconfig.getParameter<edm::ParameterSet>("vertexreco");
0091   vrec_ = new ConfigurableVertexReconstructor(vtxconfig);
0092   cout << "[CVRTest] vtxconfig=" << vtxconfig << endl;
0093 }
0094 
0095 CVRTest::~CVRTest() {
0096   if (vrec_)
0097     delete vrec_;
0098 }
0099 
0100 void CVRTest::discussPrimary(const edm::Event& iEvent) const {
0101   edm::Handle<reco::VertexCollection> retColl;
0102   iEvent.getByLabel(vertexcoll_, retColl);
0103   if (retColl->size()) {
0104     const reco::Vertex& vtx = *(retColl->begin());
0105     cout << "[CVRTest] persistent primary: " << vtx.x() << ", " << vtx.y() << ", " << vtx.z() << endl;
0106   }
0107 }
0108 
0109 void CVRTest::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0110   edm::EventNumber_t const evt = iEvent.id().event();
0111   cout << "[CVRTest] next event: " << evt << endl;
0112   edm::ESHandle<MagneticField> magneticField = iSetup.getHandle(estoken_mf);
0113   edm::ESHandle<TransientTrackBuilder> builder = iSetup.getHandle(estoken_ttk);
0114 
0115   edm::Handle<reco::TrackCollection> tks;
0116   iEvent.getByLabel(trackcoll_, tks);
0117   discussPrimary(iEvent);
0118 
0119   edm::Handle<reco::BeamSpot> bs;
0120   iEvent.getByLabel(beamspot_, bs);
0121   discussBeamSpot(*bs);
0122 
0123   vector<reco::TransientTrack> ttks;
0124   ttks = builder->build(tks);
0125   cout << "[CVRTest] got " << ttks.size() << " tracks." << endl;
0126 
0127   cout << "[CVRTest] fit w/o beamspot constraint" << endl;
0128   vector<TransientVertex> vtces = vrec_->vertices(ttks);
0129   printVertices(vtces);
0130 
0131   // cout << "[CVRTest] fit w beamspot constraint" << endl;
0132   // vector < TransientVertex > bvtces = vrec_->vertices ( ttks, *bs );
0133   // printVertices ( bvtces );
0134 }
0135 
0136 //define this as a plug-in
0137 DEFINE_FWK_MODULE(CVRTest);