CMS 3D CMS Logo

Run3ScoutingParticleToRecoPFCandidateProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
11 
15 
17 #include "fastjet/contrib/SoftKiller.hh"
18 
20 public:
23 
24  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
25  void beginStream(edm::StreamID) override {}
26  void produce(edm::Event &iEvent, edm::EventSetup const &setup) override;
27  void endStream() override {}
28 
29  void createPFCandidates(edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
30  std::unique_ptr<reco::PFCandidateCollection> &pfcands);
31  void createPFCandidatesSK(edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
32  std::unique_ptr<reco::PFCandidateCollection> &pfcands);
34  void clearVars();
35 
36 private:
40  bool use_CHS_;
42 
43  std::vector<int8_t> vertexIndex_;
44  std::vector<float> normchi2_;
45  std::vector<float> dz_;
46  std::vector<float> dxy_;
47  std::vector<float> dzsig_;
48  std::vector<float> dxysig_;
49  std::vector<int> lostInnerHits_;
50  std::vector<int> quality_;
51  std::vector<float> trkPt_;
52  std::vector<float> trkEta_;
53  std::vector<float> trkPhi_;
54 };
55 
56 //
57 // constructors and destructor
58 //
60  const edm::ParameterSet &iConfig)
61  : input_scoutingparticle_token_(consumes(iConfig.getParameter<edm::InputTag>("scoutingparticle"))),
62  particletable_token_(esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord>()),
63  use_softKiller_(iConfig.getParameter<bool>("softKiller")),
64  use_CHS_(iConfig.getParameter<bool>("CHS")) {
65  //register products
66  produces<reco::PFCandidateCollection>();
67  produces<edm::ValueMap<int>>("vertexIndex");
68  produces<edm::ValueMap<float>>("normchi2");
69  produces<edm::ValueMap<float>>("dz");
70  produces<edm::ValueMap<float>>("dxy");
71  produces<edm::ValueMap<float>>("dzsig");
72  produces<edm::ValueMap<float>>("dxysig");
73  produces<edm::ValueMap<int>>("lostInnerHits");
74  produces<edm::ValueMap<int>>("quality");
75  produces<edm::ValueMap<float>>("trkPt");
76  produces<edm::ValueMap<float>>("trkEta");
77  produces<edm::ValueMap<float>>("trkPhi");
78 }
79 
81 
83  auto m = pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId())) != nullptr
84  ? pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId()))->mass()
85  : -99.f;
86  auto q = pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId())) != nullptr
87  ? pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId()))->charge()
88  : -99.f;
89  if (m < -90 or q < -90) {
90  LogDebug("createPFCand") << "<Run3ScoutingParticleToRecoPFCandidateProducer::createPFCand>:" << std::endl
91  << "Unrecognisable pdgId - skipping particle" << std::endl;
92  return reco::PFCandidate();
93  }
94 
95  float px = scoutingparticle.pt() * cos(scoutingparticle.phi());
96  float py = scoutingparticle.pt() * sin(scoutingparticle.phi());
97  float pz = scoutingparticle.pt() * sinh(scoutingparticle.eta());
98  float p = scoutingparticle.pt() * cosh(scoutingparticle.eta());
99  float energy = std::sqrt(p * p + m * m);
101 
102  static const reco::PFCandidate dummy;
103  auto pfcand = reco::PFCandidate(q, p4, dummy.translatePdgIdToType(scoutingparticle.pdgId()));
104 
105  bool relativeTrackVars = scoutingparticle.relative_trk_vars();
106  vertexIndex_.push_back(scoutingparticle.vertex());
107  normchi2_.push_back(scoutingparticle.normchi2());
108  dz_.push_back(scoutingparticle.dz());
109  dxy_.push_back(scoutingparticle.dxy());
110  dzsig_.push_back(scoutingparticle.dzsig());
111  dxysig_.push_back(scoutingparticle.dxysig());
112  lostInnerHits_.push_back(scoutingparticle.lostInnerHits());
113  quality_.push_back(scoutingparticle.quality());
114  trkPt_.push_back(relativeTrackVars ? scoutingparticle.trk_pt() + scoutingparticle.pt() : scoutingparticle.trk_pt());
115  trkEta_.push_back(relativeTrackVars ? scoutingparticle.trk_eta() + scoutingparticle.eta()
116  : scoutingparticle.trk_eta());
117  trkPhi_.push_back(relativeTrackVars ? scoutingparticle.trk_phi() + scoutingparticle.phi()
118  : scoutingparticle.trk_phi());
119 
120  return pfcand;
121 }
122 
124  edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
125  std::unique_ptr<reco::PFCandidateCollection> &pfcands) {
126  for (unsigned int icand = 0; icand < scoutingparticleHandle->size(); ++icand) {
127  auto &scoutingparticle = (*scoutingparticleHandle)[icand];
128 
129  if (use_CHS_ and scoutingparticle.vertex() > 0)
130  continue;
131 
133  if (pfcand.energy() != 0)
134  pfcands->push_back(pfcand);
135  }
136 }
137 
139  edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
140  std::unique_ptr<reco::PFCandidateCollection> &pfcands) {
141  std::vector<fastjet::PseudoJet> fj;
142 
143  for (auto iter = scoutingparticleHandle->begin(),
144  ibegin = scoutingparticleHandle->begin(),
145  iend = scoutingparticleHandle->end();
146  iter != iend;
147  ++iter) {
148  auto m = pdTable_->particle(HepPDT::ParticleID(iter->pdgId())) != nullptr
149  ? pdTable_->particle(HepPDT::ParticleID(iter->pdgId()))->mass()
150  : -99.f;
151  if (m < -90) {
152  LogDebug("createPFCandidatesSK") << "<Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidatesSK>:"
153  << std::endl
154  << "Unrecognisable pdgId - skipping particle" << std::endl;
155  continue;
156  }
157  math::PtEtaPhiMLorentzVector p4(iter->pt(), iter->eta(), iter->phi(), m);
158  fj.push_back(fastjet::PseudoJet(p4.px(), p4.py(), p4.pz(), p4.energy()));
159  fj.back().set_user_index(iter - ibegin);
160  }
161 
162  fastjet::contrib::SoftKiller soft_killer(5, 0.4);
163  std::vector<fastjet::PseudoJet> soft_killed_particles = soft_killer(fj);
164 
165  for (auto &particle : soft_killed_particles) {
166  const Run3ScoutingParticle scoutingparticle = scoutingparticleHandle->at(particle.user_index());
168  if (pfcand.energy() != 0)
169  pfcands->push_back(pfcand);
170  }
171 }
172 
173 // ------------ method called to produce the data ------------
175  using namespace edm;
176 
177  auto pdt = setup.getHandle(particletable_token_);
178  pdTable_ = pdt.product();
179 
180  Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle;
181  iEvent.getByToken(input_scoutingparticle_token_, scoutingparticleHandle);
182 
183  auto pfcands = std::make_unique<reco::PFCandidateCollection>();
184 
185  if (use_softKiller_) {
186  createPFCandidatesSK(scoutingparticleHandle, pfcands);
187  } else {
188  createPFCandidates(scoutingparticleHandle, pfcands);
189  }
190 
192 
193  std::unique_ptr<edm::ValueMap<int>> vertexIndex_VM(new edm::ValueMap<int>());
194  edm::ValueMap<int>::Filler filler_vertexIndex(*vertexIndex_VM);
195  filler_vertexIndex.insert(oh, vertexIndex_.begin(), vertexIndex_.end());
196  filler_vertexIndex.fill();
197  iEvent.put(std::move(vertexIndex_VM), "vertexIndex");
198 
199  std::unique_ptr<edm::ValueMap<float>> normchi2_VM(new edm::ValueMap<float>());
200  edm::ValueMap<float>::Filler filler_normchi2(*normchi2_VM);
201  filler_normchi2.insert(oh, normchi2_.begin(), normchi2_.end());
202  filler_normchi2.fill();
203  iEvent.put(std::move(normchi2_VM), "normchi2");
204 
205  std::unique_ptr<edm::ValueMap<float>> dz_VM(new edm::ValueMap<float>());
206  edm::ValueMap<float>::Filler filler_dz(*dz_VM);
207  filler_dz.insert(oh, dz_.begin(), dz_.end());
208  filler_dz.fill();
209  iEvent.put(std::move(dz_VM), "dz");
210 
211  std::unique_ptr<edm::ValueMap<float>> dxy_VM(new edm::ValueMap<float>());
212  edm::ValueMap<float>::Filler filler_dxy(*dxy_VM);
213  filler_dxy.insert(oh, dxy_.begin(), dxy_.end());
214  filler_dxy.fill();
215  iEvent.put(std::move(dxy_VM), "dxy");
216 
217  std::unique_ptr<edm::ValueMap<float>> dzsig_VM(new edm::ValueMap<float>());
218  edm::ValueMap<float>::Filler filler_dzsig(*dzsig_VM);
219  filler_dzsig.insert(oh, dzsig_.begin(), dzsig_.end());
220  filler_dzsig.fill();
221  iEvent.put(std::move(dzsig_VM), "dzsig");
222 
223  std::unique_ptr<edm::ValueMap<float>> dxysig_VM(new edm::ValueMap<float>());
224  edm::ValueMap<float>::Filler filler_dxysig(*dxysig_VM);
225  filler_dxysig.insert(oh, dxysig_.begin(), dxysig_.end());
226  filler_dxysig.fill();
227  iEvent.put(std::move(dxysig_VM), "dxysig");
228 
229  std::unique_ptr<edm::ValueMap<int>> lostInnerHits_VM(new edm::ValueMap<int>());
230  edm::ValueMap<int>::Filler filler_lostInnerHits(*lostInnerHits_VM);
231  filler_lostInnerHits.insert(oh, lostInnerHits_.begin(), lostInnerHits_.end());
232  filler_lostInnerHits.fill();
233  iEvent.put(std::move(lostInnerHits_VM), "lostInnerHits");
234 
235  std::unique_ptr<edm::ValueMap<int>> quality_VM(new edm::ValueMap<int>());
236  edm::ValueMap<int>::Filler filler_quality(*quality_VM);
237  filler_quality.insert(oh, quality_.begin(), quality_.end());
238  filler_quality.fill();
239  iEvent.put(std::move(quality_VM), "quality");
240 
241  std::unique_ptr<edm::ValueMap<float>> trkPt_VM(new edm::ValueMap<float>());
242  edm::ValueMap<float>::Filler filler_trkPt(*trkPt_VM);
243  filler_trkPt.insert(oh, trkPt_.begin(), trkPt_.end());
244  filler_trkPt.fill();
245  iEvent.put(std::move(trkPt_VM), "trkPt");
246 
247  std::unique_ptr<edm::ValueMap<float>> trkEta_VM(new edm::ValueMap<float>());
248  edm::ValueMap<float>::Filler filler_trkEta(*trkEta_VM);
249  filler_trkEta.insert(oh, trkEta_.begin(), trkEta_.end());
250  filler_trkEta.fill();
251  iEvent.put(std::move(trkEta_VM), "trkEta");
252 
253  std::unique_ptr<edm::ValueMap<float>> trkPhi_VM(new edm::ValueMap<float>());
254  edm::ValueMap<float>::Filler filler_trkPhi(*trkPhi_VM);
255  filler_trkPhi.insert(oh, trkPhi_.begin(), trkPhi_.end());
256  filler_trkPhi.fill();
257  iEvent.put(std::move(trkPhi_VM), "trkPhi");
258 
259  clearVars();
260 }
261 
263  vertexIndex_.clear();
264  normchi2_.clear();
265  dz_.clear();
266  dxy_.clear();
267  dzsig_.clear();
268  dxysig_.clear();
269  lostInnerHits_.clear();
270  quality_.clear();
271  trkPt_.clear();
272  trkEta_.clear();
273  trkPhi_.clear();
274 }
275 
276 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
279  desc.add<edm::InputTag>("scoutingparticle", edm::InputTag("hltScoutingPFPacker"));
280  desc.add<bool>("softKiller", false);
281  desc.add<bool>("CHS", false);
282  descriptions.addWithDefaultLabel(desc);
283 }
284 
285 // declare this class as a framework plugin
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HepPDT::ParticleDataTable ParticleDataTable
void createPFCandidates(edm::Handle< std::vector< Run3ScoutingParticle >> scoutingparticleHandle, std::unique_ptr< reco::PFCandidateCollection > &pfcands)
void createPFCandidatesSK(edm::Handle< std::vector< Run3ScoutingParticle >> scoutingparticleHandle, std::unique_ptr< reco::PFCandidateCollection > &pfcands)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
reco::PFCandidate createPFCand(Run3ScoutingParticle scoutingparticle)
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particletable_token_
const edm::EDGetTokenT< std::vector< Run3ScoutingParticle > > input_scoutingparticle_token_
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
HLT enums.
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &iEvent, edm::EventSetup const &setup) override
#define LogDebug(id)