Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    JetVertexAssociation
0004 // Class:      JetVertexAssociation
0005 //
0006 /**\class JetVertexAssociation JetVertexAssociation.cc JetMETCorrections/JetVertexAssociation/src/JetVertexAssociation.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Natalia Ilina
0015 // Modified by Eduardo Luiggi
0016 //
0017 //         Created:  Tue Oct 31 10:52:41 CET 2006
0018 //
0019 //
0020 
0021 /**
0022   * 'JetVertexAssociation' represents the association of the jet with the signal vertex
0023   *
0024   * Parameters of the method: JV_deltaZ, JV_alpha_threshold(alpha_0 or beta_0),
0025   *                           JV_cone_size, JV_type_Algo ("1" - alpha, "2" - beta) - (the details are in CMS NOTE 2006/091),
0026   *
0027   * Output: <pair<double, bool> >.
0028   *                    The first - variable alpha(beta) for the jet,
0029   *                    the second - "true" for jet from signal vertex, "false" for jet from pile-up.
0030   **/
0031 
0032 #include <memory>
0033 #include <iostream>
0034 #include <iomanip>
0035 #include <cmath>
0036 
0037 #include "FWCore/Framework/interface/Frameworkfwd.h"
0038 #include "FWCore/Framework/interface/Event.h"
0039 
0040 #include "FWCore/PluginManager/interface/ModuleDef.h"
0041 #include "FWCore/Framework/interface/MakerMacros.h"
0042 
0043 #include "DataFormats/Common/interface/Handle.h"
0044 #include "FWCore/Framework/interface/EventSetup.h"
0045 #include "DataFormats/JetReco/interface/CaloJet.h"
0046 #include "DataFormats/VertexReco/interface/Vertex.h"
0047 
0048 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0049 #include "JetMETCorrections/JetVertexAssociation/interface/JetVertexAssociation.h"
0050 #include "JetMETCorrections/JetVertexAssociation/interface/JetVertexMain.h"
0051 
0052 using namespace std;
0053 using namespace reco;
0054 namespace cms {
0055 
0056   JetVertexAssociation::JetVertexAssociation(const edm::ParameterSet& iConfig)
0057       : m_algo(iConfig),
0058         jet_token(consumes<CaloJetCollection>(edm::InputTag(iConfig.getParameter<std::string>("JET_ALGO")))),
0059         track_token(consumes<TrackCollection>(edm::InputTag(iConfig.getParameter<std::string>("TRACK_ALGO")))),
0060         vertex_token(consumes<VertexCollection>(edm::InputTag(iConfig.getParameter<std::string>("VERTEX_ALGO")))) {
0061     produces<ResultCollection1>("Var");
0062     produces<ResultCollection2>("JetType");
0063   }
0064 
0065   void JetVertexAssociation::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0066     edm::Handle<CaloJetCollection> jets;
0067     iEvent.getByToken(jet_token, jets);
0068 
0069     edm::Handle<TrackCollection> tracks;
0070     iEvent.getByToken(track_token, tracks);
0071 
0072     edm::Handle<VertexCollection> vertexes;
0073     iEvent.getByToken(vertex_token, vertexes);
0074 
0075     double SIGNAL_V_Z = 0.;
0076     double SIGNAL_V_Z_ERROR = 0.;
0077     double ptmax = -100.;
0078 
0079     VertexCollection::const_iterator vert = vertexes->begin();
0080     if (!vertexes->empty()) {
0081       for (; vert != vertexes->end(); vert++) {
0082         SIGNAL_V_Z = vert->z();
0083         double pt = 0.;
0084         reco::Vertex::trackRef_iterator tr = vert->tracks_begin();
0085         for (; tr != vert->tracks_end(); tr++)
0086           pt += (*tr)->pt();
0087         if (pt >= ptmax) {
0088           ptmax = pt;
0089           SIGNAL_V_Z = vert->z();
0090           SIGNAL_V_Z_ERROR = vert->zError();
0091         }
0092       }
0093     }
0094 
0095     pair<double, bool> result;
0096     std::unique_ptr<ResultCollection1> result1(new ResultCollection1);
0097     std::unique_ptr<ResultCollection2> result2(new ResultCollection2);
0098 
0099     CaloJetCollection::const_iterator jet = jets->begin();
0100 
0101     if (!jets->empty()) {
0102       for (; jet != jets->end(); jet++) {
0103         result = m_algo.Main(*jet, tracks, SIGNAL_V_Z, SIGNAL_V_Z_ERROR);
0104         result1->push_back(result.first);
0105         result2->push_back(result.second);
0106       }
0107     }
0108 
0109     iEvent.put(std::move(result1), "Var");
0110     iEvent.put(std::move(result2), "JetType");
0111   }
0112 }  // namespace cms