CMS 3D CMS Logo

L1SeedConePFJetProducer.cc
Go to the documentation of this file.
1 
2 #include <vector>
3 #include <numeric>
4 
6 // FRAMEWORK HEADERS
11 
15 
16 // bitwise emulation headers
19 
21 public:
23  ~L1SeedConePFJetProducer() override;
24 
25 private:
28  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
30 
31  float coneSize;
32  unsigned nJets;
33  bool HW;
34  bool debug;
37 
38  std::vector<l1t::PFJet> processEvent_SW(std::vector<edm::Ptr<l1t::PFCandidate>>& parts) const;
39  std::vector<l1t::PFJet> processEvent_HW(std::vector<edm::Ptr<l1t::PFCandidate>>& parts) const;
40 
42 
43  static std::pair<std::vector<L1SCJetEmu::Particle>,
44  std::unordered_map<const l1t::PFCandidate*, edm::Ptr<l1t::PFCandidate>>>
46 
47  static std::vector<l1t::PFJet> convertHWToEDM(
48  std::vector<L1SCJetEmu::Jet> hwJets,
49  std::unordered_map<const l1t::PFCandidate*, edm::Ptr<l1t::PFCandidate>> constituentMap);
50 };
51 
53  : coneSize(cfg.getParameter<double>("coneSize")),
54  nJets(cfg.getParameter<unsigned>("nJets")),
55  HW(cfg.getParameter<bool>("HW")),
56  debug(cfg.getParameter<bool>("debug")),
58  l1PFToken(consumes<std::vector<l1t::PFCandidate>>(cfg.getParameter<edm::InputTag>("L1PFObjects"))) {
59  produces<l1t::PFJetCollection>();
60 }
61 
64  const edm::EventSetup& iSetup) const {
65  std::unique_ptr<l1t::PFJetCollection> newPFJetCollection(new l1t::PFJetCollection);
66 
68  iEvent.getByToken(l1PFToken, l1PFCandidates);
69 
70  std::vector<edm::Ptr<l1t::PFCandidate>> particles;
71  for (unsigned i = 0; i < (*l1PFCandidates).size(); i++) {
72  particles.push_back(edm::Ptr<l1t::PFCandidate>(l1PFCandidates, i));
73  }
74 
75  std::vector<l1t::PFJet> jets;
76  if (HW) {
78  } else {
80  }
81 
82  std::sort(jets.begin(), jets.end(), [](l1t::PFJet i, l1t::PFJet j) { return (i.pt() > j.pt()); });
83  newPFJetCollection->swap(jets);
84  iEvent.put(std::move(newPFJetCollection));
85 }
86 
88 // DESTRUCTOR
90 
92  l1t::PFCandidate seed = *parts.at(0);
93 
94  auto sumpt = [](float a, const edm::Ptr<l1t::PFCandidate>& b) { return a + b->pt(); };
95 
96  // Sum the pt
97  float pt = std::accumulate(parts.begin(), parts.end(), 0., sumpt);
98 
99  // pt weighted d eta
100  std::vector<float> pt_deta;
101  pt_deta.resize(parts.size());
102  std::transform(parts.begin(), parts.end(), pt_deta.begin(), [&seed, &pt](const edm::Ptr<l1t::PFCandidate>& part) {
103  return (part->pt() / pt) * (part->eta() - seed.eta());
104  });
105  // Accumulate the pt weighted etas. Init to the seed eta, start accumulating at begin()+1 to skip seed
106  float eta = std::accumulate(pt_deta.begin() + 1, pt_deta.end(), seed.eta());
107 
108  // pt weighted d phi
109  std::vector<float> pt_dphi;
110  pt_dphi.resize(parts.size());
111  std::transform(parts.begin(), parts.end(), pt_dphi.begin(), [&seed, &pt](const edm::Ptr<l1t::PFCandidate>& part) {
112  return (part->pt() / pt) * reco::deltaPhi(part->phi(), seed.phi());
113  });
114  // Accumulate the pt weighted phis. Init to the seed phi, start accumulating at begin()+1 to skip seed
115  float phi = std::accumulate(pt_dphi.begin() + 1, pt_dphi.end(), seed.phi());
116 
117  l1t::PFJet jet(pt, eta, phi);
118  for (auto it = parts.begin(); it != parts.end(); it++) {
119  jet.addConstituent(*it);
120  }
121 
122  return jet;
123 }
124 
126  // The floating point algorithm simulation
127  std::stable_sort(work.begin(), work.end(), [](edm::Ptr<l1t::PFCandidate> i, edm::Ptr<l1t::PFCandidate> j) {
128  return (i->pt() > j->pt());
129  });
130  std::vector<l1t::PFJet> jets;
131  jets.reserve(nJets);
132  while (!work.empty() && jets.size() < nJets) {
133  // Take the first (highest pt) candidate as a seed
135  // Get the particles within a _coneSize of the seed
136  std::vector<edm::Ptr<l1t::PFCandidate>> particlesInCone;
137  std::copy_if(
138  work.begin(), work.end(), std::back_inserter(particlesInCone), [&](const edm::Ptr<l1t::PFCandidate>& part) {
139  return reco::deltaR<l1t::PFCandidate, l1t::PFCandidate>(*seed, *part) <= coneSize;
140  });
141  jets.push_back(makeJet_SW(particlesInCone));
142  // remove the clustered particles
143  work.erase(std::remove_if(work.begin(),
144  work.end(),
145  [&](const edm::Ptr<l1t::PFCandidate>& part) {
146  return reco::deltaR<l1t::PFCandidate, l1t::PFCandidate>(*seed, *part) <= coneSize;
147  }),
148  work.end());
149  }
150  return jets;
151 }
152 
154  // The fixed point emulator
155  // Convert the EDM format to the hardware format, and call the standalone emulator
156  std::pair<std::vector<L1SCJetEmu::Particle>, std::unordered_map<const l1t::PFCandidate*, edm::Ptr<l1t::PFCandidate>>>
158  std::vector<L1SCJetEmu::Jet> jets = emulator.emulateEvent(particles.first);
159  return convertHWToEDM(jets, particles.second);
160 }
161 
162 std::pair<std::vector<L1SCJetEmu::Particle>, std::unordered_map<const l1t::PFCandidate*, edm::Ptr<l1t::PFCandidate>>>
164  std::vector<l1ct::PuppiObjEmu> hwParticles;
165  std::unordered_map<const l1t::PFCandidate*, edm::Ptr<l1t::PFCandidate>> candidateMap;
166  std::for_each(edmParticles.begin(), edmParticles.end(), [&](edm::Ptr<l1t::PFCandidate>& edmParticle) {
167  l1ct::PuppiObjEmu particle;
168  particle.initFromBits(edmParticle->encodedPuppi64());
169  particle.srcCand = edmParticle.get();
170  candidateMap.insert(std::make_pair(edmParticle.get(), edmParticle));
171  hwParticles.push_back(particle);
172  });
173  return std::make_pair(hwParticles, candidateMap);
174 }
175 
177  std::vector<L1SCJetEmu::Jet> hwJets,
178  std::unordered_map<const l1t::PFCandidate*, edm::Ptr<l1t::PFCandidate>> constituentMap) {
179  std::vector<l1t::PFJet> edmJets;
180  std::for_each(hwJets.begin(), hwJets.end(), [&](L1SCJetEmu::Jet jet) {
181  l1t::PFJet edmJet(
182  jet.floatPt(), jet.floatEta(), jet.floatPhi(), /*mass=*/0., jet.intPt(), jet.intEta(), jet.intPhi());
183  // get back the references to the constituents
184  std::vector<edm::Ptr<l1t::PFCandidate>> constituents;
185  std::for_each(jet.constituents.begin(), jet.constituents.end(), [&](auto constituent) {
186  edmJet.addConstituent(constituentMap[constituent.srcCand]);
187  });
188  l1gt::Jet gtJet = jet.toGT();
189  edmJet.setEncodedJet(jet.toGT().pack());
190  edmJets.push_back(edmJet);
191  });
192  return edmJets;
193 }
194 
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
std::vector< l1t::PFJet > processEvent_SW(std::vector< edm::Ptr< l1t::PFCandidate >> &parts) const
l1t::PFJet makeJet_SW(const std::vector< edm::Ptr< l1t::PFCandidate >> &parts) const
delete x;
Definition: CaloConfig.h:22
T get() const
get a component
const l1t::PFCandidate * srcCand
std::vector< l1t::PFJet > PFJetCollection
Definition: PFJet.h:49
int iEvent
Definition: GenABIO.cc:224
static std::vector< l1t::PFJet > convertHWToEDM(std::vector< L1SCJetEmu::Jet > hwJets, std::unordered_map< const l1t::PFCandidate *, edm::Ptr< l1t::PFCandidate >> constituentMap)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
L1SeedConePFJetProducer(const edm::ParameterSet &)
std::vector< Jet > emulateEvent(std::vector< Particle > &parts) const
static std::pair< std::vector< L1SCJetEmu::Particle >, std::unordered_map< const l1t::PFCandidate *, edm::Ptr< l1t::PFCandidate > > > convertEDMToHW(std::vector< edm::Ptr< l1t::PFCandidate >> &edmParticles)
void initFromBits(const ap_uint< BITWIDTH > &src)
Definition: puppi.h:171
#define debug
Definition: HDRShower.cc:19
std::vector< l1t::PFJet > processEvent_HW(std::vector< edm::Ptr< l1t::PFCandidate >> &parts) const
part
Definition: HCALResponse.h:20
double b
Definition: hdecay.h:118
float coneSize
///////////////// ///
HLT enums.
double a
Definition: hdecay.h:119
edm::EDGetTokenT< std::vector< l1t::PFCandidate > > l1PFToken
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
def move(src, dest)
Definition: eostools.py:511
unsigned transform(const HcalDetId &id, unsigned transformCode)