Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:14

0001 // -*- C++ -*-
0002 //
0003 // Package:    L1Trigger/L1TNtuples
0004 // Class:      L1RecoTreeProducer
0005 //
0006 /**\class L1RecoTreeProducer L1RecoTreeProducer.cc L1Trigger/L1TNtuples/src/L1RecoTreeProducer.cc
0007  Description: Produces tree containing reco quantities
0008 */
0009 
0010 // system include files
0011 #include <memory>
0012 
0013 // framework
0014 #include "FWCore/Framework/interface/Frameworkfwd.h"
0015 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0016 #include "FWCore/Framework/interface/Event.h"
0017 #include "FWCore/Framework/interface/EventSetup.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0021 #include "FWCore/Framework/interface/ESHandle.h"
0022 
0023 // ROOT output stuff
0024 #include "FWCore/ServiceRegistry/interface/Service.h"
0025 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0026 #include "TH1.h"
0027 #include "TTree.h"
0028 #include "TF1.h"
0029 
0030 // EDM formats
0031 #include "DataFormats/VertexReco/interface/Vertex.h"
0032 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0033 
0034 //local  data formats
0035 #include "L1Trigger/L1TNtuples/interface/L1AnalysisRecoVertexDataFormat.h"
0036 
0037 //
0038 // class declaration
0039 //
0040 
0041 class L1RecoTreeProducer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0042 public:
0043   explicit L1RecoTreeProducer(const edm::ParameterSet&);
0044   ~L1RecoTreeProducer() override;
0045 
0046 private:
0047   void beginJob(void) override;
0048   void analyze(const edm::Event&, const edm::EventSetup&) override;
0049   void endJob() override;
0050 
0051 public:
0052   L1Analysis::L1AnalysisRecoVertexDataFormat* vtxData_;
0053 
0054 private:
0055   // output file
0056   edm::Service<TFileService> fs_;
0057 
0058   // tree
0059   TTree* tree_;
0060 
0061   // EDM input tags
0062   edm::EDGetTokenT<reco::VertexCollection> vtxToken_;
0063 
0064   unsigned int maxVtx_;
0065 };
0066 
0067 L1RecoTreeProducer::L1RecoTreeProducer(const edm::ParameterSet& iConfig) {
0068   vtxToken_ = consumes<reco::VertexCollection>(
0069       iConfig.getUntrackedParameter("vtxToken", edm::InputTag("offlinePrimaryVertices")));
0070 
0071   maxVtx_ = iConfig.getParameter<unsigned int>("maxVtx");
0072 
0073   vtxData_ = new L1Analysis::L1AnalysisRecoVertexDataFormat();
0074 
0075   usesResource(TFileService::kSharedResource);
0076   // set up output
0077   tree_ = fs_->make<TTree>("RecoTree", "RecoTree");
0078   tree_->Branch("Vertex", "L1Analysis::L1AnalysisRecoVertexDataFormat", &vtxData_, 32000, 3);
0079 }
0080 
0081 L1RecoTreeProducer::~L1RecoTreeProducer() {
0082   // do anything here that needs to be done at desctruction time
0083   // (e.g. close files, deallocate resources etc.)
0084 }
0085 
0086 //
0087 // member functions
0088 //
0089 
0090 // ------------ method called to for each event  ------------
0091 void L1RecoTreeProducer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0092   vtxData_->Reset();
0093 
0094   // get vertices
0095   edm::Handle<reco::VertexCollection> vertices;
0096   iEvent.getByToken(vtxToken_, vertices);
0097 
0098   if (vertices.isValid()) {
0099     for (reco::VertexCollection::const_iterator it = vertices->begin();
0100          it != vertices->end() && vtxData_->nVtx < maxVtx_;
0101          ++it) {
0102       if (!it->isFake()) {
0103         vtxData_->NDoF.push_back(it->ndof());
0104         vtxData_->Z.push_back(it->z());
0105         vtxData_->Rho.push_back(it->position().rho());
0106         vtxData_->nVtx++;
0107       }
0108     }
0109     tree_->Fill();
0110   }
0111 }
0112 
0113 // ------------ method called once each job just before starting event loop  ------------
0114 void L1RecoTreeProducer::beginJob(void) {}
0115 
0116 // ------------ method called once each job just after ending the event loop  ------------
0117 void L1RecoTreeProducer::endJob() {}
0118 
0119 //define this as a plug-in
0120 DEFINE_FWK_MODULE(L1RecoTreeProducer);