CMS 3D CMS Logo

JetSubstructurePacker.cc
Go to the documentation of this file.
4 
6  : distMax_(iConfig.getParameter<double>("distMax")),
7  jetToken_(consumes<edm::View<pat::Jet>>(iConfig.getParameter<edm::InputTag>("jetSrc"))),
8  algoLabels_(iConfig.getParameter<std::vector<std::string>>("algoLabels")),
9  algoTags_(iConfig.getParameter<std::vector<edm::InputTag>>("algoTags")),
10  fixDaughters_(iConfig.getParameter<bool>("fixDaughters")) {
11  algoTokens_ =
12  edm::vector_transform(algoTags_, [this](edm::InputTag const& tag) { return consumes<edm::View<pat::Jet>>(tag); });
13  if (fixDaughters_) {
14  pf2pc_ = consumes<edm::Association<pat::PackedCandidateCollection>>(
15  iConfig.getParameter<edm::InputTag>("packedPFCandidates"));
16  pc2pf_ = consumes<edm::Association<reco::PFCandidateCollection>>(
17  iConfig.getParameter<edm::InputTag>("packedPFCandidates"));
18  }
19  //register products
20  produces<std::vector<pat::Jet>>();
21 }
22 
24 
25 // ------------ method called to produce the data ------------
27  auto outputs = std::make_unique<std::vector<pat::Jet>>();
28 
30  std::vector<edm::Handle<edm::View<pat::Jet>>> algoHandles;
31 
34  if (fixDaughters_) {
35  iEvent.getByToken(pf2pc_, pf2pc);
36  iEvent.getByToken(pc2pf_, pc2pf);
37  }
38 
39  iEvent.getByToken(jetToken_, jetHandle);
40  algoHandles.resize(algoTags_.size());
41  for (size_t i = 0; i < algoTags_.size(); ++i) {
42  iEvent.getByToken(algoTokens_[i], algoHandles[i]);
43  }
44 
45  // Loop over the input jets that will be modified.
46  for (auto const& ijet : *jetHandle) {
47  // Copy the jet.
48  outputs->push_back(ijet);
49  // Loop over the substructure collections
50  unsigned int index = 0;
51 
52  for (auto const& ialgoHandle : algoHandles) {
53  std::vector<edm::Ptr<pat::Jet>> nextSubjets;
54  float dRMin = distMax_;
55  for (auto const& jjet : *ialgoHandle) {
56  if (reco::deltaR(ijet, jjet) < dRMin) {
57  for (auto const& userfloatstr : jjet.userFloatNames()) {
58  outputs->back().addUserFloat(userfloatstr, jjet.userFloat(userfloatstr));
59  }
60  for (auto const& userintstr : jjet.userIntNames()) {
61  outputs->back().addUserInt(userintstr, jjet.userInt(userintstr));
62  }
63  for (auto const& usercandstr : jjet.userCandNames()) {
64  outputs->back().addUserCand(usercandstr, jjet.userCand(usercandstr));
65  }
66  for (size_t ida = 0; ida < jjet.numberOfDaughters(); ++ida) {
67  reco::CandidatePtr candPtr = jjet.daughterPtr(ida);
68  nextSubjets.push_back(edm::Ptr<pat::Jet>(candPtr));
69  }
70  break;
71  }
72  }
73  outputs->back().addSubjets(nextSubjets, algoLabels_[index]);
74  ++index;
75  }
76 
77  // fix daughters
78  if (fixDaughters_) {
79  std::vector<reco::CandidatePtr> daughtersInSubjets;
80  std::vector<reco::CandidatePtr> daughtersNew;
81  const std::vector<reco::CandidatePtr>& jdausPF = outputs->back().daughterPtrVector();
82  std::vector<reco::CandidatePtr> jdaus;
83  jdaus.reserve(jdausPF.size());
84  // Convert the daughters to packed candidates. This is easier than ref-navigating through PUPPI or CHS to particleFlow.
85  for (auto const& jdau : jdausPF) {
86  jdaus.push_back(edm::refToPtr((*pf2pc)[jdau]));
87  }
88 
89  for (const edm::Ptr<pat::Jet>& subjet : outputs->back().subjets()) {
90  const std::vector<reco::CandidatePtr>& sjdaus = subjet->daughterPtrVector();
91  // check that the subjet does not contain any extra constituents not contained in the jet
92  bool skipSubjet = false;
93  for (const reco::CandidatePtr& dau : sjdaus) {
94  if (std::find(jdaus.begin(), jdaus.end(), dau) == jdaus.end()) {
95  skipSubjet = true;
96  break;
97  }
98  }
99  if (skipSubjet)
100  continue;
101 
102  daughtersInSubjets.insert(daughtersInSubjets.end(), sjdaus.begin(), sjdaus.end());
103  daughtersNew.push_back(reco::CandidatePtr(subjet));
104  }
105  for (const reco::CandidatePtr& dau : jdaus) {
106  if (std::find(daughtersInSubjets.begin(), daughtersInSubjets.end(), dau) == daughtersInSubjets.end()) {
107  daughtersNew.push_back(dau);
108  }
109  }
110  outputs->back().clearDaughters();
111  for (const auto& dau : daughtersNew)
112  outputs->back().addDaughter(dau);
113  }
114  }
115 
116  iEvent.put(std::move(outputs));
117 }
118 
119 //define this as a plug-in
electrons_cff.bool
bool
Definition: electrons_cff.py:393
mps_fire.i
i
Definition: mps_fire.py:428
sistrip::View
View
Definition: ConstantsForView.h:26
edm
HLT enums.
Definition: AlignableModifier.h:19
PatBasicFWLiteJetAnalyzer_Selector_cfg.outputs
outputs
Definition: PatBasicFWLiteJetAnalyzer_Selector_cfg.py:48
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:85964
JetSubstructurePacker
Definition: JetSubstructurePacker.h:34
HLT_FULL_cff.dRMin
dRMin
Definition: HLT_FULL_cff.py:8727
JetSubstructurePacker::~JetSubstructurePacker
~JetSubstructurePacker() override
Definition: JetSubstructurePacker.cc:23
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
edm::Handle
Definition: AssociativeIterator.h:50
JetSubstructurePacker::algoTags_
std::vector< edm::InputTag > algoTags_
Definition: JetSubstructurePacker.h:48
JetSubstructurePacker::jetToken_
edm::EDGetTokenT< edm::View< pat::Jet > > jetToken_
Definition: JetSubstructurePacker.h:46
edm::refToPtr
Ptr< typename C::value_type > refToPtr(Ref< C, typename C::value_type, refhelper::FindUsingAdvance< C, typename C::value_type > > const &ref)
Definition: RefToPtr.h:18
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Jet
Definition: Jet.py:1
GlobalPosition_Frontier_DevDB_cff.tag
tag
Definition: GlobalPosition_Frontier_DevDB_cff.py:11
JetSubstructurePacker::pc2pf_
edm::EDGetTokenT< edm::Association< reco::PFCandidateCollection > > pc2pf_
Definition: JetSubstructurePacker.h:52
JetSubstructurePacker::fixDaughters_
bool fixDaughters_
Definition: JetSubstructurePacker.h:50
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::vector_transform
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
JetSubstructurePacker::algoTokens_
std::vector< edm::EDGetTokenT< edm::View< pat::Jet > > > algoTokens_
Definition: JetSubstructurePacker.h:49
edm::ParameterSet
Definition: ParameterSet.h:47
deltaR.h
JetSubstructurePacker::JetSubstructurePacker
JetSubstructurePacker(const edm::ParameterSet &)
Definition: JetSubstructurePacker.cc:5
JetSubstructurePacker::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: JetSubstructurePacker.cc:26
JetSubstructurePacker::distMax_
float distMax_
Definition: JetSubstructurePacker.h:45
iEvent
int iEvent
Definition: GenABIO.cc:224
RefToPtr.h
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
edm::EventSetup
Definition: EventSetup.h:57
pat
Definition: HeavyIon.h:7
edm::Ptr< Candidate >
JetSubstructurePacker::pf2pc_
edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > pf2pc_
Definition: JetSubstructurePacker.h:51
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
JetSubstructurePacker::algoLabels_
std::vector< std::string > algoLabels_
Definition: JetSubstructurePacker.h:47
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
JetSubstructurePacker.h