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 288 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.

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

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

Referenced by produce().

325  {
327  std::vector<int> tensor_shape(tensor.shape().dims());
328  for (int d = 0; d < tensor.shape().dims(); d++) {
329  tensor_shape.at(d) = tensor.shape().dim_size(d);
330  }
331  if (tensor_shape.size() != 4) {
332  throw cms::Exception("InvalidTensor") << "Tensor shape does not have 4 dimensions!";
333  }
334  for (int tau_idx = 0; tau_idx < tensor_shape.at(0); tau_idx++) {
335  for (int phi_idx = 0; phi_idx < tensor_shape.at(1); phi_idx++) {
336  for (int eta_idx = 0; eta_idx < tensor_shape.at(2); eta_idx++) {
337  for (int var_idx = 0; var_idx < tensor_shape.at(3); var_idx++) {
338  auto getCell = [&](NNInputs input) -> float& {
339  return getCellImpl(tensor, tau_idx, phi_idx, eta_idx, input);
340  };
341  auto nonstd_var = getCell(static_cast<NNInputs>(var_idx));
342  if (edm::isNotFinite(nonstd_var)) {
343  edm::LogWarning("InputVar") << "var is nan \nvar name= "
344  << L2TauTagNNv1::varNameMap.at(static_cast<L2TauTagNNv1::NNInputs>(var_idx))
345  << "\t var_idx = " << var_idx << "\t eta_idx = " << eta_idx
346  << "\t phi_idx = " << phi_idx << "\t tau_idx = " << tau_idx;
347  if (debugLevel > 2) {
348  edm::LogWarning("InputVar") << "other vars in same cell \n";
349  if (var_idx + 1 < tensor_shape.at(3))
350  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 1))
351  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 1));
352  if (var_idx + 2 < tensor_shape.at(3))
353  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 2))
354  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 2));
355  if (var_idx + 3 < tensor_shape.at(3))
356  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 3))
357  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 3));
358  if (var_idx + 4 < tensor_shape.at(3))
359  edm::LogWarning("InputVar") << L2TauTagNNv1::varNameMap.at(static_cast<NNInputs>(var_idx + 4))
360  << "\t = " << getCell(static_cast<NNInputs>(var_idx + 4));
361  }
362  }
363  }
364  }
365  }
366  }
367 }
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 435 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().

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

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

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

Referenced by produce().

403  {
405 
406  const int nTaus = allTaus.size();
407  for (int tau_idx = 0; tau_idx < nTaus; tau_idx++) {
408  for (int eta_idx = 0; eta_idx < L2TauTagNNv1::nCellEta; eta_idx++) {
409  for (int phi_idx = 0; phi_idx < L2TauTagNNv1::nCellPhi; phi_idx++) {
410  auto getCell = [&](NNInputs input) -> float& {
411  return getCellImpl(cellGridMatrix, tau_idx, phi_idx, eta_idx, input);
412  };
413  getCell(NNInputs::l1Tau_pt) = allTaus[tau_idx]->pt();
414  getCell(NNInputs::l1Tau_eta) = allTaus[tau_idx]->eta();
415  getCell(NNInputs::l1Tau_hwIso) = allTaus[tau_idx]->hwIso();
416  }
417  }
418  }
419 }
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 653 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().

