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 
29  delete[] VertexCand_x;
30  delete[] VertexCand_y;
31  delete[] VertexCand_z;
32  delete[] VertexCand_tracks;
33  delete[] VertexCand_chi2;
34  delete[] VertexCand_ndof;
35 }
36 
38 {
39  NVrtx = 0;
40  std::memset(VertexCand_x, '\0', kMaxVrt * sizeof(float));
41  std::memset(VertexCand_y, '\0', kMaxVrt * sizeof(float));
42  std::memset(VertexCand_z, '\0', kMaxVrt * sizeof(float));
43  std::memset(VertexCand_tracks, '\0', kMaxVrt * sizeof(int));
44  std::memset(VertexCand_chi2, '\0', kMaxVrt * sizeof(float));
45  std::memset(VertexCand_ndof, '\0', kMaxVrt * sizeof(float));
46 }
47 
48 /* Setup the analysis to put the branch-variables into the tree. */
49 void RECOVertex::setup(const edm::ParameterSet& pSet, TTree* HltTree, std::string vertexType) {
50 
51  edm::ParameterSet myHltParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
52  std::vector<std::string> parameterNames = myHltParams.getParameterNames() ;
53 
54  for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
55  iParam != parameterNames.end(); iParam++ ){
56  if ( (*iParam) == "Debug" ) _Debug = myHltParams.getParameter<bool>( *iParam );
57  }
58 
59  TString br_recoNVrt = "recoNVrt";
60  br_recoNVrt.Append(vertexType);
61  HltTree->Branch(br_recoNVrt, & NVrtx, "NVrtx/I");
62 
63  TString br_recoVrtX = "recoVrtX";
64  br_recoVrtX.Append(vertexType);
65  HltTree->Branch(br_recoVrtX, VertexCand_x, "recoVrtX[NVrtx]/F");
66 
67  TString br_recoVrtY = "recoVrtY";
68  br_recoVrtY.Append(vertexType);
69  HltTree->Branch(br_recoVrtY, VertexCand_y, "recoVrtY[NVrtx]/F");
70 
71  TString br_recoVrtZ = "recoVrtZ";
72  br_recoVrtZ.Append(vertexType);
73  HltTree->Branch(br_recoVrtZ, VertexCand_z, "recoVrtZ[NVrtx]/F");
74 
75  TString br_recoVrtNtrk = "recoVrtNtrk";
76  br_recoVrtNtrk.Append(vertexType);
77  HltTree->Branch(br_recoVrtNtrk, VertexCand_tracks, "recoVrtNtrk[NVrtx]/I");
78 
79  TString br_recoVrtChi2 = "recoVrtChi2";
80  br_recoVrtChi2.Append(vertexType);
81  HltTree->Branch(br_recoVrtChi2, VertexCand_chi2, "recoVrtChi2[NVrtx]/F");
82 
83  TString br_recoVrtNdof = "recoVrtNdof";
84  br_recoVrtNdof.Append(vertexType);
85  HltTree->Branch(br_recoVrtNdof, VertexCand_ndof, "recoVrtNdof[NVrtx]/F");
86 
87 
88 
89 }
90 
91 /* **Analyze the event** */
92 void RECOVertex::analyze(edm::Handle<reco::VertexCollection> recoVertexs, TTree* HltTree) {
93 
94  // reset the tree variables
95  clear();
96 
97  if ( recoVertexs.isValid() ) {
98  const reco::VertexCollection* vertexs = recoVertexs.product();
99  reco::VertexCollection::const_iterator vertex_i;
100 
101  size_t size = std::min(kMaxVrt, size_t(vertexs->size()) );
102  NVrtx= size;
103 
104  int nVertexCand=0;
105  if (_Debug) std::cout << "Found " << vertexs->size() << " vertices" << std::endl;
106  for (vertex_i = vertexs->begin(); vertex_i != vertexs->end(); vertex_i++){
107  if (nVertexCand>=NVrtx) break;
108  VertexCand_x[nVertexCand] = vertex_i->x();
109  VertexCand_y[nVertexCand] = vertex_i->y();
110  VertexCand_z[nVertexCand] = vertex_i->z();
111  VertexCand_tracks[nVertexCand] = vertex_i->tracksSize();
112  VertexCand_chi2[nVertexCand] = vertex_i->chi2();
113  VertexCand_ndof[nVertexCand] = vertex_i->ndof();
114  if (_Debug) {
115  std::cout << "RECOVertex -- VX, VY VZ = "
116  << VertexCand_x[nVertexCand] << " "
117  << VertexCand_y[nVertexCand] << " "
118  << VertexCand_z[nVertexCand]
119  << std::endl;
120  std::cout << "RECOVertex -- Ntracks, Chi2/Dof = "
121  << VertexCand_tracks[nVertexCand] << " "
122  << VertexCand_chi2[nVertexCand] << " / " << VertexCand_ndof[nVertexCand]
123  << std::endl;
124  }
125  nVertexCand++;
126 
127  }
128  }
129 }
130 
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
static const size_t kMaxVrt
Definition: RECOVertex.cc:12
void setup(const edm::ParameterSet &pSet, TTree *tree, std::string vertexType)
Definition: RECOVertex.cc:49
float * VertexCand_ndof
Definition: RECOVertex.h:37
void analyze(edm::Handle< reco::VertexCollection > recoVertexs, TTree *tree)
Definition: RECOVertex.cc:92
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< std::string > getParameterNames() const
bool isValid() const
Definition: HandleBase.h:75
float * VertexCand_z
Definition: RECOVertex.h:34
bool _Debug
Definition: RECOVertex.h:40
T const * product() const
Definition: Handle.h:81
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:145
void clear(void)
Definition: RECOVertex.cc:37
tuple size
Write out results.