CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
Run3ScoutingParticleToRecoPFCandidateProducer Class Reference
Inheritance diagram for Run3ScoutingParticleToRecoPFCandidateProducer:
edm::stream::EDProducer<>

Public Member Functions

void beginStream (edm::StreamID) override
 
void clearVars ()
 
reco::PFCandidate createPFCand (Run3ScoutingParticle scoutingparticle)
 
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)
 
void endStream () override
 
void produce (edm::Event &iEvent, edm::EventSetup const &setup) override
 
 Run3ScoutingParticleToRecoPFCandidateProducer (const edm::ParameterSet &)
 
 ~Run3ScoutingParticleToRecoPFCandidateProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

std::vector< float > dxy_
 
std::vector< float > dxysig_
 
std::vector< float > dz_
 
std::vector< float > dzsig_
 
const edm::EDGetTokenT< std::vector< Run3ScoutingParticle > > input_scoutingparticle_token_
 
std::vector< int > lostInnerHits_
 
std::vector< float > normchi2_
 
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecordparticletable_token_
 
const HepPDT::ParticleDataTablepdTable_
 
std::vector< int > quality_
 
std::vector< float > trkEta_
 
std::vector< float > trkPhi_
 
std::vector< float > trkPt_
 
bool use_CHS_
 
bool use_softKiller_
 
std::vector< int8_t > vertexIndex_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 19 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

Constructor & Destructor Documentation

◆ Run3ScoutingParticleToRecoPFCandidateProducer()

Run3ScoutingParticleToRecoPFCandidateProducer::Run3ScoutingParticleToRecoPFCandidateProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 59 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

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 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particletable_token_
const edm::EDGetTokenT< std::vector< Run3ScoutingParticle > > input_scoutingparticle_token_

◆ ~Run3ScoutingParticleToRecoPFCandidateProducer()

Run3ScoutingParticleToRecoPFCandidateProducer::~Run3ScoutingParticleToRecoPFCandidateProducer ( )
overridedefault

Member Function Documentation

◆ beginStream()

void Run3ScoutingParticleToRecoPFCandidateProducer::beginStream ( edm::StreamID  )
inlineoverride

Definition at line 25 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

25 {}

◆ clearVars()

void Run3ScoutingParticleToRecoPFCandidateProducer::clearVars ( )

Definition at line 262 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

References dxy_, dxysig_, dz_, dzsig_, lostInnerHits_, normchi2_, quality_, trkEta_, trkPhi_, trkPt_, and vertexIndex_.

Referenced by produce().

262  {
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 }

◆ createPFCand()

reco::PFCandidate Run3ScoutingParticleToRecoPFCandidateProducer::createPFCand ( Run3ScoutingParticle  scoutingparticle)

Definition at line 82 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

References funct::cos(), dxy_, dxysig_, dz_, dzsig_, hcalRecHitTable_cff::energy, LogDebug, lostInnerHits_, visualization-live-secondInstance_cfg::m, normchi2_, or, AlCaHLTBitMon_ParallelJobs::p, LHEGenericFilter_cfi::ParticleID, pdTable_, pfDeepBoostedJetPreprocessParams_cfi::pfcand, pfLinker_cff::PFCandidate, multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, submitPVResolutionJobs::q, quality_, HLT_2024v14_cff::relativeTrackVars, run3scouting_cff::scoutingparticle, funct::sin(), mathSSE::sqrt(), trkEta_, trkPhi_, trkPt_, and vertexIndex_.

Referenced by createPFCandidates(), and createPFCandidatesSK().

