CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes
L2TauNNProducerAlpaka Class Reference
Inheritance diagram for L2TauNNProducerAlpaka:
edm::stream::EDProducer< edm::GlobalCache< L2TauNNProducerAlpakaCacheData > >

Classes

struct  caloRecHitCollections
 
struct  InputDescTau
 

Public Types

using TracksHost = pixelTrack::TracksHostPhase1
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< L2TauNNProducerAlpakaCacheData > >
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
 

Public Member Functions

 L2TauNNProducerAlpaka (const edm::ParameterSet &, const L2TauNNProducerAlpakaCacheData *)
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< L2TauNNProducerAlpakaCacheData > >
 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 &)
 
static void globalEndJob (L2TauNNProducerAlpakaCacheData *)
 
static std::unique_ptr< L2TauNNProducerAlpakaCacheDatainitializeGlobalCache (const edm::ParameterSet &)
 

Static Public Attributes

static constexpr float dEta_width = 2 * L2TauTagNNv1::dR_max / static_cast<float>(L2TauTagNNv1::nCellEta)
 
static constexpr float dPhi_width = 2 * L2TauTagNNv1::dR_max / static_cast<float>(L2TauTagNNv1::nCellPhi)
 
static constexpr float dR2_max = L2TauTagNNv1::dR_max * L2TauTagNNv1::dR_max
 

Private Member Functions

void checknan (tensorflow::Tensor &tensor, int debugLevel)
 
void fillCaloRecHits (tensorflow::Tensor &cellGridMatrix, const std::vector< l1t::TauRef > &allTaus, const caloRecHitCollections &caloRecHits)
 
void fillL1TauVars (tensorflow::Tensor &cellGridMatrix, const std::vector< l1t::TauRef > &allTaus)
 
void fillPatatracks (tensorflow::Tensor &cellGridMatrix, const std::vector< l1t::TauRef > &allTaus, const TracksHost &patatracks_tsoa, const ZVertexHost &patavtx_soa, const reco::BeamSpot &beamspot, const MagneticField *magfi)
 
template<typename VPos , typename LVec >
std::tuple< float, float, int, int > getEtaPhiIndices (const VPos &position, const LVec &tau_p4)
 
template<typename LVec >
std::tuple< float, float, int, int > getEtaPhiIndices (float eta, float phi, const LVec &tau_p4)
 
std::vector< float > getTauScore (const tensorflow::Tensor &cellGridMatrix)
 
std::pair< float, float > impactParameter (int it, const TracksHost &patatracks_tsoa, float patatrackPhi, const reco::BeamSpot &beamspot, const MagneticField *magfi)
 
void produce (edm::Event &event, const edm::EventSetup &eventsetup) override
 
void selectGoodTracksAndVertices (const ZVertexHost &patavtx_soa, const TracksHost &patatracks_tsoa, std::vector< int > &trkGood, std::vector< int > &vtxGood)
 
void standardizeTensor (tensorflow::Tensor &tensor)
 

Private Attributes

const edm::EDGetTokenT< reco::BeamSpotbeamSpotToken_
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordbFieldToken_
 
const int debugLevel_
 
const edm::EDGetTokenT< EcalRecHitCollectionebToken_
 
const edm::EDGetTokenT< EcalRecHitCollectioneeToken_
 
const float fractionSumPt2_
 
const edm::ESGetToken< CaloGeometry, CaloGeometryRecordgeometryToken_
 
const edm::EDGetTokenT< HBHERecHitCollectionhbheToken_
 
const edm::EDGetTokenT< HORecHitCollectionhoToken_
 
std::string inputTensorName_
 
std::vector< InputDescTauL1TauDesc_
 
const L2TauNNProducerAlpakaCacheDataL2cacheData_
 
const unsigned int maxVtx_
 
const float minSumPt2_
 
std::string outputTensorName_
 
const edm::EDGetTokenT< TracksHostpataTracksToken_
 
const edm::EDGetTokenT< ZVertexHostpataVerticesToken_
 
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefstauTriggerToken_
 
const float trackChi2Max_
 
const float trackPtMax_
 
const float trackPtMin_
 

Detailed Description

Definition at line 147 of file L2TauTagNNProducerAlpaka.cc.

Member Typedef Documentation

◆ TracksHost

Definition at line 149 of file L2TauTagNNProducerAlpaka.cc.

Constructor & Destructor Documentation

◆ L2TauNNProducerAlpaka()

L2TauNNProducerAlpaka::L2TauNNProducerAlpaka ( const edm::ParameterSet cfg,
const L2TauNNProducerAlpakaCacheData cacheData 
)
explicit

Definition at line 290 of file L2TauTagNNProducerAlpaka.cc.

References looper::cfg, L2TauNNProducerAlpaka::InputDescTau::CollectionName, submitPVResolutionJobs::desc, Exception, L2TauNNProducerAlpakaCacheData::graphDef, inputTensorName_, L2TauNNProducerAlpaka::InputDescTau::inputToken_, L1TauDesc_, L2cacheData_, Skims_PA_cff::name, outputTensorName_, and AlCaHLTBitMon_QueryRunRegistry::string.

292  : debugLevel_(cfg.getParameter<int>("debugLevel")),
293  hbheToken_(consumes<HBHERecHitCollection>(cfg.getParameter<edm::InputTag>("hbheInput"))),
294  hoToken_(consumes<HORecHitCollection>(cfg.getParameter<edm::InputTag>("hoInput"))),
295  ebToken_(consumes<EcalRecHitCollection>(cfg.getParameter<edm::InputTag>("ebInput"))),
296  eeToken_(consumes<EcalRecHitCollection>(cfg.getParameter<edm::InputTag>("eeInput"))),
297  geometryToken_(esConsumes<CaloGeometry, CaloGeometryRecord>()),
298  bFieldToken_(esConsumes<MagneticField, IdealMagneticFieldRecord>()),
299  pataVerticesToken_(consumes(cfg.getParameter<edm::InputTag>("pataVertices"))),
300  pataTracksToken_(consumes(cfg.getParameter<edm::InputTag>("pataTracks"))),
301  beamSpotToken_(consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("BeamSpot"))),
302  maxVtx_(cfg.getParameter<uint>("maxVtx")),
303  fractionSumPt2_(cfg.getParameter<double>("fractionSumPt2")),
304  minSumPt2_(cfg.getParameter<double>("minSumPt2")),
305  trackPtMin_(cfg.getParameter<double>("track_pt_min")),
306  trackPtMax_(cfg.getParameter<double>("track_pt_max")),
307  trackChi2Max_(cfg.getParameter<double>("track_chi2_max")) {
308  if (cacheData->graphDef == nullptr) {
309  throw cms::Exception("InvalidCacheData") << "Invalid Cache Data.";
310  }
311  inputTensorName_ = cacheData->graphDef->node(0).name();
312  outputTensorName_ = cacheData->graphDef->node(cacheData->graphDef->node_size() - 1).name();
313  L2cacheData_ = cacheData;
314  std::vector<edm::ParameterSet> L1TauCollections = cfg.getParameter<std::vector<edm::ParameterSet>>("L1Taus");
315  L1TauDesc_.reserve(L1TauCollections.size());
316  for (const auto& l1TauInput : L1TauCollections) {
317  InputDescTau toInsert;
318  toInsert.CollectionName = l1TauInput.getParameter<std::string>("L1CollectionName");
319  toInsert.inputToken_ =
320  consumes<trigger::TriggerFilterObjectWithRefs>(l1TauInput.getParameter<edm::InputTag>("L1TauTrigger"));
321  L1TauDesc_.push_back(toInsert);
322  }
323  for (const auto& desc : L1TauDesc_)
324  produces<std::vector<float>>(desc.CollectionName);
325 }
const L2TauNNProducerAlpakaCacheData * L2cacheData_
const edm::EDGetTokenT< TracksHost > pataTracksToken_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > bFieldToken_
const edm::EDGetTokenT< HBHERecHitCollection > hbheToken_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geometryToken_
const edm::EDGetTokenT< ZVertexHost > pataVerticesToken_
const edm::EDGetTokenT< HORecHitCollection > hoToken_
std::vector< InputDescTau > L1TauDesc_
const edm::EDGetTokenT< EcalRecHitCollection > eeToken_
const edm::EDGetTokenT< EcalRecHitCollection > ebToken_

