DESY Hbb Analysis Framework
Public Member Functions | Static Public Attributes | Protected Attributes | List of all members
analysis::ntuple::Vertices Class Reference

#include <Vertices.h>

Public Member Functions

void Fill (const edm::Event &)
 
void ReadFromEvent (const edm::Event &)
 
 Vertices ()
 
 Vertices (const edm::InputTag &, TTree *)
 
 ~Vertices ()
 

Static Public Attributes

static const int maxPVs = 100
 

Protected Attributes

reco::VertexCollection candidates_
 
float chi2_ [maxPVs]
 
std::string configParameter_
 
bool fake_ [maxPVs]
 
edm::InputTag input_collection_
 
int n_
 
float ndof_ [maxPVs]
 
float rho_ [maxPVs]
 
TTree * tree_
 
float x_ [maxPVs]
 
float xe_ [maxPVs]
 
float y_ [maxPVs]
 
float ye_ [maxPVs]
 
float z_ [maxPVs]
 
float ze_ [maxPVs]
 

Detailed Description

Definition at line 45 of file Vertices.h.

Constructor & Destructor Documentation

Vertices::Vertices ( )

Definition at line 44 of file Vertices.cc.

45 {
46  // default constructor
47 }
Vertices::Vertices ( const edm::InputTag &  tag,
TTree *  tree 
)

Definition at line 49 of file Vertices.cc.

References chi2_, fake_, input_collection_, n_, ndof_, rho_, tree_, x_, xe_, y_, ye_, z_, and ze_.

50 {
51  input_collection_ = tag;
52  tree_ = tree;
53 
54  tree_->Branch("n", &this->n_, "n/I");
55  tree_->Branch("x", this->x_, "x[n]/F");
56  tree_->Branch("y", this->y_, "y[n]/F");
57  tree_->Branch("z", this->z_, "z[n]/F");
58  tree_->Branch("xe", this->xe_, "xe[n]/F");
59  tree_->Branch("ye", this->ye_, "ye[n]/F");
60  tree_->Branch("ze", this->ze_, "ze[n]/F");
61  tree_->Branch("rho", this->rho_, "rho[n]/F");
62  tree_->Branch("chi2", this->chi2_, "chi2[n]/F");
63  tree_->Branch("ndof", this->ndof_, "ndof[n]/F");
64  tree_->Branch("fake", this->fake_, "fake[n]/O");
65 
66 }
float ndof_[maxPVs]
Definition: Vertices.h:71
edm::InputTag input_collection_
Definition: Vertices.h:58
float rho_[maxPVs]
Definition: Vertices.h:72
float chi2_[maxPVs]
Definition: Vertices.h:70
Vertices::~Vertices ( )

Definition at line 68 of file Vertices.cc.

69 {
70  // do anything here that needs to be done at desctruction time
71  // (e.g. close files, deallocate resources etc.)
72 }

Member Function Documentation

void Vertices::Fill ( const edm::Event &  event)

Definition at line 112 of file Vertices.cc.

References ReadFromEvent(), and tree_.

113 {
114  this->ReadFromEvent(event);
115  tree_->Fill();
116 }
void ReadFromEvent(const edm::Event &)
Definition: Vertices.cc:80
void Vertices::ReadFromEvent ( const edm::Event &  event)

Definition at line 80 of file Vertices.cc.

References candidates_, chi2_, fake_, input_collection_, maxPVs, n_, ndof_, rho_, x_, xe_, y_, ye_, z_, and ze_.

Referenced by Fill().

81 {
82  using namespace edm;
83 
84  edm::Handle<reco::VertexCollection> handler;
85  event.getByLabel(this->input_collection_, handler);
86  candidates_ = *(handler.product());
87 
88  int n = 0;
89  for ( size_t i = 0 ; i < candidates_.size(); ++i )
90  {
91  if ( n >= maxPVs ) break;
92 
93  this->x_[n] = candidates_[i].x();
94  this->y_[n] = candidates_[i].y();
95  this->z_[n] = candidates_[i].z();
96  this->xe_[n] = candidates_[i].xError();
97  this->ye_[n] = candidates_[i].yError();
98  this->ze_[n] = candidates_[i].zError();
99  this->fake_[n] = candidates_[i].isFake();
100  this->ndof_[n] = candidates_[i].ndof();
101  this->chi2_[n] = candidates_[i].chi2();
102  this->rho_[n] = candidates_[i].position().Rho();
103 
104  ++n;
105 
106  }
107  this->n_ = n;
108 
109 
110 }
float ndof_[maxPVs]
Definition: Vertices.h:71
static const int maxPVs
Definition: Vertices.h:52
reco::VertexCollection candidates_
Definition: Vertices.h:56
edm::InputTag input_collection_
Definition: Vertices.h:58
float rho_[maxPVs]
Definition: Vertices.h:72
float chi2_[maxPVs]
Definition: Vertices.h:70

Member Data Documentation

reco::VertexCollection analysis::ntuple::Vertices::candidates_
protected

Definition at line 56 of file Vertices.h.

Referenced by ReadFromEvent().

float analysis::ntuple::Vertices::chi2_[maxPVs]
protected

Definition at line 70 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

std::string analysis::ntuple::Vertices::configParameter_
protected

Definition at line 57 of file Vertices.h.

bool analysis::ntuple::Vertices::fake_[maxPVs]
protected

Definition at line 69 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

edm::InputTag analysis::ntuple::Vertices::input_collection_
protected

Definition at line 58 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

const int analysis::ntuple::Vertices::maxPVs = 100
static

Definition at line 52 of file Vertices.h.

Referenced by ReadFromEvent().

int analysis::ntuple::Vertices::n_
protected

Definition at line 62 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

float analysis::ntuple::Vertices::ndof_[maxPVs]
protected

Definition at line 71 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

float analysis::ntuple::Vertices::rho_[maxPVs]
protected

Definition at line 72 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

TTree* analysis::ntuple::Vertices::tree_
protected

Definition at line 74 of file Vertices.h.

Referenced by Fill(), and Vertices().

float analysis::ntuple::Vertices::x_[maxPVs]
protected

Definition at line 63 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

float analysis::ntuple::Vertices::xe_[maxPVs]
protected

Definition at line 66 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

float analysis::ntuple::Vertices::y_[maxPVs]
protected

Definition at line 64 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

float analysis::ntuple::Vertices::ye_[maxPVs]
protected

Definition at line 67 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

float analysis::ntuple::Vertices::z_[maxPVs]
protected

Definition at line 65 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().

float analysis::ntuple::Vertices::ze_[maxPVs]
protected

Definition at line 68 of file Vertices.h.

Referenced by ReadFromEvent(), and Vertices().


The documentation for this class was generated from the following files: