Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:18:27

0001 // \class JetTracksAssociatorAtCaloFace JetTracksAssociatorAtCaloFace.cc
0002 //
0003 // Original Author:  Andrea Rizzi
0004 //         Created:  Wed Apr 12 11:12:49 CEST 2006
0005 // Accommodated for Jet Package by: Fedor Ratnikov Jul. 30, 2007
0006 //
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 
0009 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0010 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "DataFormats/JetReco/interface/Jet.h"
0013 #include "DataFormats/JetReco/interface/CaloJet.h"
0014 #include "DataFormats/TrackReco/interface/Track.h"
0015 
0016 #include "JetTracksAssociatorAtCaloFace.h"
0017 
0018 JetTracksAssociatorAtCaloFace::JetTracksAssociatorAtCaloFace(const edm::ParameterSet& fConfig)
0019     : mGeometry(esConsumes()), dR_(fConfig.getParameter<double>("coneSize")) {
0020   mJets = consumes<edm::View<reco::Jet> >(fConfig.getParameter<edm::InputTag>("jets"));
0021   mExtrapolations =
0022       consumes<std::vector<reco::TrackExtrapolation> >(fConfig.getParameter<edm::InputTag>("extrapolations")),
0023 
0024   produces<reco::JetTracksAssociation::Container>();
0025 }
0026 
0027 void JetTracksAssociatorAtCaloFace::produce(edm::Event& fEvent, const edm::EventSetup& fSetup) {
0028   // Get geometry
0029   auto const& geo = fSetup.getData(mGeometry);
0030 
0031   // get stuff from Event
0032   edm::Handle<edm::View<reco::Jet> > jets_h;
0033   fEvent.getByToken(mJets, jets_h);
0034   edm::Handle<std::vector<reco::TrackExtrapolation> > extrapolations_h;
0035   fEvent.getByToken(mExtrapolations, extrapolations_h);
0036 
0037   auto jetTracks = std::make_unique<reco::JetTracksAssociation::Container>(reco::JetRefBaseProd(jets_h));
0038 
0039   // Check to make sure we have inputs
0040   if (jets_h->empty()) {
0041     // store output regardless the size of the inputs
0042     fEvent.put(std::move(jetTracks));
0043     return;
0044   }
0045   // Check to make sure the inputs are calo jets
0046   reco::CaloJet const* caloJet0 = dynamic_cast<reco::CaloJet const*>(&(jets_h->at(0)));
0047   // Disallowed non-CaloJet inputs
0048   if (caloJet0 == nullptr) {
0049     throw cms::Exception("InvalidInput") << " Jet-track association is only defined for CaloJets.";
0050   }
0051 
0052   // format inputs
0053   std::vector<edm::RefToBase<reco::Jet> > allJets;
0054   allJets.reserve(jets_h->size());
0055   for (unsigned i = 0; i < jets_h->size(); ++i)
0056     allJets.push_back(jets_h->refAt(i));
0057   mAssociator.produce(&*jetTracks, allJets, *extrapolations_h, geo, dR_);
0058 
0059   // store output
0060   fEvent.put(std::move(jetTracks));
0061 }