82  {
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 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:23
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
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
#define LogDebug(id)

◆ createPFCandidates()

void Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidates ( edm::Handle< std::vector< Run3ScoutingParticle >>  scoutingparticleHandle,
std::unique_ptr< reco::PFCandidateCollection > &  pfcands 
)

Definition at line 123 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

References createPFCand(), pfDeepBoostedJetPreprocessParams_cfi::pfcand, HLT_2024v14_cff::pfcands, run3scouting_cff::scoutingparticle, and use_CHS_.

Referenced by produce().

125  {
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 }
reco::PFCandidate createPFCand(Run3ScoutingParticle scoutingparticle)

◆ createPFCandidatesSK()

void Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidatesSK ( edm::Handle< std::vector< Run3ScoutingParticle >>  scoutingparticleHandle,
std::unique_ptr< reco::PFCandidateCollection > &  pfcands 
)

Definition at line 138 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

References createPFCand(), LogDebug, visualization-live-secondInstance_cfg::m, LHEGenericFilter_cfi::ParticleID, pdTable_, pfDeepBoostedJetPreprocessParams_cfi::pfcand, HLT_2024v14_cff::pfcands, and run3scouting_cff::scoutingparticle.

Referenced by produce().

140  {
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 }
reco::PFCandidate createPFCand(Run3ScoutingParticle scoutingparticle)
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
#define LogDebug(id)

◆ endStream()

void Run3ScoutingParticleToRecoPFCandidateProducer::endStream ( )
inlineoverride

Definition at line 27 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

27 {}

◆ fillDescriptions()

void Run3ScoutingParticleToRecoPFCandidateProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 277 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

277  {
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 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

◆ produce()

void Run3ScoutingParticleToRecoPFCandidateProducer::produce ( edm::Event iEvent,
edm::EventSetup const &  setup 
)
override

Definition at line 174 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

References clearVars(), createPFCandidates(), createPFCandidatesSK(), dxy_, dxysig_, dz_, dzsig_, edm::helper::Filler< Map >::fill(), iEvent, input_scoutingparticle_token_, edm::helper::Filler< Map >::insert(), lostInnerHits_, eostools::move(), normchi2_, particletable_token_, pdTable_, HLT_2024v14_cff::pfcands, quality_, singleTopDQM_cfi::setup, trkEta_, trkPhi_, trkPt_, use_softKiller_, and vertexIndex_.

174  {
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 }
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)
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particletable_token_
const edm::EDGetTokenT< std::vector< Run3ScoutingParticle > > input_scoutingparticle_token_
HLT enums.
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ dxy_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::dxy_
private

◆ dxysig_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::dxysig_
private

◆ dz_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::dz_
private

◆ dzsig_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::dzsig_
private

◆ input_scoutingparticle_token_

const edm::EDGetTokenT<std::vector<Run3ScoutingParticle> > Run3ScoutingParticleToRecoPFCandidateProducer::input_scoutingparticle_token_
private

Definition at line 37 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

Referenced by produce().

◆ lostInnerHits_

std::vector<int> Run3ScoutingParticleToRecoPFCandidateProducer::lostInnerHits_
private

◆ normchi2_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::normchi2_
private

◆ particletable_token_

const edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> Run3ScoutingParticleToRecoPFCandidateProducer::particletable_token_
private

Definition at line 38 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

Referenced by produce().

◆ pdTable_

const HepPDT::ParticleDataTable* Run3ScoutingParticleToRecoPFCandidateProducer::pdTable_
private

◆ quality_

std::vector<int> Run3ScoutingParticleToRecoPFCandidateProducer::quality_
private

◆ trkEta_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::trkEta_
private

◆ trkPhi_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::trkPhi_
private

◆ trkPt_

std::vector<float> Run3ScoutingParticleToRecoPFCandidateProducer::trkPt_
private

◆ use_CHS_

bool Run3ScoutingParticleToRecoPFCandidateProducer::use_CHS_
private

Definition at line 40 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

Referenced by createPFCandidates().

◆ use_softKiller_

bool Run3ScoutingParticleToRecoPFCandidateProducer::use_softKiller_
private

Definition at line 39 of file Run3ScoutingParticleToRecoPFCandidateProducer.cc.

Referenced by produce().

◆ vertexIndex_

std::vector<int8_t> Run3ScoutingParticleToRecoPFCandidateProducer::vertexIndex_
private