658  {
660  using patatrackHelpers = TracksUtilities<pixelTopology::Phase1>;
661  float deta, dphi;
662  int eta_idx = 0;
663  int phi_idx = 0;
664  int tau_idx = 0;
665 
666  auto getCell = [&](NNInputs input) -> float& {
667  return getCellImpl(cellGridMatrix, tau_idx, phi_idx, eta_idx, input);
668  };
669 
670  std::vector<int> trkGood;
671  std::vector<int> vtxGood;
672 
673  selectGoodTracksAndVertices(patavtx_soa, patatracks_tsoa, trkGood, vtxGood);
674 
675  const int nTaus = allTaus.size();
676  for (tau_idx = 0; tau_idx < nTaus; tau_idx++) {
677  const float tauEta = allTaus[tau_idx]->eta();
678  const float tauPhi = allTaus[tau_idx]->phi();
679 
680  for (const auto it : trkGood) {
681  const float patatrackPt = patatracks_tsoa.const_view()[it].pt();
682  if (patatrackPt <= 0)
683  continue;
684  const float patatrackPhi = reco::phi(patatracks_tsoa.const_view(), it);
685  const float patatrackEta = patatracks_tsoa.const_view()[it].eta();
686  const float patatrackCharge = reco::charge(patatracks_tsoa.const_view(), it);
687  const float patatrackChi2OverNdof = patatracks_tsoa.view()[it].chi2();
688  const auto nHits = patatrackHelpers::nHits(patatracks_tsoa.const_view(), it);
689  if (nHits <= 0)
690  continue;
691  const int patatrackNdof = 2 * std::min(6, nHits) - 5;
692 
693  const int vtx_idx_assTrk = patavtx_soa.view()[it].idv();
694  if (reco::deltaR2(patatrackEta, patatrackPhi, tauEta, tauPhi) < dR2_max) {
695  std::tie(deta, dphi, eta_idx, phi_idx) =
696  getEtaPhiIndices(patatrackEta, patatrackPhi, allTaus[tau_idx]->polarP4());
697  getCell(NNInputs::PatatrackPtSum) += patatrackPt;
698  getCell(NNInputs::PatatrackSize) += 1.;
699  getCell(NNInputs::PatatrackChargeSum) += patatrackCharge;
700  getCell(NNInputs::PatatrackDeltaEta) += deta * patatrackPt;
701  getCell(NNInputs::PatatrackDeltaPhi) += dphi * patatrackPt;
702  getCell(NNInputs::PatatrackChi2OverNdof) += patatrackChi2OverNdof * patatrackPt;
703  getCell(NNInputs::PatatrackNdof) += patatrackNdof * patatrackPt;
704  std::pair<float, float> impactParameters = impactParameter(it, patatracks_tsoa, patatrackPhi, beamspot, magfi);
705  getCell(NNInputs::PatatrackDxy) += impactParameters.first * patatrackPt;
706  getCell(NNInputs::PatatrackDz) += impactParameters.second * patatrackPt;
707  if ((std::find(vtxGood.begin(), vtxGood.end(), vtx_idx_assTrk) != vtxGood.end())) {
708  getCell(NNInputs::PatatrackPtSumWithVertex) += patatrackPt;
709  getCell(NNInputs::PatatrackSizeWithVertex) += 1.;
710  }
711  }
712  }
713 
714  // normalize to sum and define stdDev
715  for (eta_idx = 0; eta_idx < L2TauTagNNv1::nCellEta; eta_idx++) {
716  for (phi_idx = 0; phi_idx < L2TauTagNNv1::nCellPhi; phi_idx++) {
717  getCell(NNInputs::nVertices) = vtxGood.size();
718  if (getCell(NNInputs::PatatrackPtSum) > 0.) {
719  getCell(NNInputs::PatatrackDeltaEta) /= getCell(NNInputs::PatatrackPtSum);
720  getCell(NNInputs::PatatrackDeltaPhi) /= getCell(NNInputs::PatatrackPtSum);
721  getCell(NNInputs::PatatrackChi2OverNdof) /= getCell(NNInputs::PatatrackPtSum);
722  getCell(NNInputs::PatatrackNdof) /= getCell(NNInputs::PatatrackPtSum);
723  getCell(NNInputs::PatatrackDxy) /= getCell(NNInputs::PatatrackPtSum);
724  getCell(NNInputs::PatatrackDz) /= getCell(NNInputs::PatatrackPtSum);
725  }
726  }
727  }
728  }
729 }
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 431 of file L2TauTagNNProducerAlpaka.cc.

References position.

Referenced by fillCaloRecHits(), and fillPatatracks().

431  {
432  return getEtaPhiIndices(position.eta(), position.phi(), tau_p4);
433 }
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 422 of file L2TauTagNNProducerAlpaka.cc.

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

422  {
423  const float deta = eta - tau_p4.eta();
424  const float dphi = reco::deltaPhi(phi, tau_p4.phi());
425  const int eta_idx = static_cast<int>(floor((deta + L2TauTagNNv1::dR_max) / dEta_width));
426  const int phi_idx = static_cast<int>(floor((dphi + L2TauTagNNv1::dR_max) / dPhi_width));
427  return std::make_tuple(deta, dphi, eta_idx, phi_idx);
428 }
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 731 of file L2TauTagNNProducerAlpaka.cc.

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

