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 //
19 //
20 
32 #include <memory>
33 #include <iostream>
34 #include <iomanip>
35 #include <cmath>
36 
37 
41 
44 
49 
53 
54 using namespace std;
55 using namespace reco;
56 namespace cms{
57 
58  JetVertexAssociation::JetVertexAssociation(const edm::ParameterSet& iConfig): m_algo(iConfig),
59  jet_token(consumes<CaloJetCollection>(edm::InputTag(iConfig.getParameter<std::string>("JET_ALGO")))),
60  track_token(consumes<TrackCollection>(edm::InputTag(iConfig.getParameter<std::string>("TRACK_ALGO")))),
61  vertex_token(consumes<VertexCollection>(edm::InputTag(iConfig.getParameter<std::string>("VERTEX_ALGO")))) {
62 
63 
64 
65  produces<ResultCollection1>("Var");
66  produces<ResultCollection2>("JetType");
67 
68 
69  }
70 
72 
74  iEvent.getByToken(jet_token, jets);
75 
77  iEvent.getByToken(track_token, tracks);
78 
80  iEvent.getByToken(vertex_token, vertexes);
81 
82  double SIGNAL_V_Z = 0.;
83  double SIGNAL_V_Z_ERROR = 0.;
84  double ptmax = -100.;
85 
86  VertexCollection::const_iterator vert = vertexes->begin ();
87  if(vertexes->size() > 0 ) {
88  for (; vert != vertexes->end (); vert++) {
89 
90  SIGNAL_V_Z = vert->z();
91  double pt = 0.;
92  reco::Vertex::trackRef_iterator tr = vert->tracks_begin();
93  for (; tr != vert->tracks_end(); tr++) pt += (*tr)->pt();
94  if( pt >= ptmax ){
95 
96  ptmax = pt;
97  SIGNAL_V_Z = vert->z();
98  SIGNAL_V_Z_ERROR = vert->zError();
99 
100  }
101 
102  }
103  }
104 
105  pair<double, bool> result;
106  std::auto_ptr<ResultCollection1> result1 (new ResultCollection1) ;
107  std::auto_ptr<ResultCollection2> result2 (new ResultCollection2) ;
108 
109  CaloJetCollection::const_iterator jet = jets->begin ();
110 
111  if(jets->size() > 0 ) {
112  for (; jet != jets->end (); jet++) {
113  result = m_algo.Main(*jet, tracks, SIGNAL_V_Z, SIGNAL_V_Z_ERROR);
114  result1->push_back(result.first);
115  result2->push_back(result.second);
116 
117  }
118  }
119 
120  iEvent.put(result1, "Var");
121  iEvent.put(result2, "JetType");
122 
123  }
124 }
std::pair< double, bool > Main(const reco::CaloJet &jet, edm::Handle< reco::TrackCollection > tracks, double SIGNAL_V_Z, double SIGNAL_V_Z_Error)
edm::EDGetTokenT< reco::VertexCollection > vertex_token
void produce(edm::Event &e, const edm::EventSetup &c)
edm::EDGetTokenT< reco::TrackCollection > track_token
std::vector< double > ResultCollection1
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::CaloJetCollection > jet_token
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
tuple result
Definition: query.py:137
std::vector< bool > ResultCollection2
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:37
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects