CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
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  // set up output
76  tree_ = fs_->make<TTree>("RecoTree", "RecoTree");
77  tree_->Branch("Vertex", "L1Analysis::L1AnalysisRecoVertexDataFormat", &vtxData_, 32000, 3);
78 }
79 
81  // do anything here that needs to be done at desctruction time
82  // (e.g. close files, deallocate resources etc.)
83 }
84 
85 //
86 // member functions
87 //
88 
89 // ------------ method called to for each event ------------
91  vtxData_->Reset();
92 
93  // get vertices
95  iEvent.getByToken(vtxToken_, vertices);
96 
97  if (vertices.isValid()) {
98  for (reco::VertexCollection::const_iterator it = vertices->begin();
99  it != vertices->end() && vtxData_->nVtx < maxVtx_;
100  ++it) {
101  if (!it->isFake()) {
102  vtxData_->NDoF.push_back(it->ndof());
103  vtxData_->Z.push_back(it->z());
104  vtxData_->Rho.push_back(it->position().rho());
105  vtxData_->nVtx++;
106  }
107  }
108  tree_->Fill();
109  }
110 }
111 
112 // ------------ method called once each job just before starting event loop ------------
114 
115 // ------------ method called once each job just after ending the event loop ------------
117 
118 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
void analyze(const edm::Event &, const edm::EventSetup &) override
L1Analysis::L1AnalysisRecoVertexDataFormat * vtxData_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
L1RecoTreeProducer(const edm::ParameterSet &)
void endJob() override
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
int iEvent
Definition: GenABIO.cc:224
~L1RecoTreeProducer() override
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void beginJob(void) override
edm::Service< TFileService > fs_