CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetVertexAssociation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: JetVertexAssociation
4 // Class: JetVertexAssociation
5 //
13 //
14 // Original Author: Natalia Ilina
15 // Modified by Eduardo Luiggi
16 //
17 // Created: Tue Oct 31 10:52:41 CET 2006
18 // $Id: JetVertexAssociation.cc,v 1.8 2011/05/20 17:17:28 wmtan Exp $
19 //
20 //
21 
33 #include <memory>
34 #include <iostream>
35 #include <iomanip>
36 #include <cmath>
37 
38 
42 
45 
53 
57 
58 using namespace std;
59 using namespace reco;
60 namespace cms{
61 
62  JetVertexAssociation::JetVertexAssociation(const edm::ParameterSet& iConfig): m_algo(iConfig),
63  jet_algo(iConfig.getParameter<std::string>("JET_ALGO")),
64  track_algo(iConfig.getParameter<std::string>("TRACK_ALGO")),
65  vertex_algo(iConfig.getParameter<std::string>("VERTEX_ALGO")) {
66 
67 
68 
69  produces<ResultCollection1>("Var");
70  produces<ResultCollection2>("JetType");
71 
72 
73  }
74 
76 
78  iEvent.getByLabel(jet_algo, jets);
79 
81  iEvent.getByLabel(track_algo, tracks);
82 
84  iEvent.getByLabel(vertex_algo, vertexes);
85 
86  double SIGNAL_V_Z = 0.;
87  double SIGNAL_V_Z_ERROR = 0.;
88  double ptmax = -100.;
89 
90  VertexCollection::const_iterator vert = vertexes->begin ();
91  if(vertexes->size() > 0 ) {
92  for (; vert != vertexes->end (); vert++) {
93 
94  SIGNAL_V_Z = vert->z();
95  double pt = 0.;
96  reco::Vertex::trackRef_iterator tr = vert->tracks_begin();
97  for (; tr != vert->tracks_end(); tr++) pt += (*tr)->pt();
98  if( pt >= ptmax ){
99 
100  ptmax = pt;
101  SIGNAL_V_Z = vert->z();
102  SIGNAL_V_Z_ERROR = vert->zError();
103 
104  }
105 
106  }
107  }
108 
109  pair<double, bool> result;
110  std::auto_ptr<ResultCollection1> result1 (new ResultCollection1) ;
111  std::auto_ptr<ResultCollection2> result2 (new ResultCollection2) ;
112 
113  CaloJetCollection::const_iterator jet = jets->begin ();
114 
115  if(jets->size() > 0 ) {
116  for (; jet != jets->end (); jet++) {
117  result = m_algo.Main(*jet, tracks, SIGNAL_V_Z, SIGNAL_V_Z_ERROR);
118  result1->push_back(result.first);
119  result2->push_back(result.second);
120 
121  }
122  }
123 
124  iEvent.put(result1, "Var");
125  iEvent.put(result2, "JetType");
126 
127  }
128 }
std::pair< double, bool > Main(const reco::CaloJet &jet, edm::Handle< reco::TrackCollection > tracks, double SIGNAL_V_Z, double SIGNAL_V_Z_Error)
void produce(edm::Event &e, const edm::EventSetup &c)
std::vector< double > ResultCollection1
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
vector< PseudoJet > jets
tuple result
Definition: query.py:137
std::vector< bool > ResultCollection2
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple tracks
Definition: testEve_cfg.py:39
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38