Member Function Documentation

◆ checknan()

void L2TauNNProducerAlpaka::checknan ( tensorflow::Tensor &  tensor,
int  debugLevel 
)
private

Definition at line 327 of file L2TauTagNNProducerAlpaka.cc.

References ztail::d, HLT_2024v10_cff::debugLevel, Exception, input, edm::isNotFinite(), and L2TauTagNNv1::varNameMap.

Referenced by produce().

327  {
329  std::vector<int> tensor_shape(tensor.shape().dims());
330  for (int d = 0; d < tensor.shape().dims(); d++) {
331  tensor_shape.at(d) = tensor.shape().dim_size(d);
332  }
333  if (tensor_shape.size() != 4) {
334  throw cms::Exception("InvalidTensor") << "Tensor shape does not have 4 dimensions!";
335  }
336  for (int tau_idx = 0; tau_idx < tensor_shape.at(0); tau_idx++) {
337  for (int phi_idx = 0; phi_idx < tensor_shape.at(1); phi_idx++) {
338  for (int eta_idx = 0; eta_idx < tensor_shape.at(2); eta_idx++) {
339  for (int var_idx = 0; var_idx < tensor_shape.at(3); var_idx++) {
340  auto getCell = [&](NNInputs input) -> float& {
341  return getCellImpl(tensor, tau_idx, phi_idx, eta_idx, input);
342  };
343  auto nonstd_var = getCell(static_cast<NNInputs>(var_idx));
344  if (edm::isNotFinite(nonstd_var)) {
345  edm::LogWarning("InputVar") << "var is nan \nvar name= "
346  << L2TauTagNNv1::varNameMap.at(static_cast<L2TauTagNNv1::NNInputs>(var_idx))
347  << "\t var_idx = " << var_idx << "\t eta_idx = " << eta_idx
348  << "\t phi_idx = " << phi_idx << "\t tau_idx = " << tau_idx;
349  if (debugLevel > 2) {
350  edm::LogWarning("InputVar") << "other vars in same cell \n";
351  if (var_idx + 1 < tensor_shape.at(3))
352  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 1))
353  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 1));
354  if (var_idx + 2 < tensor_shape.at(3))
355  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 2))
356  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 2));
357  if (var_idx + 3 < tensor_shape.at(3))
358  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 3))
359  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 3));
360  if (var_idx + 4 < tensor_shape.at(3))
361  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 4))
362  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 4));
363  }
364  }
365  }
366  }
367  }
368  }
369 }
const std::map< NNInputs, std::string > varNameMap
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
static std::string const input
Definition: EdmProvDump.cc:50
d
Definition: ztail.py:151
Log< level::Warning, false > LogWarning

◆ fillCaloRecHits()

void L2TauNNProducerAlpaka::fillCaloRecHits ( tensorflow::Tensor &  cellGridMatrix,
const std::vector< l1t::TauRef > &  allTaus,
const caloRecHitCollections caloRecHits 
)
private

Definition at line 437 of file L2TauTagNNProducerAlpaka.cc.

References reco::deltaR2(), dR2_max, L2TauNNProducerAlpaka::caloRecHitCollections::eb, L2TauNNProducerAlpaka::caloRecHitCollections::ee, L2TauNNProducerAlpaka::caloRecHitCollections::geometry, getEtaPhiIndices(), CaloGeometry::getGeometry(), L2TauNNProducerAlpaka::caloRecHitCollections::hbhe, L2TauNNProducerAlpaka::caloRecHitCollections::ho, input, L2TauTagNNv1::nCellEta, L2TauTagNNv1::nCellPhi, and position.

Referenced by produce().

