Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-31 23:02:24

0001 #include "FWCore/Framework/interface/MakerMacros.h"
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "DataFormats/Common/interface/Handle.h"
0007 #include "FWCore/ServiceRegistry/interface/Service.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 
0011 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0012 #include "DataFormats/TrackReco/interface/Track.h"
0013 
0014 #include "SimDataFormats/Track/interface/SimTrack.h"
0015 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
0016 
0017 #include <iostream>
0018 
0019 #include "TH1F.h"
0020 #include "TFile.h"
0021 #include "TObjArray.h"
0022 
0023 using namespace std;
0024 
0025 class PixelTrackTest : public edm::one::EDAnalyzer<> {
0026 public:
0027   explicit PixelTrackTest(const edm::ParameterSet& conf);
0028   ~PixelTrackTest();
0029   virtual void beginJob() {}
0030   virtual void analyze(const edm::Event& ev, const edm::EventSetup& es);
0031   virtual void endJob() {}
0032 
0033 private:
0034   string collectionLabel;
0035   void myprint(const reco::Track& track) const;
0036   TObjArray hList;
0037 };
0038 
0039 PixelTrackTest::PixelTrackTest(const edm::ParameterSet& conf) {
0040   hList.Add(new TH1F("h_Pt", "h_Pt", 100, -1.2, 1.2));
0041   hList.Add(new TH1F("h_tip", "h_tip", 100, -0.03, 0.03));
0042   hList.SetOwner();
0043   collectionLabel = conf.getParameter<std::string>("TrackCollection");
0044   edm::LogInfo("PixelTrackTest") << " CTOR";
0045 }
0046 
0047 PixelTrackTest::~PixelTrackTest() {
0048   std::string fileName = collectionLabel + ".root";
0049   TFile f(fileName.c_str(), "RECREATE");
0050   hList.Write();
0051   f.Close();
0052 
0053   edm::LogInfo("PixelTrackTest") << " DTOR";
0054 }
0055 
0056 void PixelTrackTest::analyze(const edm::Event& ev, const edm::EventSetup& es) {
0057   using namespace edm;
0058   using namespace std;
0059   using namespace reco;
0060 
0061   edm::LogPrint("PixelTrackTest") << "*** PixelTrackTest, analyze event: " << ev.id();
0062   Handle<SimTrackContainer> simTrks;
0063   ev.getByLabel("g4SimHits", simTrks);
0064   //  std::edm::LogPrint("PixelTrackTest") << "simtrks " << simTrks->size() << std::endl;
0065 
0066   float pt_gen = 0.0;
0067   typedef SimTrackContainer::const_iterator IP;
0068   for (IP p = simTrks->begin(); p != simTrks->end(); p++) {
0069     if ((*p).noVertex())
0070       continue;
0071     if ((*p).type() == -99)
0072       continue;
0073     if ((*p).vertIndex() != 0)
0074       continue;
0075     if ((*p).momentum().Pt() > pt_gen)
0076       pt_gen = (*p).momentum().Pt();
0077   }
0078   //  edm::LogPrint("PixelTrackTest") << "pt_gen: " << pt_gen ;
0079   if (pt_gen < 0.9)
0080     return;
0081 
0082   typedef reco::TrackCollection::const_iterator IT;
0083 
0084   edm::Handle<reco::TrackCollection> trackCollection;
0085   ev.getByLabel(collectionLabel, trackCollection);
0086   const reco::TrackCollection tracks = *(trackCollection.product());
0087   edm::LogPrint("PixelTrackTest") << "Number of tracks: " << tracks.size() << " tracks" << std::endl;
0088   for (IT it = tracks.begin(); it != tracks.end(); it++) {
0089     //math::XYZVector mom_rec = (*it).momentum();
0090     float pt_rec = (*it).pt();
0091     myprint(*it);
0092     static_cast<TH1*>(hList.FindObject("h_Pt"))->Fill((pt_gen - pt_rec) / pt_gen);
0093     static_cast<TH1*>(hList.FindObject("h_tip"))->Fill((*it).d0());
0094   }
0095 
0096   edm::LogPrint("PixelTrackTest") << "------------------------------------------------";
0097 }
0098 
0099 void PixelTrackTest::myprint(const reco::Track& track) const {
0100   edm::LogPrint("PixelTrackTest") << "--- RECONSTRUCTED TRACK: ";
0101   edm::LogPrint("PixelTrackTest") << "\tmomentum: " << track.momentum() << "\tPT: " << track.pt();
0102   edm::LogPrint("PixelTrackTest") << "\tvertex: " << track.vertex() << "\t zip: " << track.dz() << "+/-"
0103                                   << track.dzError() << "\t tip: " << track.d0() << "+/-" << track.d0Error();
0104   edm::LogPrint("PixelTrackTest") << "\t chi2: " << track.chi2();
0105   edm::LogPrint("PixelTrackTest") << "\tcharge: " << track.charge();
0106 }
0107 
0108 DEFINE_FWK_MODULE(PixelTrackTest);