Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:23

0001 #include <iostream>
0002 #include <sstream>
0003 #include <istream>
0004 #include <fstream>
0005 #include <iomanip>
0006 #include <cstdlib>
0007 #include <cstring>
0008 
0009 #include "RECOVertex.h"
0010 #include "HLTMessages.h"
0011 
0012 static const size_t kMaxVrt = 50;
0013 
0014 RECOVertex::RECOVertex() {
0015   //set parameter defaults
0016   _Debug = false;
0017 
0018   NVrtx = 0;
0019   VertexCand_x = new float[kMaxVrt];
0020   VertexCand_y = new float[kMaxVrt];
0021   VertexCand_z = new float[kMaxVrt];
0022   VertexCand_tracks = new int[kMaxVrt];
0023   VertexCand_chi2 = new float[kMaxVrt];
0024   VertexCand_ndof = new float[kMaxVrt];
0025 }
0026 
0027 RECOVertex::~RECOVertex() {
0028   delete[] VertexCand_x;
0029   delete[] VertexCand_y;
0030   delete[] VertexCand_z;
0031   delete[] VertexCand_tracks;
0032   delete[] VertexCand_chi2;
0033   delete[] VertexCand_ndof;
0034 }
0035 
0036 void RECOVertex::clear() {
0037   NVrtx = 0;
0038   std::memset(VertexCand_x, '\0', kMaxVrt * sizeof(float));
0039   std::memset(VertexCand_y, '\0', kMaxVrt * sizeof(float));
0040   std::memset(VertexCand_z, '\0', kMaxVrt * sizeof(float));
0041   std::memset(VertexCand_tracks, '\0', kMaxVrt * sizeof(int));
0042   std::memset(VertexCand_chi2, '\0', kMaxVrt * sizeof(float));
0043   std::memset(VertexCand_ndof, '\0', kMaxVrt * sizeof(float));
0044 }
0045 
0046 /*  Setup the analysis to put the branch-variables into the tree. */
0047 void RECOVertex::setup(const edm::ParameterSet& pSet, TTree* HltTree, std::string vertexType) {
0048   edm::ParameterSet myHltParams = pSet.getParameter<edm::ParameterSet>("RunParameters");
0049   std::vector<std::string> parameterNames = myHltParams.getParameterNames();
0050 
0051   for (auto& parameterName : parameterNames) {
0052     if (parameterName == "Debug")
0053       _Debug = myHltParams.getParameter<bool>(parameterName);
0054   }
0055 
0056   TString br_recoNVrt = "recoNVrt";
0057   br_recoNVrt.Append(vertexType);
0058   HltTree->Branch(br_recoNVrt, &NVrtx, "NVrtx/I");
0059 
0060   TString br_recoVrtX = "recoVrtX";
0061   br_recoVrtX.Append(vertexType);
0062   HltTree->Branch(br_recoVrtX, VertexCand_x, "recoVrtX[NVrtx]/F");
0063 
0064   TString br_recoVrtY = "recoVrtY";
0065   br_recoVrtY.Append(vertexType);
0066   HltTree->Branch(br_recoVrtY, VertexCand_y, "recoVrtY[NVrtx]/F");
0067 
0068   TString br_recoVrtZ = "recoVrtZ";
0069   br_recoVrtZ.Append(vertexType);
0070   HltTree->Branch(br_recoVrtZ, VertexCand_z, "recoVrtZ[NVrtx]/F");
0071 
0072   TString br_recoVrtNtrk = "recoVrtNtrk";
0073   br_recoVrtNtrk.Append(vertexType);
0074   HltTree->Branch(br_recoVrtNtrk, VertexCand_tracks, "recoVrtNtrk[NVrtx]/I");
0075 
0076   TString br_recoVrtChi2 = "recoVrtChi2";
0077   br_recoVrtChi2.Append(vertexType);
0078   HltTree->Branch(br_recoVrtChi2, VertexCand_chi2, "recoVrtChi2[NVrtx]/F");
0079 
0080   TString br_recoVrtNdof = "recoVrtNdof";
0081   br_recoVrtNdof.Append(vertexType);
0082   HltTree->Branch(br_recoVrtNdof, VertexCand_ndof, "recoVrtNdof[NVrtx]/F");
0083 }
0084 
0085 /* **Analyze the event** */
0086 void RECOVertex::analyze(edm::Handle<reco::VertexCollection> recoVertexs, TTree* HltTree) {
0087   // reset the tree variables
0088   clear();
0089 
0090   if (recoVertexs.isValid()) {
0091     const reco::VertexCollection* vertexs = recoVertexs.product();
0092     reco::VertexCollection::const_iterator vertex_i;
0093 
0094     size_t size = std::min(kMaxVrt, size_t(vertexs->size()));
0095     NVrtx = size;
0096 
0097     int nVertexCand = 0;
0098     if (_Debug)
0099       std::cout << "Found " << vertexs->size() << " vertices" << std::endl;
0100     for (vertex_i = vertexs->begin(); vertex_i != vertexs->end(); vertex_i++) {
0101       if (nVertexCand >= NVrtx)
0102         break;
0103       VertexCand_x[nVertexCand] = vertex_i->x();
0104       VertexCand_y[nVertexCand] = vertex_i->y();
0105       VertexCand_z[nVertexCand] = vertex_i->z();
0106       VertexCand_tracks[nVertexCand] = vertex_i->tracksSize();
0107       VertexCand_chi2[nVertexCand] = vertex_i->chi2();
0108       VertexCand_ndof[nVertexCand] = vertex_i->ndof();
0109       if (_Debug) {
0110         std::cout << "RECOVertex -- VX, VY VZ   = " << VertexCand_x[nVertexCand] << " " << VertexCand_y[nVertexCand]
0111                   << " " << VertexCand_z[nVertexCand] << std::endl;
0112         std::cout << "RECOVertex -- Ntracks, Chi2/Dof   = " << VertexCand_tracks[nVertexCand] << " "
0113                   << VertexCand_chi2[nVertexCand] << " / " << VertexCand_ndof[nVertexCand] << std::endl;
0114       }
0115       nVertexCand++;
0116     }
0117   }
0118 }