439  {
441 
442  const int nTaus = allTaus.size();
443  float deta, dphi;
444  int eta_idx = 0;
445  int phi_idx = 0;
446  int tau_idx = 0;
447 
448  auto getCell = [&](NNInputs input) -> float& {
449  return getCellImpl(cellGridMatrix, tau_idx, phi_idx, eta_idx, input);
450  };
451  for (tau_idx = 0; tau_idx < nTaus; tau_idx++) {
452  // calorechit_EE
453  for (const auto& caloRecHit_ee : *caloRecHits.ee) {
454  if (caloRecHit_ee.energy() <= 0)
455  continue;
456  const auto& position = caloRecHits.geometry->getGeometry(caloRecHit_ee.id())->getPosition();
457  const float eeCalEn = caloRecHit_ee.energy();
458  const float eeCalChi2 = caloRecHit_ee.chi2();
459  if (reco::deltaR2(position, allTaus[tau_idx]->polarP4()) < dR2_max) {
460  std::tie(deta, dphi, eta_idx, phi_idx) = getEtaPhiIndices(position, allTaus[tau_idx]->polarP4());
461  getCell(NNInputs::EcalEnergySum) += eeCalEn;
462  getCell(NNInputs::EcalSize) += 1.;
463  getCell(NNInputs::EcalEnergyStdDev) += eeCalEn * eeCalEn;
464  getCell(NNInputs::EcalDeltaEta) += deta * eeCalEn;
465  getCell(NNInputs::EcalDeltaPhi) += dphi * eeCalEn;
466  if (eeCalChi2 >= 0) {
467  getCell(NNInputs::EcalChi2) += eeCalChi2 * eeCalEn;
468  getCell(NNInputs::EcalEnergySumForPositiveChi2) += eeCalEn;
469  getCell(NNInputs::EcalSizeForPositiveChi2) += 1.;
470  }
471  }
472  }
473 
474  // calorechit_EB
475  for (const auto& caloRecHit_eb : *caloRecHits.eb) {
476  if (caloRecHit_eb.energy() <= 0)
477  continue;
478  const auto& position = caloRecHits.geometry->getGeometry(caloRecHit_eb.id())->getPosition();
479  const float ebCalEn = caloRecHit_eb.energy();
480  const float ebCalChi2 = caloRecHit_eb.chi2();
481  if (reco::deltaR2(position, allTaus[tau_idx]->polarP4()) < dR2_max) {
482  std::tie(deta, dphi, eta_idx, phi_idx) = getEtaPhiIndices(position, allTaus[tau_idx]->polarP4());
483  getCell(NNInputs::EcalEnergySum) += ebCalEn;
484  getCell(NNInputs::EcalSize) += 1.;
485  getCell(NNInputs::EcalEnergyStdDev) += ebCalEn * ebCalEn;
486  getCell(NNInputs::EcalDeltaEta) += deta * ebCalEn;
487  getCell(NNInputs::EcalDeltaPhi) += dphi * ebCalEn;
488  if (ebCalChi2 >= 0) {
489  getCell(NNInputs::EcalChi2) += ebCalChi2 * ebCalEn;
490  getCell(NNInputs::EcalEnergySumForPositiveChi2) += ebCalEn;
491  getCell(NNInputs::EcalSizeForPositiveChi2) += 1.;
492  }
493  }
494  }
495 
496  // calorechit_HBHE
497  for (const auto& caloRecHit_hbhe : *caloRecHits.hbhe) {
498  if (caloRecHit_hbhe.energy() <= 0)
499  continue;
500  const auto& position = caloRecHits.geometry->getGeometry(caloRecHit_hbhe.id())->getPosition();
501  const float hbheCalEn = caloRecHit_hbhe.energy();
502  const float hbheCalChi2 = caloRecHit_hbhe.chi2();
503  if (reco::deltaR2(position, allTaus[tau_idx]->polarP4()) < dR2_max) {
504  std::tie(deta, dphi, eta_idx, phi_idx) = getEtaPhiIndices(position, allTaus[tau_idx]->polarP4());
505  getCell(NNInputs::HcalEnergySum) += hbheCalEn;
506  getCell(NNInputs::HcalEnergyStdDev) += hbheCalEn * hbheCalEn;
507  getCell(NNInputs::HcalSize) += 1.;
508  getCell(NNInputs::HcalDeltaEta) += deta * hbheCalEn;
509  getCell(NNInputs::HcalDeltaPhi) += dphi * hbheCalEn;
510  if (hbheCalChi2 >= 0) {
511  getCell(NNInputs::HcalChi2) += hbheCalChi2 * hbheCalEn;
512  getCell(NNInputs::HcalEnergySumForPositiveChi2) += hbheCalEn;
513  getCell(NNInputs::HcalSizeForPositiveChi2) += 1.;
514  }
515  }
516  }
517 
518  // calorechit_HO
519  for (const auto& caloRecHit_ho : *caloRecHits.ho) {
520  if (caloRecHit_ho.energy() <= 0)
521  continue;
522  const auto& position = caloRecHits.geometry->getGeometry(caloRecHit_ho.id())->getPosition();
523  const float hoCalEn = caloRecHit_ho.energy();
524  if (reco::deltaR2(position, allTaus[tau_idx]->polarP4()) < dR2_max) {
525  std::tie(deta, dphi, eta_idx, phi_idx) = getEtaPhiIndices(position, allTaus[tau_idx]->polarP4());
526  getCell(NNInputs::HcalEnergySum) += hoCalEn;
527  getCell(NNInputs::HcalEnergyStdDev) += hoCalEn * hoCalEn;
528  getCell(NNInputs::HcalSize) += 1.;
529  getCell(NNInputs::HcalDeltaEta) += deta * hoCalEn;
530  getCell(NNInputs::HcalDeltaPhi) += dphi * hoCalEn;
531  }
532  }
533 
534  // normalize to sum and define stdDev
535  for (eta_idx = 0; eta_idx < L2TauTagNNv1::nCellEta; eta_idx++) {
536  for (phi_idx = 0; phi_idx < L2TauTagNNv1::nCellPhi; phi_idx++) {
537  /* normalize eCal vars*/
538  if (getCell(NNInputs::EcalEnergySum) > 0.) {
539  getCell(NNInputs::EcalDeltaEta) /= getCell(NNInputs::EcalEnergySum);
540  getCell(NNInputs::EcalDeltaPhi) /= getCell(NNInputs::EcalEnergySum);
541  }
542  if (getCell(NNInputs::EcalEnergySumForPositiveChi2) > 0.) {
543  getCell(NNInputs::EcalChi2) /= getCell(NNInputs::EcalEnergySumForPositiveChi2);
544  }
545  if (getCell(NNInputs::EcalSize) > 1.) {
546  // (stdDev - (enSum*enSum)/size) / (size-1)
547  getCell(NNInputs::EcalEnergyStdDev) =
548  (getCell(NNInputs::EcalEnergyStdDev) -
549  (getCell(NNInputs::EcalEnergySum) * getCell(NNInputs::EcalEnergySum)) / getCell(NNInputs::EcalSize)) /
550  (getCell(NNInputs::EcalSize) - 1);
551  } else {
552  getCell(NNInputs::EcalEnergyStdDev) = 0.;
553  }
554  /* normalize hCal Vars */
555  if (getCell(NNInputs::HcalEnergySum) > 0.) {
556  getCell(NNInputs::HcalDeltaEta) /= getCell(NNInputs::HcalEnergySum);
557  getCell(NNInputs::HcalDeltaPhi) /= getCell(NNInputs::HcalEnergySum);
558  }
559  if (getCell(NNInputs::HcalEnergySumForPositiveChi2) > 0.) {
560  getCell(NNInputs::HcalChi2) /= getCell(NNInputs::HcalEnergySumForPositiveChi2);
561  }
562  if (getCell(NNInputs::HcalSize) > 1.) {
563  // (stdDev - (enSum*enSum)/size) / (size-1)
564  getCell(NNInputs::HcalEnergyStdDev) =
565  (getCell(NNInputs::HcalEnergyStdDev) -
566  (getCell(NNInputs::HcalEnergySum) * getCell(NNInputs::HcalEnergySum)) / getCell(NNInputs::HcalSize)) /
567  (getCell(NNInputs::HcalSize) - 1);
568  } else {
569  getCell(NNInputs::HcalEnergyStdDev) = 0.;
570  }
571  }
572  }
573  }
574 }
static std::string const input
Definition: EdmProvDump.cc:50
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static constexpr float dR2_max
std::tuple< float, float, int, int > getEtaPhiIndices(const VPos &position, const LVec &tau_p4)
constexpr int nCellEta
static int position[264][3]
Definition: ReadPGInfo.cc:289
constexpr int nCellPhi

◆ fillDescriptions()

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

Definition at line 258 of file L2TauTagNNProducerAlpaka.cc.

References edm::ParameterSetDescription::add(), edm::ParameterSet::addParameter(), edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, AlCaHLTBitMon_QueryRunRegistry::string, and parallelization::uint.

