CMS 3D CMS Logo

L1RecoTreeProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1Trigger/L1TNtuples
4 // Class: L1RecoTreeProducer
5 //
10 // system include files
11 #include <memory>
12 
13 // framework
22 
23 // ROOT output stuff
26 #include "TH1.h"
27 #include "TTree.h"
28 #include "TF1.h"
29 
30 // EDM formats
33 
34 //local data formats
36 
37 //
38 // class declaration
39 //
40 
41 class L1RecoTreeProducer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
42 public:
43  explicit L1RecoTreeProducer(const edm::ParameterSet&);
44  ~L1RecoTreeProducer() override;
45 
46 private:
47  void beginJob(void) override;
48  void analyze(const edm::Event&, const edm::EventSetup&) override;
49  void endJob() override;
50 
51 public:
53 
54 private:
55  // output file
57 
58  // tree
59  TTree* tree_;
60 
61  // EDM input tags
63 
64  unsigned int maxVtx_;
65 };
66 
68  vtxToken_ = consumes<reco::VertexCollection>(
69  iConfig.getUntrackedParameter("vtxToken", edm::InputTag("offlinePrimaryVertices")));
70 
71  maxVtx_ = iConfig.getParameter<unsigned int>("maxVtx");
72 
74 
75  usesResource(TFileService::kSharedResource);
76  // set up output
77  tree_ = fs_->make<TTree>("RecoTree", "RecoTree");
78  tree_->Branch("Vertex", "L1Analysis::L1AnalysisRecoVertexDataFormat", &vtxData_, 32000, 3);
79 }
80 
82  // do anything here that needs to be done at desctruction time
83  // (e.g. close files, deallocate resources etc.)
84 }
85 
86 //
87 // member functions
88 //
89 
90 // ------------ method called to for each event ------------
92  vtxData_->Reset();
93 
94  // get vertices
96  iEvent.getByToken(vtxToken_, vertices);
97 
98  if (vertices.isValid()) {
99  for (reco::VertexCollection::const_iterator it = vertices->begin();
100  it != vertices->end() && vtxData_->nVtx < maxVtx_;
101  ++it) {
102  if (!it->isFake()) {
103  vtxData_->NDoF.push_back(it->ndof());
104  vtxData_->Z.push_back(it->z());
105  vtxData_->Rho.push_back(it->position().rho());
106  vtxData_->nVtx++;
107  }
108  }
109  tree_->Fill();
110  }
111 }
112 
113 // ------------ method called once each job just before starting event loop ------------
115 
116 // ------------ method called once each job just after ending the event loop ------------
118 
119 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void analyze(const edm::Event &, const edm::EventSetup &) override
L1Analysis::L1AnalysisRecoVertexDataFormat * vtxData_
L1RecoTreeProducer(const edm::ParameterSet &)
void endJob() override
T getUntrackedParameter(std::string const &, T const &) const
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~L1RecoTreeProducer() override
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
void beginJob(void) override
edm::Service< TFileService > fs_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64