CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetTracksAssociationToTrackRefs.cc
Go to the documentation of this file.
7 
13 
15 
16 #include <unordered_set>
17 
24 public:
26 
27  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
28 
29  virtual void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
30 
31 private:
35  const double ptMin_;
36 };
37 
38 
40  associationToken_(consumes<reco::JetTracksAssociation::Container>(iConfig.getParameter<edm::InputTag>("association"))),
41  jetToken_(consumes<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("jets"))),
42  correctorToken_(consumes<reco::JetCorrector>(iConfig.getParameter<edm::InputTag>("corrector"))),
43  ptMin_(iConfig.getParameter<double>("correctedPtMin"))
44 {
45  produces<reco::TrackRefVector> ();
46 }
47 
50  desc.add<edm::InputTag>("association", edm::InputTag("ak4JetTracksAssociatorAtVertexPF"));
51  desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
52  desc.add<edm::InputTag>("corrector", edm::InputTag("ak4PFL1FastL2L3Corrector"));
53  desc.add<double>("correctedPtMin", 0);
54  descriptions.add("jetTracksAssociationToTrackRefs", desc);
55 }
56 
59  iEvent.getByToken(associationToken_, h_assoc);
60  const reco::JetTracksAssociation::Container& association = *h_assoc;
61 
63  iEvent.getByToken(jetToken_, h_jets);
64  const edm::View<reco::Jet>& jets = *h_jets;
65 
67  iEvent.getByToken(correctorToken_, h_corrector);
68  const reco::JetCorrector& corrector = *h_corrector;
69 
70  auto ret = std::make_unique<reco::TrackRefVector>();
71  std::unordered_set<reco::TrackRefVector::key_type> alreadyAdded;
72 
73  // Exctract tracks only for jets passing certain pT threshold after
74  // correction
75  for(size_t i=0; i<jets.size(); ++i) {
76  edm::RefToBase<reco::Jet> jetRef = jets.refAt(i);
77  const reco::Jet& jet = *jetRef;
78 
79  auto p4 = jet.p4();
80 
81  // Energy correction in the most general way
82  if(!corrector.vectorialCorrection()) {
83  double scale = 1;
84  if(!corrector.refRequired()) {
85  scale = corrector.correction(jet);
86  }
87  else {
88  scale = corrector.correction(jet, jetRef);
89  }
90  p4 = p4*scale;
91  }
92  else {
93  corrector.correction(jet, jetRef, p4);
94  }
95 
96  if(p4.pt() <= ptMin_)
97  continue;
98 
99  for(const auto& trackRef: association[jetRef]) {
100  if(alreadyAdded.find(trackRef.key()) == alreadyAdded.end()) {
101  ret->push_back(trackRef);
102  alreadyAdded.insert(trackRef.key());
103  }
104  }
105  }
106 
107  iEvent.put(std::move(ret));
108 }
109 
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
tuple ret
prodAgent to be discontinued
virtual void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Base class for all types of Jets.
Definition: Jet.h:20
bool refRequired() const
if correction needs the jet reference
Definition: JetCorrector.h:70
size_type size() const
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
JetTracksAssociationToTrackRefs(const edm::ParameterSet &iConfig)
edm::EDGetTokenT< reco::JetCorrector > correctorToken_
RefToBase< value_type > refAt(size_type i) const
int iEvent
Definition: GenABIO.cc:230
tuple corrector
Definition: mvaPFMET_cff.py:86
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
double p4[4]
Definition: TauolaWrapper.h:92
vector< PseudoJet > jets
def move
Definition: eostools.py:510
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::JetTracksAssociation::Container > associationToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Association between jets and float value.
bool vectorialCorrection() const
if vectorial correction is provided
Definition: JetCorrector.h:75
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99