258  {
260  desc.add<int>("debugLevel", 0)->setComment("set debug level for printing out info");
261  edm::ParameterSetDescription l1TausPset;
262  l1TausPset.add<std::string>("L1CollectionName", "DoubleTau")->setComment("Name of collections");
263  l1TausPset.add<edm::InputTag>("L1TauTrigger", edm::InputTag("hltL1sDoubleTauBigOR"))
264  ->setComment("Which trigger should the L1 Taus collection pass");
265  edm::ParameterSet l1TausPSetDefault;
266  l1TausPSetDefault.addParameter<std::string>("L1CollectionName", "DoubleTau");
267  l1TausPSetDefault.addParameter<edm::InputTag>("L1TauTrigger", edm::InputTag("hltL1sDoubleTauBigOR"));
268  desc.addVPSet("L1Taus", l1TausPset, {l1TausPSetDefault});
269  desc.add<edm::InputTag>("hbheInput", edm::InputTag("hltHbhereco"))->setComment("HBHE recHit collection");
270  desc.add<edm::InputTag>("hoInput", edm::InputTag("hltHoreco"))->setComment("HO recHit Collection");
271  desc.add<edm::InputTag>("ebInput", edm::InputTag("hltEcalRecHit:EcalRecHitsEB"))->setComment("EB recHit Collection");
272  desc.add<edm::InputTag>("eeInput", edm::InputTag("hltEcalRecHit:EcalRecHitsEE"))->setComment("EE recHit Collection");
273  desc.add<edm::InputTag>("pataVertices", edm::InputTag("hltPixelVerticesSoA"))
274  ->setComment("patatrack vertices collection");
275  desc.add<edm::InputTag>("pataTracks", edm::InputTag("hltPixelTracksSoA"))->setComment("patatrack collection");
276  desc.add<edm::InputTag>("BeamSpot", edm::InputTag("hltOnlineBeamSpot"))->setComment("BeamSpot Collection");
277  desc.add<uint>("maxVtx", 100)->setComment("max output collection size (number of accepted vertices)");
278  desc.add<double>("fractionSumPt2", 0.3)->setComment("threshold on sumPt2 fraction of the leading vertex");
279  desc.add<double>("minSumPt2", 0.)->setComment("min sumPt2");
280  desc.add<double>("track_pt_min", 1.0)->setComment("min track p_T");
281  desc.add<double>("track_pt_max", 10.0)->setComment("max track p_T");
282  desc.add<double>("track_chi2_max", 99999.)->setComment("max track chi2");
283  desc.add<std::string>("graphPath", "RecoTauTag/TrainingFiles/data/L2TauNNTag/L2TauTag_Run3v1.pb")
284  ->setComment("path to the saved CNN");
285  desc.add<std::string>("normalizationDict", "RecoTauTag/TrainingFiles/data/L2TauNNTag/NormalizationDict.json")
286  ->setComment("path to the dictionary for variable standardization");
287  descriptions.addWithDefaultLabel(desc);
288 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:136
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ fillL1TauVars()

void L2TauNNProducerAlpaka::fillL1TauVars ( tensorflow::Tensor &  cellGridMatrix,
const std::vector< l1t::TauRef > &  allTaus 
)
private

Definition at line 405 of file L2TauTagNNProducerAlpaka.cc.

References input, L2TauTagNNv1::nCellEta, and L2TauTagNNv1::nCellPhi.

Referenced by produce().

405  {
407 
408  const int nTaus = allTaus.size();
409  for (int tau_idx = 0; tau_idx < nTaus; tau_idx++) {
410  for (int eta_idx = 0; eta_idx < L2TauTagNNv1::nCellEta; eta_idx++) {
411  for (int phi_idx = 0; phi_idx < L2TauTagNNv1::nCellPhi; phi_idx++) {
412  auto getCell = [&](NNInputs input) -> float& {
413  return getCellImpl(cellGridMatrix, tau_idx, phi_idx, eta_idx, input);
414  };
415  getCell(NNInputs::l1Tau_pt) = allTaus[tau_idx]->pt();
416  getCell(NNInputs::l1Tau_eta) = allTaus[tau_idx]->eta();
417  getCell(NNInputs::l1Tau_hwIso) = allTaus[tau_idx]->hwIso();
418  }
419  }
420  }
421 }
static std::string const input
Definition: EdmProvDump.cc:50
constexpr int nCellEta
constexpr int nCellPhi

◆ fillPatatracks()

void L2TauNNProducerAlpaka::fillPatatracks ( tensorflow::Tensor &  cellGridMatrix,
const std::vector< l1t::TauRef > &  allTaus,
const TracksHost patatracks_tsoa,
const ZVertexHost patavtx_soa,
const reco::BeamSpot beamspot,
const MagneticField magfi 
)
private

Definition at line 655 of file L2TauTagNNProducerAlpaka.cc.

References reco::charge(), PortableHostCollection< T >::const_view(), reco::deltaR2(), dR2_max, spr::find(), getEtaPhiIndices(), impactParameter(), input, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, SiStripPI::min, L2TauTagNNv1::nCellEta, L2TauTagNNv1::nCellPhi, nHits, TrackCollections2monitor_cff::nVertices, reco::phi(), selectGoodTracksAndVertices(), and PortableHostCollection< T >::view().

Referenced by produce().

660  {
662  using patatrackHelpers = TracksUtilities<pixelTopology::Phase1>;
663  float deta, dphi;
664  int eta_idx = 0;
665  int phi_idx = 0;
666  int tau_idx = 0;
667 
668  auto getCell = [&](NNInputs input) -> float& {
669  return getCellImpl(cellGridMatrix, tau_idx, phi_idx, eta_idx, input);
670  };
671 
672  std::vector<int> trkGood;
673  std::vector<int> vtxGood;
674 
675  selectGoodTracksAndVertices(patavtx_soa, patatracks_tsoa, trkGood, vtxGood);
676 
677  const int nTaus = allTaus.size();
678  for (tau_idx = 0; tau_idx < nTaus; tau_idx++) {
679  const float tauEta = allTaus[tau_idx]->eta();
680  const float tauPhi = allTaus[tau_idx]->phi();
681 
682  for (const auto it : trkGood) {
683  const float patatrackPt = patatracks_tsoa.const_view()[it].pt();
684  if (patatrackPt <= 0)
685  continue;
686  const float patatrackPhi = reco::phi(patatracks_tsoa.const_view(), it);
687  const float patatrackEta = patatracks_tsoa.const_view()[it].eta();
688  const float patatrackCharge = reco::charge(patatracks_tsoa.const_view(), it);
689  const float patatrackChi2OverNdof = patatracks_tsoa.view()[it].chi2();
690  const auto nHits = patatrackHelpers::nHits(patatracks_tsoa.const_view(), it);
691  if (nHits <= 0)
692  continue;
693  const int patatrackNdof = 2 * std::min(6, nHits) - 5;
694 
695  const int vtx_idx_assTrk = patavtx_soa.view()[it].idv();
696  if (reco::deltaR2(patatrackEta, patatrackPhi, tauEta, tauPhi) < dR2_max) {
697  std::tie(deta, dphi, eta_idx, phi_idx) =
698  getEtaPhiIndices(patatrackEta, patatrackPhi, allTaus[tau_idx]->polarP4());
699  getCell(NNInputs::PatatrackPtSum) += patatrackPt;
700  getCell(NNInputs::PatatrackSize) += 1.;
701  getCell(NNInputs::PatatrackChargeSum) += patatrackCharge;
702  getCell(NNInputs::PatatrackDeltaEta) += deta * patatrackPt;
703  getCell(NNInputs::PatatrackDeltaPhi) += dphi * patatrackPt;
704  getCell(NNInputs::PatatrackChi2OverNdof) += patatrackChi2OverNdof * patatrackPt;
705  getCell(NNInputs::PatatrackNdof) += patatrackNdof * patatrackPt;
706  std::pair<float, float> impactParameters = impactParameter(it, patatracks_tsoa, patatrackPhi, beamspot, magfi);
707  getCell(NNInputs::PatatrackDxy) += impactParameters.first * patatrackPt;
708  getCell(NNInputs::PatatrackDz) += impactParameters.second * patatrackPt;
709  if ((std::find(vtxGood.begin(), vtxGood.end(), vtx_idx_assTrk) != vtxGood.end())) {
710  getCell(NNInputs::PatatrackPtSumWithVertex) += patatrackPt;
711  getCell(NNInputs::PatatrackSizeWithVertex) += 1.;
712  }
713  }
714  }
715 
716  // normalize to sum and define stdDev
717  for (eta_idx = 0; eta_idx < L2TauTagNNv1::nCellEta; eta_idx++) {
718  for (phi_idx = 0; phi_idx < L2TauTagNNv1::nCellPhi; phi_idx++) {
719  getCell(NNInputs::nVertices) = vtxGood.size();
720  if (getCell(NNInputs::PatatrackPtSum) > 0.) {
721  getCell(NNInputs::PatatrackDeltaEta) /= getCell(NNInputs::PatatrackPtSum);
722  getCell(NNInputs::PatatrackDeltaPhi) /= getCell(NNInputs::PatatrackPtSum);
723  getCell(NNInputs::PatatrackChi2OverNdof) /= getCell(NNInputs::PatatrackPtSum);
724  getCell(NNInputs::PatatrackNdof) /= getCell(NNInputs::PatatrackPtSum);
725  getCell(NNInputs::PatatrackDxy) /= getCell(NNInputs::PatatrackPtSum);
726  getCell(NNInputs::PatatrackDz) /= getCell(NNInputs::PatatrackPtSum);
727  }
728  }
729  }
730  }
731 }
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float charge(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:73
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
static std::string const input
Definition: EdmProvDump.cc:50
ConstView const & const_view() const
void selectGoodTracksAndVertices(const ZVertexHost &patavtx_soa, const TracksHost &patatracks_tsoa, std::vector< int > &trkGood, std::vector< int > &vtxGood)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float phi(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:80
std::pair< float, float > impactParameter(int it, const TracksHost &patatracks_tsoa, float patatrackPhi, const reco::BeamSpot &beamspot, const MagneticField *magfi)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static constexpr float dR2_max
std::tuple< float, float, int, int > getEtaPhiIndices(const VPos &position, const LVec &tau_p4)
constexpr int nCellEta
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
constexpr int nCellPhi

◆ getEtaPhiIndices() [1/2]

template<typename VPos , typename LVec >
std::tuple< float, float, int, int > L2TauNNProducerAlpaka::getEtaPhiIndices ( const VPos &  position,
const LVec &  tau_p4 
)
private

Definition at line 433 of file L2TauTagNNProducerAlpaka.cc.

References position.

Referenced by fillCaloRecHits(), and fillPatatracks().

433  {
434  return getEtaPhiIndices(position.eta(), position.phi(), tau_p4);
435 }
std::tuple< float, float, int, int > getEtaPhiIndices(const VPos &position, const LVec &tau_p4)
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ getEtaPhiIndices() [2/2]

template<typename LVec >
std::tuple< float, float, int, int > L2TauNNProducerAlpaka::getEtaPhiIndices ( float  eta,
float  phi,
const LVec &  tau_p4 
)
private

Definition at line 424 of file L2TauTagNNProducerAlpaka.cc.

References reco::deltaPhi(), dEta_width, dPhi_width, L2TauTagNNv1::dR_max, PVValHelper::eta, and phi.

424  {
425  const float deta = eta - tau_p4.eta();
426  const float dphi = reco::deltaPhi(phi, tau_p4.phi());
427  const int eta_idx = static_cast<int>(floor((deta + L2TauTagNNv1::dR_max) / dEta_width));
428  const int phi_idx = static_cast<int>(floor((dphi + L2TauTagNNv1::dR_max) / dPhi_width));
429  return std::make_tuple(deta, dphi, eta_idx, phi_idx);
430 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
constexpr float dR_max
static constexpr float dEta_width
static constexpr float dPhi_width

◆ getTauScore()

std::vector< float > L2TauNNProducerAlpaka::getTauScore ( const tensorflow::Tensor &  cellGridMatrix)
private

Definition at line 733 of file L2TauTagNNProducerAlpaka.cc.

References L2cacheData_, tensorflow::run(), and L2TauNNProducerAlpakaCacheData::session.

Referenced by produce().

733  {
734  std::vector<tensorflow::Tensor> pred_tensor;
735  tensorflow::run(L2cacheData_->session, {{inputTensorName_, cellGridMatrix}}, {outputTensorName_}, &pred_tensor);
736  const int nTau = cellGridMatrix.shape().dim_size(0);
737  std::vector<float> pred_vector(nTau);
738  for (int tau_idx = 0; tau_idx < nTau; ++tau_idx) {
739  pred_vector[tau_idx] = pred_tensor[0].matrix<float>()(tau_idx, 0);
740  }
741 
742  return pred_vector;
743 }
const L2TauNNProducerAlpakaCacheData * L2cacheData_
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:259

◆ globalEndJob()

void L2TauNNProducerAlpaka::globalEndJob ( L2TauNNProducerAlpakaCacheData cacheData)
static

Definition at line 252 of file L2TauTagNNProducerAlpaka.cc.

References tensorflow::closeSession(), L2TauNNProducerAlpakaCacheData::graphDef, and L2TauNNProducerAlpakaCacheData::session.

252  {
253  if (cacheData->graphDef != nullptr) {
254  delete cacheData->graphDef;
255  }
256  tensorflow::closeSession(cacheData->session);
257 }
bool closeSession(Session *&session)
Definition: TensorFlow.cc:234

◆ impactParameter()

std::pair< float, float > L2TauNNProducerAlpaka::impactParameter ( int  it,
const TracksHost patatracks_tsoa,
float  patatrackPhi,
const reco::BeamSpot beamspot,
const MagneticField magfi 
)
private

Definition at line 624 of file L2TauTagNNProducerAlpaka.cc.

References LocalTrajectoryParameters::charge(), TracksUtilities< TrackerTraits >::copyToDense(), funct::cos(), f, runTauDisplay::gp, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, M_PI_2, LocalTrajectoryParameters::momentum(), phi, LocalTrajectoryParameters::position(), createTree::pp, funct::sin(), riemannFit::transformToPerigeePlane(), and PortableHostCollection< T >::view().

Referenced by fillPatatracks().

628  {
629  /* dxy and dz */
630  riemannFit::Vector5d ipar, opar;
631  riemannFit::Matrix5d icov, ocov;
632  TracksUtilities<pixelTopology::Phase1>::copyToDense(patatracks_tsoa.view(), ipar, icov, it);
633  riemannFit::transformToPerigeePlane(ipar, icov, opar, ocov);
634  LocalTrajectoryParameters lpar(opar(0), opar(1), opar(2), opar(3), opar(4), 1.);
635  float sp = std::sin(patatrackPhi);
636  float cp = std::cos(patatrackPhi);
637  Surface::RotationType Rotation(sp, -cp, 0, 0, 0, -1.f, cp, sp, 0);
638  GlobalPoint BeamSpotPoint(beamspot.x0(), beamspot.y0(), beamspot.z0());
639  Plane impPointPlane(BeamSpotPoint, Rotation);
641  impPointPlane.toGlobal(lpar.position()), impPointPlane.toGlobal(lpar.momentum()), lpar.charge(), magfi);
642  GlobalPoint vv = gp.position();
643  math::XYZPoint pos(vv.x(), vv.y(), vv.z());
644  GlobalVector pp = gp.momentum();
645  math::XYZVector mom(pp.x(), pp.y(), pp.z());
646  auto lambda = M_PI_2 - pp.theta();
647  auto phi = pp.phi();
648  float patatrackDxy = -vv.x() * std::sin(phi) + vv.y() * std::cos(phi);
649  float patatrackDz =
650  (vv.z() * std::cos(lambda) - (vv.x() * std::cos(phi) + vv.y() * std::sin(phi)) * std::sin(lambda)) /
651  std::cos(lambda);
652  return std::make_pair(patatrackDxy, patatrackDz);
653 }
static constexpr __host__ __device__ void copyToDense(const TrackSoAConstView &tracks, V5 &v, M5 &cov, int32_t i)
Eigen::Matrix< double, 5, 5 > Matrix5d
Definition: FitResult.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Eigen::Matrix< double, 5, 1 > Vector5d
Definition: FitResult.h:13
#define M_PI_2
void transformToPerigeePlane(VI5 const &ip, MI5 const &icov, VO5 &op, MO5 &ocov)
Definition: FitUtils.h:60
Definition: Plane.h:16
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double f[11][100]
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
SOARotation< float > Rotation

◆ initializeGlobalCache()

std::unique_ptr< L2TauNNProducerAlpakaCacheData > L2TauNNProducerAlpaka::initializeGlobalCache ( const edm::ParameterSet cfg)
static

Definition at line 226 of file L2TauTagNNProducerAlpaka.cc.

References looper::cfg, tensorflow::createSession(), edm::FileInPath::fullPath(), HLT_2024v10_cff::graphPath, submitPVResolutionJobs::key, tensorflow::loadGraphDef(), normDictElement::max, normDictElement::mean, normDictElement::min, HLT_2024v10_cff::normalizationDict, L2TauTagNNv1::nVars, tensorflow::setLogging(), normDictElement::std, AlCaHLTBitMon_QueryRunRegistry::string, heppy_batch::val, trigObjTnPSource_cfi::var, and L2TauTagNNv1::varNameMap.

227  {
228  std::unique_ptr<L2TauNNProducerAlpakaCacheData> cacheData = std::make_unique<L2TauNNProducerAlpakaCacheData>();
229  cacheData->normVec.reserve(L2TauTagNNv1::nVars);
230 
231  auto const graphPath = edm::FileInPath(cfg.getParameter<std::string>("graphPath")).fullPath();
232 
233  cacheData->graphDef = tensorflow::loadGraphDef(graphPath);
234  cacheData->session = tensorflow::createSession(cacheData->graphDef);
235 
237 
238  boost::property_tree::ptree loadPtreeRoot;
239  auto const normalizationDict = edm::FileInPath(cfg.getParameter<std::string>("normalizationDict")).fullPath();
240  boost::property_tree::read_json(normalizationDict, loadPtreeRoot);
241  for (const auto& [key, val] : L2TauTagNNv1::varNameMap) {
242  boost::property_tree::ptree var = loadPtreeRoot.get_child(val);
243  normDictElement current_element;
244  current_element.mean = var.get_child("mean").get_value<float>();
245  current_element.std = var.get_child("std").get_value<float>();
246  current_element.min = var.get_child("min").get_value<float>();
247  current_element.max = var.get_child("max").get_value<float>();
248  cacheData->normVec.push_back(current_element);
249  }
250  return cacheData;
251 }
const std::map< NNInputs, std::string > varNameMap
std::string fullPath() const
Definition: FileInPath.cc:161
GraphDef * loadGraphDef(const std::string &pbFile)
Definition: TensorFlow.cc:120
constexpr int nVars
key
prepare the HTCondor submission files and eventually submit them
Session * createSession()
Definition: TensorFlow.cc:137
void setLogging(const std::string &level="3")
Definition: TensorFlow.cc:81

◆ produce()

void L2TauNNProducerAlpaka::produce ( edm::Event event,
const edm::EventSetup eventsetup 
)
overrideprivate

Definition at line 745 of file L2TauTagNNProducerAlpaka.cc.

References beamSpotToken_, bFieldToken_, checknan(), debugLevel_, HLT_2024v10_cff::distance, L2TauNNProducerAlpaka::caloRecHitCollections::eb, ebToken_, L2TauNNProducerAlpaka::caloRecHitCollections::ee, eeToken_, fillCaloRecHits(), fillL1TauVars(), fillPatatracks(), spr::find(), relativeConstraints::geometry, L2TauNNProducerAlpaka::caloRecHitCollections::geometry, geometryToken_, edm::EventSetup::getHandle(), getTauScore(), photonIsolationHIProducer_cfi::hbhe, L2TauNNProducerAlpaka::caloRecHitCollections::hbhe, hbheToken_, photonIsolationHIProducer_cfi::ho, L2TauNNProducerAlpaka::caloRecHitCollections::ho, hoToken_, L1TauDesc_, eostools::move(), L2TauTagNNv1::nCellEta, L2TauTagNNv1::nCellPhi, L2TauTagNNv1::nVars, pataTracksToken_, pataVerticesToken_, DiDispStaMuonMonitor_cfi::pt, standardizeTensor(), and trigger::TriggerL1Tau.

745  {
746  std::vector<std::vector<size_t>> TauCollectionMap(L1TauDesc_.size());
747  l1t::TauVectorRef allTaus;
748 
749  for (size_t inp_idx = 0; inp_idx < L1TauDesc_.size(); inp_idx++) {
750  l1t::TauVectorRef l1Taus;
751  auto const& l1TriggeredTaus = event.get(L1TauDesc_[inp_idx].inputToken_);
752  l1TriggeredTaus.getObjects(trigger::TriggerL1Tau, l1Taus);
753  TauCollectionMap.at(inp_idx).resize(l1Taus.size());
754 
755  for (size_t l1_idx = 0; l1_idx < l1Taus.size(); l1_idx++) {
756  size_t tau_idx;
757  const auto iter = std::find(allTaus.begin(), allTaus.end(), l1Taus[l1_idx]);
758  if (iter != allTaus.end()) {
759  tau_idx = std::distance(allTaus.begin(), iter);
760  } else {
761  allTaus.push_back(l1Taus[l1_idx]);
762  tau_idx = allTaus.size() - 1;
763  }
764  TauCollectionMap.at(inp_idx).at(l1_idx) = tau_idx;
765  }
766  }
767  const auto ebCal = event.getHandle(ebToken_);
768  const auto eeCal = event.getHandle(eeToken_);
769  const auto hbhe = event.getHandle(hbheToken_);
770  const auto ho = event.getHandle(hoToken_);
771  auto const& patatracks_SoA = event.get(pataTracksToken_);
772  auto const& vertices_SoA = event.get(pataVerticesToken_);
773  const auto bsHandle = event.getHandle(beamSpotToken_);
774 
775  auto const fieldESH = eventsetup.getHandle(bFieldToken_);
776  auto const geometry = eventsetup.getHandle(geometryToken_);
777 
778  caloRecHitCollections caloRecHits;
779  caloRecHits.hbhe = &*hbhe;
780  caloRecHits.ho = &*ho;
781  caloRecHits.eb = &*ebCal;
782  caloRecHits.ee = &*eeCal;
783  caloRecHits.geometry = &*geometry;
784 
785  const int nTaus = allTaus.size();
786  tensorflow::Tensor cellGridMatrix(tensorflow::DT_FLOAT,
788  const int n_inputs = nTaus * L2TauTagNNv1::nCellEta * L2TauTagNNv1::nCellPhi * L2TauTagNNv1::nVars;
789  for (int input_idx = 0; input_idx < n_inputs; ++input_idx) {
790  cellGridMatrix.flat<float>()(input_idx) = 0;
791  }
792  fillL1TauVars(cellGridMatrix, allTaus);
793 
794  fillCaloRecHits(cellGridMatrix, allTaus, caloRecHits);
795 
796  fillPatatracks(cellGridMatrix, allTaus, patatracks_SoA, vertices_SoA, *bsHandle, fieldESH.product());
797 
798  standardizeTensor(cellGridMatrix);
799 
800  if (debugLevel_ > 0) {
801  checknan(cellGridMatrix, debugLevel_);
802  }
803 
804  std::vector<float> tau_score = getTauScore(cellGridMatrix);
805 
806  for (size_t inp_idx = 0; inp_idx < L1TauDesc_.size(); inp_idx++) {
807  const size_t nTau = TauCollectionMap[inp_idx].size();
808  auto tau_tags = std::make_unique<std::vector<float>>(nTau);
809  for (size_t tau_pos = 0; tau_pos < nTau; ++tau_pos) {
810  const auto tau_idx = TauCollectionMap[inp_idx][tau_pos];
811  if (debugLevel_ > 0) {
812  edm::LogInfo("DebugInfo") << event.id().event() << " \t " << (allTaus[tau_idx])->pt() << " \t "
813  << tau_score.at(tau_idx) << std::endl;
814  }
815  (*tau_tags)[tau_pos] = tau_score.at(tau_idx);
816  }
817  event.put(std::move(tau_tags), L1TauDesc_[inp_idx].CollectionName);
818  }
819 }
const edm::EDGetTokenT< TracksHost > pataTracksToken_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > bFieldToken_
const edm::EDGetTokenT< HBHERecHitCollection > hbheToken_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
void fillPatatracks(tensorflow::Tensor &cellGridMatrix, const std::vector< l1t::TauRef > &allTaus, const TracksHost &patatracks_tsoa, const ZVertexHost &patavtx_soa, const reco::BeamSpot &beamspot, const MagneticField *magfi)
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geometryToken_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void fillCaloRecHits(tensorflow::Tensor &cellGridMatrix, const std::vector< l1t::TauRef > &allTaus, const caloRecHitCollections &caloRecHits)
void checknan(tensorflow::Tensor &tensor, int debugLevel)
const edm::EDGetTokenT< ZVertexHost > pataVerticesToken_
constexpr int nVars
const edm::EDGetTokenT< HORecHitCollection > hoToken_
std::vector< TauRef > TauVectorRef
Definition: Tau.h:14
std::vector< InputDescTau > L1TauDesc_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
Log< level::Info, false > LogInfo
const edm::EDGetTokenT< EcalRecHitCollection > eeToken_
const edm::EDGetTokenT< EcalRecHitCollection > ebToken_
constexpr int nCellEta
void fillL1TauVars(tensorflow::Tensor &cellGridMatrix, const std::vector< l1t::TauRef > &allTaus)
std::vector< float > getTauScore(const tensorflow::Tensor &cellGridMatrix)
def move(src, dest)
Definition: eostools.py:511
void standardizeTensor(tensorflow::Tensor &tensor)
constexpr int nCellPhi

◆ selectGoodTracksAndVertices()

void L2TauNNProducerAlpaka::selectGoodTracksAndVertices ( const ZVertexHost patavtx_soa,
const TracksHost patatracks_tsoa,
std::vector< int > &  trkGood,
std::vector< int > &  vtxGood 
)
private

Definition at line 576 of file L2TauTagNNProducerAlpaka.cc.

References cms::cuda::assert(), PortableHostCollection< T >::const_view(), fractionSumPt2_, dqmiolumiharvest::j, pixelTrack::loose, DMR_cfg::maxTracks, maxVtx_, SiStripPI::min, minSumPt2_, nHits, quality, trackChi2Max_, trackPtMax_, trackPtMin_, and PortableHostCollection< T >::view().

Referenced by fillPatatracks().

579  {
580  using patatrackHelpers = TracksUtilities<pixelTopology::Phase1>;
581  const auto maxTracks = patatracks_tsoa.view().metadata().size();
582  const int nv = patavtx_soa.view().nvFinal();
583  trkGood.clear();
584  trkGood.reserve(maxTracks);
585  vtxGood.clear();
586  vtxGood.reserve(nv);
587  auto const* quality = patatracks_tsoa.view().quality();
588 
589  // No need to sort either as the algorithms is just using the max (not even the location, just the max value of pt2sum).
590  std::vector<float> pTSquaredSum(nv, 0);
591  std::vector<int> nTrkAssociated(nv, 0);
592 
593  for (int32_t trk_idx = 0; trk_idx < maxTracks; ++trk_idx) {
594  auto nHits = patatrackHelpers::nHits(patatracks_tsoa.view(), trk_idx);
595  if (nHits == 0) {
596  break;
597  }
598  int vtx_ass_to_track = patavtx_soa.view()[trk_idx].idv();
599  if (vtx_ass_to_track >= 0 && vtx_ass_to_track < nv) {
600  auto patatrackPt = patatracks_tsoa.view()[trk_idx].pt();
601  ++nTrkAssociated[vtx_ass_to_track];
602  if (patatrackPt >= trackPtMin_ && patatracks_tsoa.const_view()[trk_idx].chi2() <= trackChi2Max_) {
603  patatrackPt = std::min(patatrackPt, trackPtMax_);
604  pTSquaredSum[vtx_ass_to_track] += patatrackPt * patatrackPt;
605  }
606  }
607  if (nHits > 0 and quality[trk_idx] >= pixelTrack::Quality::loose) {
608  trkGood.push_back(trk_idx);
609  }
610  }
611  if (nv > 0) {
612  const auto minFOM_fromFrac = (*std::max_element(pTSquaredSum.begin(), pTSquaredSum.end())) * fractionSumPt2_;
613  for (int j = nv - 1; j >= 0 && vtxGood.size() < maxVtx_; --j) {
614  auto vtx_idx = patavtx_soa.view()[j].sortInd();
615  assert(vtx_idx < nv);
616  if (nTrkAssociated[vtx_idx] >= 2 && pTSquaredSum[vtx_idx] >= minFOM_fromFrac &&
617  pTSquaredSum[vtx_idx] > minSumPt2_) {
618  vtxGood.push_back(vtx_idx);
619  }
620  }
621  }
622 }
assert(be >=bs)
ConstView const & const_view() const
string quality
maxTracks
Definition: DMR_cfg.py:158
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits

◆ standardizeTensor()

void L2TauNNProducerAlpaka::standardizeTensor ( tensorflow::Tensor &  tensor)
private

Definition at line 371 of file L2TauTagNNProducerAlpaka.cc.

References ztail::d, Exception, input, L2cacheData_, SiStripPI::max, SiStripPI::mean, SiStripPI::min, and L2TauNNProducerAlpakaCacheData::normVec.

Referenced by produce().

371  {
373  std::vector<int> tensor_shape(tensor.shape().dims());
374  for (int d = 0; d < tensor.shape().dims(); d++) {
375  tensor_shape.at(d) = tensor.shape().dim_size(d);
376  }
377  if (tensor_shape.size() != 4) {
378  throw cms::Exception("InvalidTensor") << "Tensor shape does not have 4 dimensions!";
379  }
380  for (int tau_idx = 0; tau_idx < tensor_shape.at(0); tau_idx++) {
381  for (int phi_idx = 0; phi_idx < tensor_shape.at(1); phi_idx++) {
382  for (int eta_idx = 0; eta_idx < tensor_shape.at(2); eta_idx++) {
383  for (int var_idx = 0; var_idx < tensor_shape.at(3); var_idx++) {
384  auto getCell = [&](NNInputs input) -> float& {
385  return getCellImpl(tensor, tau_idx, phi_idx, eta_idx, input);
386  };
387  float mean = L2cacheData_->normVec.at(var_idx).mean;
388  float std = L2cacheData_->normVec.at(var_idx).std;
389  float min = L2cacheData_->normVec.at(var_idx).min;
390  float max = L2cacheData_->normVec.at(var_idx).max;
391  float nonstd_var = getCell(static_cast<NNInputs>(var_idx));
392  float std_var = static_cast<float>((nonstd_var - mean) / std);
393  if (std_var > max) {
394  std_var = static_cast<float>(max);
395  } else if (std_var < min) {
396  std_var = static_cast<float>(min);
397  }
398  getCell(static_cast<NNInputs>(var_idx)) = std_var;
399  }
400  }
401  }
402  }
403 }
const L2TauNNProducerAlpakaCacheData * L2cacheData_
static std::string const input
Definition: EdmProvDump.cc:50
d
Definition: ztail.py:151
std::vector< normDictElement > normVec

Member Data Documentation

◆ beamSpotToken_

const edm::EDGetTokenT<reco::BeamSpot> L2TauNNProducerAlpaka::beamSpotToken_
private

Definition at line 214 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ bFieldToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> L2TauNNProducerAlpaka::bFieldToken_
private

Definition at line 211 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ debugLevel_

const int L2TauNNProducerAlpaka::debugLevel_
private

Definition at line 203 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ dEta_width

constexpr float L2TauNNProducerAlpaka::dEta_width = 2 * L2TauTagNNv1::dR_max / static_cast<float>(L2TauTagNNv1::nCellEta)
static

Definition at line 165 of file L2TauTagNNProducerAlpaka.cc.

Referenced by getEtaPhiIndices().

◆ dPhi_width

constexpr float L2TauNNProducerAlpaka::dPhi_width = 2 * L2TauTagNNv1::dR_max / static_cast<float>(L2TauTagNNv1::nCellPhi)
static

Definition at line 166 of file L2TauTagNNProducerAlpaka.cc.

Referenced by getEtaPhiIndices().

◆ dR2_max

constexpr float L2TauNNProducerAlpaka::dR2_max = L2TauTagNNv1::dR_max * L2TauTagNNv1::dR_max
static

Definition at line 164 of file L2TauTagNNProducerAlpaka.cc.

Referenced by fillCaloRecHits(), and fillPatatracks().

◆ ebToken_

const edm::EDGetTokenT<EcalRecHitCollection> L2TauNNProducerAlpaka::ebToken_
private

Definition at line 208 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ eeToken_

const edm::EDGetTokenT<EcalRecHitCollection> L2TauNNProducerAlpaka::eeToken_
private

Definition at line 209 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ fractionSumPt2_

const float L2TauNNProducerAlpaka::fractionSumPt2_
private

Definition at line 216 of file L2TauTagNNProducerAlpaka.cc.

Referenced by selectGoodTracksAndVertices().

◆ geometryToken_

const edm::ESGetToken<CaloGeometry, CaloGeometryRecord> L2TauNNProducerAlpaka::geometryToken_
private

Definition at line 210 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ hbheToken_

const edm::EDGetTokenT<HBHERecHitCollection> L2TauNNProducerAlpaka::hbheToken_
private

Definition at line 206 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ hoToken_

const edm::EDGetTokenT<HORecHitCollection> L2TauNNProducerAlpaka::hoToken_
private

Definition at line 207 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ inputTensorName_

std::string L2TauNNProducerAlpaka::inputTensorName_
private

Definition at line 221 of file L2TauTagNNProducerAlpaka.cc.

Referenced by L2TauNNProducerAlpaka().

◆ L1TauDesc_

std::vector<InputDescTau> L2TauNNProducerAlpaka::L1TauDesc_
private

Definition at line 205 of file L2TauTagNNProducerAlpaka.cc.

Referenced by L2TauNNProducerAlpaka(), and produce().

◆ L2cacheData_

const L2TauNNProducerAlpakaCacheData* L2TauNNProducerAlpaka::L2cacheData_
private

◆ maxVtx_

const unsigned int L2TauNNProducerAlpaka::maxVtx_
private

Definition at line 215 of file L2TauTagNNProducerAlpaka.cc.

Referenced by selectGoodTracksAndVertices().

◆ minSumPt2_

const float L2TauNNProducerAlpaka::minSumPt2_
private

Definition at line 217 of file L2TauTagNNProducerAlpaka.cc.

Referenced by selectGoodTracksAndVertices().

◆ outputTensorName_

std::string L2TauNNProducerAlpaka::outputTensorName_
private

Definition at line 222 of file L2TauTagNNProducerAlpaka.cc.

Referenced by L2TauNNProducerAlpaka().

◆ pataTracksToken_

const edm::EDGetTokenT<TracksHost> L2TauNNProducerAlpaka::pataTracksToken_
private

Definition at line 213 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ pataVerticesToken_

const edm::EDGetTokenT<ZVertexHost> L2TauNNProducerAlpaka::pataVerticesToken_
private

Definition at line 212 of file L2TauTagNNProducerAlpaka.cc.

Referenced by produce().

◆ tauTriggerToken_

const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> L2TauNNProducerAlpaka::tauTriggerToken_
private

Definition at line 204 of file L2TauTagNNProducerAlpaka.cc.

◆ trackChi2Max_

const float L2TauNNProducerAlpaka::trackChi2Max_
private

Definition at line 220 of file L2TauTagNNProducerAlpaka.cc.

Referenced by selectGoodTracksAndVertices().

◆ trackPtMax_

const float L2TauNNProducerAlpaka::trackPtMax_
private

Definition at line 219 of file L2TauTagNNProducerAlpaka.cc.

Referenced by selectGoodTracksAndVertices().

◆ trackPtMin_

const float L2TauNNProducerAlpaka::trackPtMin_
private

Definition at line 218 of file L2TauTagNNProducerAlpaka.cc.

Referenced by selectGoodTracksAndVertices().