Referenced by produce().

731  {
732  const int nTau = cellGridMatrix.shape().dim_size(0);
733  if (nTau == 0) {
734  return std::vector<float>();
735  } else {
736  std::vector<tensorflow::Tensor> pred_tensor;
737  tensorflow::run(L2cacheData_->session, {{inputTensorName_, cellGridMatrix}}, {outputTensorName_}, &pred_tensor);
738  std::vector<float> pred_vector(nTau);
739  for (int tau_idx = 0; tau_idx < nTau; ++tau_idx) {
740  pred_vector[tau_idx] = pred_tensor[0].matrix<float>()(tau_idx, 0);
741  }
742 
743  return pred_vector;
744  }
745 }
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:271

◆ globalEndJob()

void L2TauNNProducerAlpaka::globalEndJob ( L2TauNNProducerAlpakaCacheData cacheData)
static

Definition at line 250 of file L2TauTagNNProducerAlpaka.cc.

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

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

◆ 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 622 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().

626  {
627  /* dxy and dz */
628  riemannFit::Vector5d ipar, opar;
629  riemannFit::Matrix5d icov, ocov;
630  TracksUtilities<pixelTopology::Phase1>::copyToDense(patatracks_tsoa.view(), ipar, icov, it);
631  riemannFit::transformToPerigeePlane(ipar, icov, opar, ocov);
632  LocalTrajectoryParameters lpar(opar(0), opar(1), opar(2), opar(3), opar(4), 1.);
633  float sp = std::sin(patatrackPhi);
634  float cp = std::cos(patatrackPhi);
635  Surface::RotationType Rotation(sp, -cp, 0, 0, 0, -1.f, cp, sp, 0);
636  GlobalPoint BeamSpotPoint(beamspot.x0(), beamspot.y0(), beamspot.z0());
637  Plane impPointPlane(BeamSpotPoint, Rotation);
639  impPointPlane.toGlobal(lpar.position()), impPointPlane.toGlobal(lpar.momentum()), lpar.charge(), magfi);
640  GlobalPoint vv = gp.position();
641  math::XYZPoint pos(vv.x(), vv.y(), vv.z());
642  GlobalVector pp = gp.momentum();
643  math::XYZVector mom(pp.x(), pp.y(), pp.z());
644  auto lambda = M_PI_2 - pp.theta();
645  auto phi = pp.phi();
646  float patatrackDxy = -vv.x() * std::sin(phi) + vv.y() * std::cos(phi);
647  float patatrackDz =
648  (vv.z() * std::cos(lambda) - (vv.x() * std::cos(phi) + vv.y() * std::sin(phi)) * std::sin(lambda)) /
649  std::cos(lambda);
650  return std::make_pair(patatrackDxy, patatrackDz);
651 }
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_2024v14_cff::graphPath, submitPVResolutionJobs::key, tensorflow::loadGraphDef(), normDictElement::max, normDictElement::mean, normDictElement::min, HLT_2024v14_cff::normalizationDict, L2TauTagNNv1::nVars, 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 
236  boost::property_tree::ptree loadPtreeRoot;
237  auto const normalizationDict = edm::FileInPath(cfg.getParameter<std::string>("normalizationDict")).fullPath();
238  boost::property_tree::read_json(normalizationDict, loadPtreeRoot);
239  for (const auto& [key, val] : L2TauTagNNv1::varNameMap) {
240  boost::property_tree::ptree var = loadPtreeRoot.get_child(val);
241  normDictElement current_element;
242  current_element.mean = var.get_child("mean").get_value<float>();
243  current_element.std = var.get_child("std").get_value<float>();
244  current_element.min = var.get_child("min").get_value<float>();
245  current_element.max = var.get_child("max").get_value<float>();
246  cacheData->normVec.push_back(current_element);
247  }
248  return cacheData;
249 }
const std::map< NNInputs, std::string > varNameMap
std::string fullPath() const
Definition: FileInPath.cc:161
GraphDef * loadGraphDef(const std::string &pbFile)
Definition: TensorFlow.cc:119
constexpr int nVars
key
prepare the HTCondor submission files and eventually submit them
Session * createSession()
Definition: TensorFlow.cc:136

◆ produce()

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

