CMS 3D CMS Logo

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  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
30 
31 private:
35  const double ptMin_;
36 };
37 
39  : associationToken_(
40  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  produces<reco::TrackRefVector>();
45 }
46 
49  desc.add<edm::InputTag>("association", edm::InputTag("ak4JetTracksAssociatorAtVertexPF"));
50  desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
51  desc.add<edm::InputTag>("corrector", edm::InputTag("ak4PFL1FastL2L3Corrector"));
52  desc.add<double>("correctedPtMin", 0);
53  descriptions.add("jetTracksAssociationToTrackRefs", desc);
54 }
55 
58  iEvent.getByToken(associationToken_, h_assoc);
60 
62  iEvent.getByToken(jetToken_, h_jets);
63  const edm::View<reco::Jet>& jets = *h_jets;
64 
66  iEvent.getByToken(correctorToken_, h_corrector);
67  const reco::JetCorrector& corrector = *h_corrector;
68 
69  auto ret = std::make_unique<reco::TrackRefVector>();
70  std::unordered_set<reco::TrackRefVector::key_type> alreadyAdded;
71 
72  // Exctract tracks only for jets passing certain pT threshold after
73  // correction
74  for (size_t i = 0; i < jets.size(); ++i) {
75  edm::RefToBase<reco::Jet> jetRef = jets.refAt(i);
76  const reco::Jet& jet = *jetRef;
77 
78  auto p4 = jet.p4();
79 
80  // Energy correction in the most general way
81  if (!corrector.vectorialCorrection()) {
82  double scale = 1;
83  if (!corrector.refRequired()) {
84  scale = corrector.correction(jet);
85  } else {
86  scale = corrector.correction(jet, jetRef);
87  }
88  p4 = p4 * scale;
89  } else {
90  corrector.correction(jet, jetRef, p4);
91  }
92 
93  if (p4.pt() <= ptMin_)
94  continue;
95 
96  for (const auto& trackRef : association[jetRef]) {
97  if (alreadyAdded.find(trackRef.key()) == alreadyAdded.end()) {
98  ret->push_back(trackRef);
99  alreadyAdded.insert(trackRef.key());
100  }
101  }
102  }
103 
104  iEvent.put(std::move(ret));
105 }
106 
ConfigurationDescriptions.h
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:373
edm::StreamID
Definition: StreamID.h:30
JetTracksAssociation.h
mps_fire.i
i
Definition: mps_fire.py:428
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
L1EGammaCrystalsEmulatorProducer_cfi.scale
scale
Definition: L1EGammaCrystalsEmulatorProducer_cfi.py:10
sistrip::View
View
Definition: ConstantsForView.h:26
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::JetCorrector
Definition: JetCorrector.h:33
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
hgcal::association
std::tuple< layerClusterToCaloParticle, caloParticleToLayerCluster > association
Definition: LayerClusterAssociatorByEnergyScoreImpl.h:44
JetTracksAssociation
Association between jets and float value.
Jet.h
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle
Definition: AssociativeIterator.h:50
JetTracksAssociationToTrackRefs::ptMin_
const double ptMin_
Definition: JetTracksAssociationToTrackRefs.cc:35
JetCorrector
Definition: JetCorrector.h:19
reco::JetExtendedAssociation::Container
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
Definition: JetExtendedAssociation.h:29
JetTracksAssociationToTrackRefs::JetTracksAssociationToTrackRefs
JetTracksAssociationToTrackRefs(const edm::ParameterSet &iConfig)
Definition: JetTracksAssociationToTrackRefs.cc:38
pfClustersFromHGC3DClusters_cfi.corrector
corrector
Definition: pfClustersFromHGC3DClusters_cfi.py:5
JetTracksAssociationToTrackRefs::correctorToken_
edm::EDGetTokenT< reco::JetCorrector > correctorToken_
Definition: JetTracksAssociationToTrackRefs.cc:34
MakerMacros.h
Track.h
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Jet
Definition: Jet.py:1
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
JetTracksAssociationToTrackRefs::produce
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
Definition: JetTracksAssociationToTrackRefs.cc:56
JetTracksAssociationToTrackRefs::jetToken_
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition: JetTracksAssociationToTrackRefs.cc:33
JetTracksAssociationToTrackRefs
Definition: JetTracksAssociationToTrackRefs.cc:23
ParameterSetDescription.h
JetCorrector.h
edm::global::EDProducer
Definition: EDProducer.h:32
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
edm::AssociationVector
Definition: AssociationVector.h:67
edm::View
Definition: CaloClusterFwd.h:14
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
JetTracksAssociationToTrackRefs::associationToken_
edm::EDGetTokenT< reco::JetTracksAssociation::Container > associationToken_
Definition: JetTracksAssociationToTrackRefs.cc:32
iEvent
int iEvent
Definition: GenABIO.cc:224
p4
double p4[4]
Definition: TauolaWrapper.h:92
JetTracksAssociationToTrackRefs::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: JetTracksAssociationToTrackRefs.cc:47
edm::EventSetup
Definition: EventSetup.h:57
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
metsig::jet
Definition: SignAlgoResolutions.h:47
edm::RefToBase< reco::Jet >
ParameterSet.h
EDProducer.h
edm::Event
Definition: Event.h:73
EDProductfwd.h
edm::InputTag
Definition: InputTag.h:15