Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:28

0001 // -*- C++ -*-
0002 //
0003 // Package:    IPAnalyzer
0004 // Class:      IPAnalyzer
0005 //
0006 /**\class IPAnalyzer IPAnalyzer.cc RecoBTag/IPAnalyzer/src/IPAnalyzer.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Andrea Rizzi
0015 //         Created:  Wed Apr 12 11:12:49 CEST 2006
0016 //
0017 //
0018 
0019 // user include files
0020 #include "FWCore/Framework/interface/Frameworkfwd.h"
0021 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0022 #include "FWCore/Framework/interface/Event.h"
0023 #include "FWCore/Framework/interface/MakerMacros.h"
0024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0025 #include "FWCore/Utilities/interface/InputTag.h"
0026 
0027 #include "DataFormats/Common/interface/Ref.h"
0028 #include "DataFormats/JetReco/interface/Jet.h"
0029 #include "DataFormats/JetReco/interface/CaloJet.h"
0030 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
0031 #include "DataFormats/TrackReco/interface/Track.h"
0032 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
0033 
0034 #include "DataFormats/Math/interface/Vector3D.h"
0035 
0036 // Math
0037 #include "Math/GenVector/VectorUtil.h"
0038 #include "Math/GenVector/PxPyPzE4D.h"
0039 #include "DataFormats/VertexReco/interface/Vertex.h"
0040 
0041 // system include files
0042 #include <string>
0043 #include <iostream>
0044 
0045 using namespace std;
0046 using namespace reco;
0047 
0048 //
0049 // class decleration
0050 //
0051 
0052 class IPAnalyzer : public edm::one::EDAnalyzer<> {
0053 public:
0054   explicit IPAnalyzer(const edm::ParameterSet&);
0055   ~IPAnalyzer() {}
0056 
0057   virtual void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup);
0058 
0059 private:
0060   edm::InputTag m_assoc;
0061   edm::InputTag m_jets;
0062   edm::EDGetTokenT<reco::TrackIPTagInfoCollection> token_ipassoc;
0063 };
0064 
0065 //
0066 // constructors and destructor
0067 //
0068 IPAnalyzer::IPAnalyzer(const edm::ParameterSet& iConfig) {
0069   m_jets = iConfig.getParameter<edm::InputTag>("jets");
0070   m_assoc = iConfig.getParameter<edm::InputTag>("association");
0071   token_ipassoc = consumes<reco::TrackIPTagInfoCollection>(iConfig.getParameter<edm::InputTag>("ipassociation"));
0072 }
0073 
0074 void IPAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0075   using namespace edm;
0076   using namespace reco;
0077 
0078   Handle<TrackIPTagInfoCollection> ipHandle;
0079   iEvent.getByToken(token_ipassoc, ipHandle);
0080   const TrackIPTagInfoCollection& ip = *(ipHandle.product());
0081   cout << "Found " << ip.size() << " TagInfo" << endl;
0082 
0083   cout << boolalpha;
0084   cout << fixed;
0085 
0086   TrackIPTagInfoCollection::const_iterator it = ip.begin();
0087   for (; it != ip.end(); it++) {
0088     cout << "Jet pt: " << it->jet()->pt() << endl;
0089     cout << "Tot tracks: " << it->tracks().size() << endl;
0090     TrackRefVector selTracks = it->selectedTracks();
0091     int n = selTracks.size();
0092     cout << "Sel tracks: " << n << endl;
0093     // false      cout << " Pt  \t d len \t jet dist \t p3d \t p2d\t ip3d \t ip2d " << endl;
0094     GlobalPoint pv(
0095         it->primaryVertex()->position().x(), it->primaryVertex()->position().y(), it->primaryVertex()->position().z());
0096     cout << pv << " vs " << it->primaryVertex()->position() << endl;
0097     for (int j = 0; j < n; j++) {
0098       btag::TrackIPData data = it->impactParameterData()[j];
0099       cout << selTracks[j]->pt() << "\t";
0100       cout << it->probabilities(0)[j] << "\t";
0101       cout << it->probabilities(1)[j] << "\t";
0102       cout << data.ip3d.value() << "\t";
0103       cout << data.ip3d.significance() << "\t";
0104       cout << data.distanceToJetAxis.value() << "\t";
0105       cout << data.distanceToJetAxis.significance() << "\t";
0106       cout << data.distanceToGhostTrack.value() << "\t";
0107       cout << data.distanceToGhostTrack.significance() << "\t";
0108       cout << data.closestToJetAxis << "\t";
0109       cout << (data.closestToJetAxis - pv).mag() << "\t";
0110       cout << data.closestToGhostTrack << "\t";
0111       cout << (data.closestToGhostTrack - pv).mag() << "\t";
0112       cout << data.ip2d.value() << "\t";
0113       cout << data.ip2d.significance() << endl;
0114     }
0115   }
0116 }
0117 
0118 //define this as a plug-in
0119 DEFINE_FWK_MODULE(IPAnalyzer);