CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RECOVertex.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <sstream>
3 #include <istream>
4 #include <fstream>
5 #include <iomanip>
6 #include <stdlib.h>
7 #include <string.h>
8 
10 #include "HLTMessages.h"
11 
12 static const size_t kMaxVrt = 50;
13 
15 
16  //set parameter defaults
17  _Debug=false;
18 
19  NVrtx = 0;
20  VertexCand_x = new float[kMaxVrt];
21  VertexCand_y = new float[kMaxVrt];
22  VertexCand_z = new float[kMaxVrt];
23  VertexCand_tracks = new int[kMaxVrt];
24  VertexCand_chi2 = new float[kMaxVrt];
25  VertexCand_ndof = new float[kMaxVrt];
26 
27 }
28 
30 
31 }
32 
34 {
35  NVrtx = 0;
36  std::memset(VertexCand_x, '\0', kMaxVrt * sizeof(float));
37  std::memset(VertexCand_y, '\0', kMaxVrt * sizeof(float));
38  std::memset(VertexCand_z, '\0', kMaxVrt * sizeof(float));
39  std::memset(VertexCand_tracks, '\0', kMaxVrt * sizeof(int));
40  std::memset(VertexCand_chi2, '\0', kMaxVrt * sizeof(float));
41  std::memset(VertexCand_ndof, '\0', kMaxVrt * sizeof(float));
42 }
43 
44 /* Setup the analysis to put the branch-variables into the tree. */
45 void RECOVertex::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
46 
47  edm::ParameterSet myHltParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
48  std::vector<std::string> parameterNames = myHltParams.getParameterNames() ;
49 
50  for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
51  iParam != parameterNames.end(); iParam++ ){
52  if ( (*iParam) == "Debug" ) _Debug = myHltParams.getParameter<bool>( *iParam );
53  }
54 
55  HltTree->Branch("recoNVrt", & NVrtx, "NVrtx/I");
56  HltTree->Branch("recoVrtX", VertexCand_x, "recoVrtX[NVrtx]/F");
57  HltTree->Branch("recoVrtY", VertexCand_y, "recoVrtY[NVrtx]/F");
58  HltTree->Branch("recoVrtZ", VertexCand_z, "recoVrtZ[NVrtx]/F");
59  HltTree->Branch("recoVrtNtrk", VertexCand_tracks, "recoVrtNtrk[NVrtx]/I");
60  HltTree->Branch("recoVrtChi2", VertexCand_chi2, "recoVrtChi2[NVrtx]/F");
61  HltTree->Branch("recoVrtNdof", VertexCand_ndof, "recoVrtNdof[NVrtx]/F");
62 
63 }
64 
65 /* **Analyze the event** */
66 void RECOVertex::analyze(edm::Handle<reco::VertexCollection> recoVertexs, TTree* HltTree) {
67 
68  // reset the tree variables
69  clear();
70 
71  if ( recoVertexs.isValid() ) {
72  const reco::VertexCollection* vertexs = recoVertexs.product();
73  reco::VertexCollection::const_iterator vertex_i;
74 
75  size_t size = std::min(kMaxVrt, size_t(vertexs->size()) );
76  NVrtx= size;
77 
78  int nVertexCand=0;
79  if (_Debug) std::cout << "Found " << vertexs->size() << " vertices" << std::endl;
80  for (vertex_i = vertexs->begin(); vertex_i != vertexs->end(); vertex_i++){
81  if (nVertexCand>=NVrtx) break;
82  VertexCand_x[nVertexCand] = vertex_i->x();
83  VertexCand_y[nVertexCand] = vertex_i->y();
84  VertexCand_z[nVertexCand] = vertex_i->z();
85  VertexCand_tracks[nVertexCand] = vertex_i->tracksSize();
86  VertexCand_chi2[nVertexCand] = vertex_i->chi2();
87  VertexCand_ndof[nVertexCand] = vertex_i->ndof();
88  if (_Debug) {
89  std::cout << "RECOVertex -- VX, VY VZ = "
90  << VertexCand_x[nVertexCand] << " "
91  << VertexCand_y[nVertexCand] << " "
92  << VertexCand_z[nVertexCand]
93  << std::endl;
94  std::cout << "RECOVertex -- Ntracks, Chi2/Dof = "
95  << VertexCand_tracks[nVertexCand] << " "
96  << VertexCand_chi2[nVertexCand] << " / " << VertexCand_ndof[nVertexCand]
97  << std::endl;
98  }
99  nVertexCand++;
100 
101  }
102  }
103 }
104 
T getParameter(std::string const &) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
float * VertexCand_y
Definition: RECOVertex.h:34
#define min(a, b)
Definition: mlp_lapack.h:161
void setup(const edm::ParameterSet &pSet, TTree *tree)
Definition: RECOVertex.cc:45
static const size_t kMaxVrt
Definition: RECOVertex.cc:12
float * VertexCand_ndof
Definition: RECOVertex.h:37
void analyze(edm::Handle< reco::VertexCollection > recoVertexs, TTree *tree)
Definition: RECOVertex.cc:66
std::vector< std::string > getParameterNames() const
bool isValid() const
Definition: HandleBase.h:76
float * VertexCand_z
Definition: RECOVertex.h:34
bool _Debug
Definition: RECOVertex.h:40
T const * product() const
Definition: Handle.h:74
float * VertexCand_chi2
Definition: RECOVertex.h:36
int * VertexCand_tracks
Definition: RECOVertex.h:35
float * VertexCand_x
Definition: RECOVertex.h:34
tuple cout
Definition: gather_cfg.py:41
void clear(void)
Definition: RECOVertex.cc:33
tuple size
Write out results.