Definition at line 747 of file L2TauTagNNProducerAlpaka.cc.

References beamSpotToken_, bFieldToken_, checknan(), debugLevel_, HLT_2024v14_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.

747  {
748  std::vector<std::vector<size_t>> TauCollectionMap(L1TauDesc_.size());
749  l1t::TauVectorRef allTaus;
750 
751  for (size_t inp_idx = 0; inp_idx < L1TauDesc_.size(); inp_idx++) {
752  l1t::TauVectorRef l1Taus;
753  auto const& l1TriggeredTaus = event.get(L1TauDesc_[inp_idx].inputToken_);
754  l1TriggeredTaus.getObjects(trigger::TriggerL1Tau, l1Taus);
755  TauCollectionMap.at(inp_idx).resize(l1Taus.size());
756 
757  for (size_t l1_idx = 0; l1_idx < l1Taus.size(); l1_idx++) {
758  size_t tau_idx;
759  const auto iter = std::find(allTaus.begin(), allTaus.end(), l1Taus[l1_idx]);
760  if (iter != allTaus.end()) {
761  tau_idx = std::distance(allTaus.begin(), iter);
762  } else {
763  allTaus.push_back(l1Taus[l1_idx]);
764  tau_idx = allTaus.size() - 1;
765  }
766  TauCollectionMap.at(inp_idx).at(l1_idx) = tau_idx;
767  }
768  }
769  const auto ebCal = event.getHandle(ebToken_);
770  const auto eeCal = event.getHandle(eeToken_);
771  const auto hbhe = event.getHandle(hbheToken_);
772  const auto ho = event.getHandle(hoToken_);
773  auto const& patatracks_SoA = event.get(pataTracksToken_);
774  auto const& vertices_SoA = event.get(pataVerticesToken_);
775  const auto bsHandle = event.getHandle(beamSpotToken_);
776 
777  auto const fieldESH = eventsetup.getHandle(bFieldToken_);
778  auto const geometry = eventsetup.getHandle(geometryToken_);
779 
780  caloRecHitCollections caloRecHits;
781  caloRecHits.hbhe = &*hbhe;
782  caloRecHits.ho = &*ho;
783  caloRecHits.eb = &*ebCal;
784  caloRecHits.ee = &*eeCal;
785  caloRecHits.geometry = &*geometry;
786 
787  const int nTaus = allTaus.size();
788  tensorflow::Tensor cellGridMatrix(tensorflow::DT_FLOAT,
790  const int n_inputs = nTaus * L2TauTagNNv1::nCellEta * L2TauTagNNv1::nCellPhi * L2TauTagNNv1::nVars;
791  for (int input_idx = 0; input_idx < n_inputs; ++input_idx) {
792  cellGridMatrix.flat<float>()(input_idx) = 0;
793  }
794  fillL1TauVars(cellGridMatrix, allTaus);
795 
796  fillCaloRecHits(cellGridMatrix, allTaus, caloRecHits);
797 
798  fillPatatracks(cellGridMatrix, allTaus, patatracks_SoA, vertices_SoA, *bsHandle, fieldESH.product());
799 
800  standardizeTensor(cellGridMatrix);
801 
802  if (debugLevel_ > 0) {
803  checknan(cellGridMatrix, debugLevel_);
804  }
805 
806  std::vector<float> tau_score = getTauScore(cellGridMatrix);
807 
808  for (size_t inp_idx = 0; inp_idx < L1TauDesc_.size(); inp_idx++) {
809  const size_t nTau = TauCollectionMap[inp_idx].size();
810  auto tau_tags = std::make_unique<std::vector<float>>(nTau);
811  for (size_t tau_pos = 0; tau_pos < nTau; ++tau_pos) {
812  const auto tau_idx = TauCollectionMap[inp_idx][tau_pos];
813  if (debugLevel_ > 0) {
814  edm::LogInfo("DebugInfo") << event.id().event() << " \t " << (allTaus[tau_idx])->pt() << " \t "
815  << tau_score.at(tau_idx) << std::endl;
816  }
817  (*tau_tags)[tau_pos] = tau_score.at(tau_idx);
818  }
819  event.put(std::move(tau_tags), L1TauDesc_[inp_idx].CollectionName);
820  }
821 }
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 574 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().

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

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

Referenced by produce().

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