Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-20 01:53:30

0001 /* -*- C++ -*-
0002 
0003 Package: L1CaloTrigger
0004 Class: L1NNCaloTauProducer
0005 Frinedly name: The TauMinator
0006 
0007 \class L1NNCaloTauProducer L1NNCaloTauProducer.cc
0008 
0009 Description: 
0010 Perform reconstruction and identification of tau 
0011 candidates at L1 Trigger with a CNN.
0012 
0013 Implementation:
0014 Take as input the HCAL TPs, the ECAL TPs from
0015 l1tEGammaClusterEmuProducer, and the HGCAL TPs
0016 from l1tHGCalTowerProducer and l1tHGCalBackEndLayer2Producer.
0017 Proceed to clustering of trigger towers in NxM
0018 clusters, match to HGcal 3D clusters in the endcap.
0019 Finally apply the CNNs.
0020 
0021 Original Author: Jona Motta
0022 Created: Tue May 30th 2023
0023 
0024 */
0025 
0026 #include <iostream>
0027 #include <vector>
0028 #include <cmath>
0029 
0030 #include "boost/property_tree/ptree.hpp"
0031 #include "boost/property_tree/json_parser.hpp"
0032 
0033 #include "FWCore/Framework/interface/Frameworkfwd.h"
0034 #include "FWCore/Framework/interface/stream/EDProducer.h"
0035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0036 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0037 #include "FWCore/ServiceRegistry/interface/Service.h"
0038 #include "FWCore/Framework/interface/Event.h"
0039 #include "FWCore/Framework/interface/ESHandle.h"
0040 #include "FWCore/Framework/interface/MakerMacros.h"
0041 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0042 
0043 #include "DataFormats/L1TCalorimeterPhase2/interface/CaloTower.h"
0044 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0045 #include "DataFormats/L1THGCal/interface/HGCalMulticluster.h"
0046 #include "DataFormats/L1TParticleFlow/interface/PFCluster.h"
0047 #include "DataFormats/L1THGCal/interface/HGCalTower.h"
0048 #include "DataFormats/Math/interface/deltaPhi.h"
0049 #include "DataFormats/L1Trigger/interface/Tau.h"
0050 
0051 #include "CalibFormats/CaloTPG/interface/CaloTPGTranscoder.h"
0052 #include "CalibFormats/CaloTPG/interface/CaloTPGRecord.h"
0053 
0054 #include "L1Trigger/L1THGCal/interface/backend/HGCalTriggerClusterIdentificationBase.h"
0055 #include "L1Trigger/Phase2L1ParticleFlow/interface/HGC3DClusterEgID.h"
0056 #include "L1Trigger/L1TCalorimeter/interface/CaloTools.h"
0057 
0058 #include "PhysicsTools/TensorFlow/interface/TensorFlow.h"
0059 
0060 struct NNmodels_GlobalCache {
0061   std::string CNNmodel_CB_path;
0062   std::string DNNident_CB_path;
0063   std::string DNNcalib_CB_path;
0064 
0065   std::string CNNmodel_CE_path;
0066   std::string DNNident_CE_path;
0067   std::string DNNcalib_CE_path;
0068   std::string FeatScaler_CE_path;
0069   boost::property_tree::ptree FeatScaler_CE;
0070 
0071   tensorflow::GraphDef* CNNmodel_CB;
0072   tensorflow::GraphDef* DNNident_CB;
0073   tensorflow::GraphDef* DNNcalib_CB;
0074 
0075   tensorflow::Session* CNNmodel_CBsession;
0076   tensorflow::Session* DNNident_CBsession;
0077   tensorflow::Session* DNNcalib_CBsession;
0078 
0079   tensorflow::GraphDef* CNNmodel_CE;
0080   tensorflow::GraphDef* DNNident_CE;
0081   tensorflow::GraphDef* DNNcalib_CE;
0082 
0083   tensorflow::Session* CNNmodel_CEsession;
0084   tensorflow::Session* DNNident_CEsession;
0085   tensorflow::Session* DNNcalib_CEsession;
0086 };
0087 
0088 class L1NNCaloTauProducer : public edm::stream::EDProducer<edm::GlobalCache<NNmodels_GlobalCache>> {
0089 public:
0090   explicit L1NNCaloTauProducer(const edm::ParameterSet&, const NNmodels_GlobalCache*);
0091 
0092   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0093   static std::unique_ptr<NNmodels_GlobalCache> initializeGlobalCache(const edm::ParameterSet&);
0094   static void globalEndJob(const NNmodels_GlobalCache*) { /*do nothing*/ }
0095 
0096 private:
0097   //----edm control---
0098   void produce(edm::Event&, const edm::EventSetup&) override;
0099 
0100   //----private functions----
0101   int tower_dIPhi(int& iPhi_1, int& iPhi_2) const;
0102   int tower_dIEta(int& iEta_1, int& iEta_2) const;
0103   int endcap_iphi(float& phi) const;
0104   int endcap_ieta(float& eta) const;
0105   float inputQuantizer(float inputF, float LSB, int nbits);
0106   float inputScaler(float inputF, std::string feature);
0107 
0108   //----tokens and handles----
0109   edm::EDGetTokenT<l1tp2::CaloTowerCollection> l1TowersToken;
0110   edm::Handle<l1tp2::CaloTowerCollection> l1CaloTowerHandle;
0111 
0112   edm::EDGetToken hgcalTowersToken;
0113   edm::Handle<l1t::HGCalTowerBxCollection> hgcalTowersHandle;
0114 
0115   edm::EDGetTokenT<l1t::HGCalMulticlusterBxCollection> HGClusterToken;
0116   edm::Handle<l1t::HGCalMulticlusterBxCollection> HGClusterHandle;
0117 
0118   //----private variables----
0119   enum class UseEmInterp { No, EmOnly, AllKeepHad, AllKeepTot };
0120   UseEmInterp scenario;
0121   StringCutObjectSelector<l1t::HGCalMulticluster> preEmId;
0122   l1tpf::HGC3DClusterEgID VsPuId;
0123 
0124   double EcalEtMinForClustering;
0125   double HcalEtMinForClustering;
0126   double EtMinForSeeding;
0127   double EtaRestriction;
0128   double CB_CE_split;
0129 
0130   double IdWp90_CB;
0131   double IdWp95_CB;
0132   double IdWp99_CB;
0133 
0134   double IdWp90_CE;
0135   double IdWp95_CE;
0136   double IdWp99_CE;
0137 
0138   bool DEBUG;
0139 
0140   // hardoced dimensions of the tower clusters
0141   const int seedIdx = 22;
0142   const int IEta_dim = 5;
0143   const int IPhi_dim = 9;
0144   const float Eta_dim = 0.2;
0145   const float Phi_dim = 0.4;
0146   const float Eta_dim_seed = 0.35;
0147   const float Phi_dim_seed = 0.7;
0148   const float Eta_limit = 2.83;
0149 
0150   // classes of objects used only in this producer
0151   class SimpleTowerHit {
0152   public:
0153     float towerEta = -99.;
0154     float towerPhi = -99.;
0155     float towerEm = 0.;
0156     float towerHad = 0.;
0157     float l1egTowerEt = 0.;
0158     float towerEt = 0.;
0159     int towerIeta = -99;
0160     int towerIphi = -99;
0161     bool isBarrel = true;
0162     bool stale = false;
0163     bool stale4seed = false;
0164   };
0165 
0166   class SimpleTowerCluster {
0167   public:
0168     bool barrelSeeded = false;
0169     int seedIeta = -99;
0170     int seedIphi = -99;
0171     float seedEta = -99.;
0172     float seedPhi = -99.;
0173     float rawEt = 0.;
0174     float IDscore = -99.;
0175     float calibPt = -99.;
0176 
0177     std::vector<SimpleTowerHit> towerHits;
0178 
0179     void InitHits(int N, int M) { towerHits.resize(N * M); }
0180   };
0181 
0182   class SimpleHGCluster {
0183   public:
0184     float pt = -99.;
0185     float eta = -99.;
0186     float phi = -99.;
0187     float showerlength = -99.;
0188     float coreshowerlength = -99.;
0189     float spptot = -99.;
0190     float szz = -99.;
0191     float srrtot = -99.;
0192     float meanz = -99.;
0193     bool stale = false;
0194   };
0195 };
0196 
0197 /*
0198 ████████ ██   ██ ██████     ████████  █████  ██   ██ ███    ███ ██ ███    ██  █████  ████████  ██████  ██████  
0199    ██    ██   ██ ██            ██    ██   ██ ██   ██ ████  ████ ██ ████   ██ ██   ██    ██    ██    ██ ██   ██ 
0200    ██    ███████ █████         ██    ███████ ██   ██ ██ ████ ██ ██ ██ ██  ██ ███████    ██    ██    ██ ██████  
0201    ██    ██   ██ ██            ██    ██   ██ ██   ██ ██  ██  ██ ██ ██  ██ ██ ██   ██    ██    ██    ██ ██   ██ 
0202    ██    ██   ██ ██████        ██    ██   ██ ███████ ██      ██ ██ ██   ████ ██   ██    ██     ██████  ██    ██
0203 */
0204 
0205 std::unique_ptr<NNmodels_GlobalCache> L1NNCaloTauProducer::initializeGlobalCache(const edm::ParameterSet& iConfig) {
0206   edm::LogInfo("Initialization") << "Init NN models Global Cache " << std::endl;
0207 
0208   std::unique_ptr<NNmodels_GlobalCache> GlobalCache(new NNmodels_GlobalCache);
0209 
0210   GlobalCache->CNNmodel_CB_path = iConfig.getParameter<std::string>("CNNmodel_CB_path");
0211   GlobalCache->DNNident_CB_path = iConfig.getParameter<std::string>("DNNident_CB_path");
0212   GlobalCache->DNNcalib_CB_path = iConfig.getParameter<std::string>("DNNcalib_CB_path");
0213   GlobalCache->CNNmodel_CE_path = iConfig.getParameter<std::string>("CNNmodel_CE_path");
0214   GlobalCache->DNNident_CE_path = iConfig.getParameter<std::string>("DNNident_CE_path");
0215   GlobalCache->DNNcalib_CE_path = iConfig.getParameter<std::string>("DNNcalib_CE_path");
0216   GlobalCache->FeatScaler_CE_path = iConfig.getParameter<std::string>("FeatScaler_CE_path");
0217 
0218   // Create sessions for Tensorflow inferece
0219   (GlobalCache->CNNmodel_CB) = tensorflow::loadGraphDef(edm::FileInPath((GlobalCache->CNNmodel_CB_path)).fullPath());
0220   (GlobalCache->CNNmodel_CBsession) = tensorflow::createSession((GlobalCache->CNNmodel_CB));
0221 
0222   (GlobalCache->DNNident_CB) = tensorflow::loadGraphDef(edm::FileInPath((GlobalCache->DNNident_CB_path)).fullPath());
0223   (GlobalCache->DNNident_CBsession) = tensorflow::createSession((GlobalCache->DNNident_CB));
0224 
0225   (GlobalCache->DNNcalib_CB) = tensorflow::loadGraphDef(edm::FileInPath((GlobalCache->DNNcalib_CB_path)).fullPath());
0226   (GlobalCache->DNNcalib_CBsession) = tensorflow::createSession((GlobalCache->DNNcalib_CB));
0227 
0228   (GlobalCache->CNNmodel_CE) = tensorflow::loadGraphDef(edm::FileInPath((GlobalCache->CNNmodel_CE_path)).fullPath());
0229   (GlobalCache->CNNmodel_CEsession) = tensorflow::createSession((GlobalCache->CNNmodel_CE));
0230 
0231   (GlobalCache->DNNident_CE) = tensorflow::loadGraphDef(edm::FileInPath((GlobalCache->DNNident_CE_path)).fullPath());
0232   (GlobalCache->DNNident_CEsession) = tensorflow::createSession((GlobalCache->DNNident_CE));
0233 
0234   (GlobalCache->DNNcalib_CE) = tensorflow::loadGraphDef(edm::FileInPath((GlobalCache->DNNcalib_CE_path)).fullPath());
0235   (GlobalCache->DNNcalib_CEsession) = tensorflow::createSession((GlobalCache->DNNcalib_CE));
0236 
0237   // Read features scaler
0238   boost::property_tree::read_json(edm::FileInPath((GlobalCache->FeatScaler_CE_path)).fullPath(),
0239                                   (GlobalCache->FeatScaler_CE));
0240 
0241   return GlobalCache;
0242 }
0243 
0244 // ----Constructor and Destructor -----
0245 L1NNCaloTauProducer::L1NNCaloTauProducer(const edm::ParameterSet& iConfig, const NNmodels_GlobalCache* globalCache)
0246     : l1TowersToken(consumes<l1tp2::CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("l1CaloTowers"))),
0247       hgcalTowersToken(consumes<l1t::HGCalTowerBxCollection>(iConfig.getParameter<edm::InputTag>("hgcalTowers"))),
0248 
0249       HGClusterToken(
0250           consumes<l1t::HGCalMulticlusterBxCollection>(iConfig.getParameter<edm::InputTag>("HgcalClusters"))),
0251       scenario(UseEmInterp::No),
0252       preEmId(iConfig.getParameter<std::string>("preEmId")),
0253       VsPuId(iConfig.getParameter<edm::ParameterSet>("VsPuId")),
0254 
0255       EcalEtMinForClustering(iConfig.getParameter<double>("EcalEtMinForClustering")),
0256       HcalEtMinForClustering(iConfig.getParameter<double>("HcalEtMinForClustering")),
0257       EtMinForSeeding(iConfig.getParameter<double>("EtMinForSeeding")),
0258       EtaRestriction(iConfig.getParameter<double>("EtaRestriction")),
0259       CB_CE_split(iConfig.getParameter<double>("CB_CE_split")),
0260 
0261       IdWp90_CB(iConfig.getParameter<double>("IdWp90_CB")),
0262       IdWp95_CB(iConfig.getParameter<double>("IdWp95_CB")),
0263       IdWp99_CB(iConfig.getParameter<double>("IdWp99_CB")),
0264 
0265       IdWp90_CE(iConfig.getParameter<double>("IdWp90_CE")),
0266       IdWp95_CE(iConfig.getParameter<double>("IdWp95_CE")),
0267       IdWp99_CE(iConfig.getParameter<double>("IdWp99_CE")),
0268 
0269       DEBUG(iConfig.getParameter<bool>("DEBUG")) {
0270   // Initialize HGCAL BDTs
0271   if (!VsPuId.method().empty()) {
0272     VsPuId.prepareTMVA();
0273   }
0274 
0275   // Create produced outputs
0276   produces<BXVector<l1t::Tau>>("L1NNCaloTauCollectionBXV");
0277 
0278   // Settings output
0279   edm::LogInfo("Settings") << "EtaRestriction = " << EtaRestriction << " , CB_CE_split = " << CB_CE_split
0280                            << " , EtMinForSeeding = " << EtMinForSeeding
0281                            << " , HcalTpEtMin = " << HcalEtMinForClustering
0282                            << " , EcalTpEtMin = " << EcalEtMinForClustering << std::endl;
0283 }
0284 
0285 void L1NNCaloTauProducer::produce(edm::Event& iEvent, const edm::EventSetup& eSetup) {
0286   // Output collection
0287   std::unique_ptr<BXVector<l1t::Tau>> L1NNCaloTauCollectionBXV(new l1t::TauBxCollection);
0288 
0289   // Create and Fill collection of all calotowers and their attributes
0290   std::vector<SimpleTowerHit> l1CaloTowers;
0291 
0292   iEvent.getByToken(l1TowersToken, l1CaloTowerHandle);
0293   int warnings = 0;
0294   for (auto& hit : *l1CaloTowerHandle.product()) {
0295     // Skip this weird towers and store warning
0296     if (hit.towerIEta() == -1016 && hit.towerIPhi() == -962) {
0297       warnings += 1;
0298       continue;
0299     }
0300 
0301     SimpleTowerHit l1Hit;
0302     l1Hit.isBarrel = true;
0303     l1Hit.l1egTowerEt = hit.l1egTowerEt();
0304     l1Hit.towerEta = hit.towerEta();
0305     l1Hit.towerPhi = hit.towerPhi();
0306     l1Hit.towerEm = hit.ecalTowerEt();
0307     l1Hit.towerHad = hit.hcalTowerEt();
0308     l1Hit.towerEt = l1Hit.towerEm + l1Hit.towerHad + l1Hit.l1egTowerEt;
0309     l1Hit.towerIeta = hit.towerIEta();
0310     l1Hit.towerIphi = hit.towerIPhi();
0311 
0312     l1CaloTowers.push_back(l1Hit);
0313   }
0314   if (warnings != 0 && DEBUG) {
0315     edm::LogWarning("BrokenTowers") << " ** WARNING : FOUND " << warnings
0316                                     << " TOWERS WITH towerIeta=-1016 AND towerIphi=-962" << std::endl;
0317   }
0318 
0319   iEvent.getByToken(hgcalTowersToken, hgcalTowersHandle);
0320   for (auto& hit : *hgcalTowersHandle.product()) {
0321     SimpleTowerHit l1Hit;
0322     l1Hit.isBarrel = false;
0323     l1Hit.l1egTowerEt = 0.0;
0324     l1Hit.towerEta = hit.eta();
0325     l1Hit.towerPhi = hit.phi();
0326     l1Hit.towerEm = hit.etEm();
0327     l1Hit.towerHad = hit.etHad();
0328     l1Hit.towerEt = l1Hit.towerEm + l1Hit.towerHad;
0329     l1Hit.towerIeta = endcap_ieta(l1Hit.towerEta);  // computed and filled but not used
0330     l1Hit.towerIphi = endcap_iphi(l1Hit.towerPhi);  // computed and filled but not used
0331 
0332     l1CaloTowers.push_back(l1Hit);
0333   }
0334 
0335   // Sort the ECAL+HCAL+L1EGs tower sums based on total ET
0336   std::sort(begin(l1CaloTowers), end(l1CaloTowers), [](const SimpleTowerHit& a, SimpleTowerHit& b) {
0337     return a.towerEt > b.towerEt;
0338   });
0339 
0340   // Create and Fill the collection of 3D clusters and their attributes
0341   std::vector<SimpleHGCluster> AllHGClusters;
0342   iEvent.getByToken(HGClusterToken, HGClusterHandle);
0343 
0344   for (auto cl3dIt = HGClusterHandle->begin(0); cl3dIt != HGClusterHandle->end(0); ++cl3dIt) {
0345     auto& cl3d = *cl3dIt;
0346 
0347     // Implement cl3d PU ID as done in
0348     // https://github.com/cms-sw/cmssw/blob/master/L1Trigger/Phase2L1ParticleFlow/plugins/PFClusterProducerFromHGC3DClusters.cc#L120
0349     bool isEM = preEmId(*cl3dIt);
0350     // FIXME: I suspect the PFCluster and the code handling its energy are not actually needed
0351     l1t::PFCluster cluster(cl3d.pt(), cl3d.eta(), cl3d.phi(), cl3d.hOverE());
0352     if (scenario == UseEmInterp::EmOnly)  // for emID objs, use EM interp as pT and set H = 0
0353     {
0354       if (isEM) {
0355         float pt_new = cl3d.iPt(l1t::HGCalMulticluster::EnergyInterpretation::EM);
0356         float hoe_new = 0.;
0357         cluster = l1t::PFCluster(pt_new, cl3d.eta(), cl3d.phi(), hoe_new, isEM);
0358       }
0359     } else if (scenario == UseEmInterp::AllKeepHad)  // for all objs, replace EM part with EM interp, preserve H
0360     {
0361       float had_old = cl3d.pt() - cluster.emEt();
0362       float em_new = cl3d.iPt(l1t::HGCalMulticluster::EnergyInterpretation::EM);
0363       float pt_new = had_old + em_new;
0364       float hoe_new = em_new > 0 ? (had_old / em_new) : -1;
0365       cluster = l1t::PFCluster(pt_new, cl3d.eta(), cl3d.phi(), hoe_new, isEM);
0366     } else if (scenario == UseEmInterp::AllKeepTot)  // for all objs, replace EM part with EM interp, preserve pT
0367     {
0368       float em_new = cl3d.iPt(l1t::HGCalMulticluster::EnergyInterpretation::EM);
0369       float hoe_new = em_new > 0 ? (cl3d.pt() / em_new - 1) : -1;
0370       cluster = l1t::PFCluster(cl3d.pt(), cl3d.eta(), cl3d.phi(), hoe_new, isEM);
0371     }
0372     float mvaOut = -1;
0373     if (!VsPuId.method().empty()) {
0374       bool id = VsPuId.passID(*cl3dIt, mvaOut);
0375       if (!id) {
0376         continue;
0377       }  // skip cl3d if it does not pass puid
0378     }
0379 
0380     SimpleHGCluster HGCluster;
0381     HGCluster.pt = cl3d.pt();
0382     HGCluster.eta = cl3d.eta();
0383     HGCluster.phi = cl3d.phi();
0384     HGCluster.showerlength = cl3d.showerLength();
0385     HGCluster.coreshowerlength = cl3d.coreShowerLength();
0386     HGCluster.spptot = cl3d.sigmaPhiPhiTot();
0387     HGCluster.szz = cl3d.sigmaZZ();
0388     HGCluster.srrtot = cl3d.sigmaRRTot();
0389     HGCluster.meanz = cl3d.zBarycenter();
0390 
0391     AllHGClusters.push_back(HGCluster);
0392   }
0393 
0394   // Order the collection in pt (the input to the GCT will be pt ordered)
0395   std::sort(begin(AllHGClusters), end(AllHGClusters), [](const SimpleHGCluster& a, SimpleHGCluster& b) {
0396     return a.pt > b.pt;
0397   });
0398 
0399   // Make NxM TowerClusters and HGClusters collections for TauMinator
0400   std::vector<SimpleTowerCluster> l1TowerClustersNxM_CB;
0401   std::vector<SimpleTowerCluster> l1TowerClustersNxM_CE;
0402   std::vector<SimpleHGCluster> HGClusters;
0403 
0404   // Supporting collection of endcap clusters before cl3d matching
0405   std::vector<SimpleTowerCluster> AllL1TowerClustersNxM_CE;
0406 
0407   bool caloTauSeedingFinished = false;
0408   // Loop for seeding of clNxM objects
0409   while (!caloTauSeedingFinished) {
0410     SimpleTowerCluster clNxM;
0411     clNxM.InitHits(IEta_dim, IPhi_dim);
0412     bool seeded = false;
0413 
0414     for (auto& l1CaloTower : l1CaloTowers) {
0415       // Skip seeding in towers that would make the cluster extend in HF
0416       // Skip l1CaloTowers which are already used by this clusters' mask
0417       if (abs(l1CaloTower.towerEta) > Eta_limit || abs(l1CaloTower.towerEta) > EtaRestriction ||
0418           l1CaloTower.stale4seed) {
0419         continue;
0420       }
0421 
0422       // If not seded do the seeding
0423       if (!seeded) {
0424         // The leading unused tower has ET < min, stop jet clustering
0425         if (l1CaloTower.towerEt < EtMinForSeeding) {
0426           caloTauSeedingFinished = true;
0427           continue;
0428         }
0429 
0430         clNxM.seedIeta = l1CaloTower.towerIeta;
0431         clNxM.seedIphi = l1CaloTower.towerIphi;
0432         clNxM.seedEta = l1CaloTower.towerEta;
0433         clNxM.seedPhi = l1CaloTower.towerPhi;
0434         if (l1CaloTower.isBarrel) {
0435           clNxM.barrelSeeded = true;
0436         }
0437 
0438         clNxM.rawEt += l1CaloTower.towerEt;
0439         clNxM.towerHits[seedIdx] = l1CaloTower;
0440         l1CaloTower.stale4seed = true;
0441         l1CaloTower.stale = true;
0442         seeded = true;
0443 
0444         continue;
0445       }
0446 
0447       int d_iEta = 99;
0448       int d_iPhi = 99;
0449       float d_Eta = 99.;
0450       float d_Phi = 99.;
0451       // Ese iEta/iPhi comparisons in the barrel and eta/phi in HGCal
0452       if (clNxM.barrelSeeded && l1CaloTower.isBarrel) {
0453         d_iEta = tower_dIEta(l1CaloTower.towerIeta, clNxM.seedIeta);
0454         d_iPhi = tower_dIPhi(l1CaloTower.towerIphi, clNxM.seedIphi);
0455       } else {
0456         d_Eta = l1CaloTower.towerEta - clNxM.seedEta;
0457         d_Phi = reco::deltaPhi(l1CaloTower.towerPhi, clNxM.seedPhi);
0458       }
0459 
0460       // Stale tower for seeding if it would lead to overalp between clusters
0461       if ((abs(d_iEta) <= IEta_dim - 1 && abs(d_iPhi) <= IPhi_dim - 1) ||
0462           (abs(d_Eta) < Eta_dim_seed && abs(d_Phi) < Phi_dim_seed)) {
0463         l1CaloTower.stale4seed = true;
0464       }
0465 
0466     }  // End for loop over TPs
0467 
0468     // Pushback seeds split in barrel and endcap
0469     if (seeded) {
0470       if (abs(clNxM.seedEta) < CB_CE_split) {
0471         l1TowerClustersNxM_CB.push_back(clNxM);
0472       } else {
0473         AllL1TowerClustersNxM_CE.push_back(clNxM);
0474       }
0475     }
0476 
0477   }  // End while loop of TowerClusters seeding
0478 
0479   // Loop for barrel NxM TowerClusters clustering starting from the seeds
0480   for (auto& clNxM : l1TowerClustersNxM_CB) {
0481     for (auto& l1CaloTower : l1CaloTowers) {
0482       // Skip l1CaloTowers which are already used
0483       if (l1CaloTower.stale) {
0484         continue;
0485       }
0486 
0487       int d_iEta = 99;
0488       int d_iPhi = 99;
0489       float d_Eta = 99.;
0490       float d_Phi = 99.;
0491       int hitIdx = 99.;
0492       // Use iEta/iPhi comparisons in the barrel and use eta/phi in HGCal
0493       if (l1CaloTower.isBarrel) {
0494         d_iEta = tower_dIEta(l1CaloTower.towerIeta, clNxM.seedIeta);
0495         d_iPhi = tower_dIPhi(l1CaloTower.towerIphi, clNxM.seedIphi);
0496 
0497         hitIdx = d_iEta * IPhi_dim + d_iPhi + seedIdx;
0498       } else {
0499         d_Eta = l1CaloTower.towerEta - clNxM.seedEta;
0500         d_Phi = reco::deltaPhi(l1CaloTower.towerPhi, clNxM.seedPhi);
0501 
0502         int dieta = d_Eta / 0.0807;  // minimal difference in endcap is 0.0808
0503         int diphi = d_Phi / 0.0872;
0504         hitIdx = dieta * IPhi_dim + diphi + seedIdx;
0505       }
0506 
0507       // Cluster all towers in a NxM towers mask
0508       if ((abs(d_iEta) <= (IEta_dim - 1) / 2 && abs(d_iPhi) <= (IPhi_dim - 1) / 2) ||
0509           (abs(d_Eta) < Eta_dim && abs(d_Phi) < Phi_dim)) {
0510         clNxM.rawEt += l1CaloTower.towerEt;
0511         clNxM.towerHits[hitIdx] = l1CaloTower;
0512         l1CaloTower.stale = true;
0513       }
0514 
0515     }  // End for loop over TPs
0516 
0517   }  // End while loop of barrel TowerClusters creation
0518 
0519   // In the endcap cross-loop over clNxM and cl3d to match them
0520   // (we can do it before full clustering just using the seed info)
0521   for (auto& clNxM : AllL1TowerClustersNxM_CE) {
0522     bool matched = false;
0523     for (auto& HGCluster : AllHGClusters) {
0524       // In case the clNxM or HGCluster have already been matched just continue through the list to the end
0525       // only use cl3ds above 4GeV
0526       if (matched || HGCluster.stale || HGCluster.pt < 4) {
0527         continue;
0528       }
0529 
0530       float d_Eta = HGCluster.eta - clNxM.seedEta;
0531       float d_Phi = reco::deltaPhi(HGCluster.phi, clNxM.seedPhi);
0532       float d_R2 = pow(d_Eta, 2) + pow(d_Phi, 2);
0533 
0534       if (d_R2 < 0.25) {
0535         HGCluster.stale = true;
0536         HGClusters.push_back(HGCluster);
0537         l1TowerClustersNxM_CE.push_back(clNxM);
0538         matched = true;
0539       }
0540 
0541     }  // End for loop over cl3ds
0542 
0543   }  // End for loop over clNxM
0544 
0545   // Loop for endcap matched NxM TowerClusters clustering starting from the seeds just found
0546   for (auto& clNxM : l1TowerClustersNxM_CE) {
0547     for (auto& l1CaloTower : l1CaloTowers) {
0548       // Skip l1CaloTowers which are already used
0549       if (l1CaloTower.stale) {
0550         continue;
0551       }
0552 
0553       int d_iEta = 99;
0554       int d_iPhi = 99;
0555       float d_Eta = 99.;
0556       float d_Phi = 99.;
0557       int hitIdx = 99.;
0558       // Use iEta/iPhi comparisons in the endcap and use eta/phi in HGCal
0559       if (l1CaloTower.isBarrel) {
0560         d_iEta = tower_dIEta(l1CaloTower.towerIeta, clNxM.seedIeta);
0561         d_iPhi = tower_dIPhi(l1CaloTower.towerIphi, clNxM.seedIphi);
0562 
0563         hitIdx = d_iEta * IPhi_dim + d_iPhi + seedIdx;
0564       } else {
0565         d_Eta = l1CaloTower.towerEta - clNxM.seedEta;
0566         d_Phi = reco::deltaPhi(l1CaloTower.towerPhi, clNxM.seedPhi);
0567 
0568         int dieta = d_Eta / 0.0807;  // minimal difference in endcap is 0.0808
0569         int diphi = d_Phi / 0.0872;
0570         hitIdx = dieta * IPhi_dim + diphi + seedIdx;
0571       }
0572 
0573       // Cluster all towers in a NxM towers mask
0574       if ((abs(d_iEta) <= (IEta_dim - 1) / 2 && abs(d_iPhi) <= (IPhi_dim - 1) / 2) ||
0575           (abs(d_Eta) < Eta_dim && abs(d_Phi) < Phi_dim)) {
0576         clNxM.rawEt += l1CaloTower.towerEt;
0577         clNxM.towerHits[hitIdx] = l1CaloTower;
0578         l1CaloTower.stale = true;
0579       }
0580 
0581     }  // End for loop over TPs
0582 
0583   }  // End while loop of endcap TowerClusters creation
0584 
0585   // Barrel TauMinator application
0586   int batchSize_CB = (int)(l1TowerClustersNxM_CB.size());
0587   tensorflow::TensorShape imageShape_CB({batchSize_CB, IEta_dim, IPhi_dim, 2});
0588   tensorflow::TensorShape positionShape_CB({batchSize_CB, 2});
0589   tensorflow::Tensor TowerClusterImage_CB(tensorflow::DT_FLOAT, imageShape_CB);
0590   tensorflow::Tensor TowerClusterPosition_CB(tensorflow::DT_FLOAT, positionShape_CB);
0591 
0592   int clIdx = 0;
0593   for (auto& clNxM : l1TowerClustersNxM_CB) {
0594     // Fill inputs for Tensorflow inference
0595     for (int eta = 0; eta < IEta_dim; ++eta) {
0596       for (int phi = 0; phi < IPhi_dim; ++phi) {
0597         int towerIdx = eta * IPhi_dim + phi;
0598         TowerClusterImage_CB.tensor<float, 4>()(clIdx, eta, phi, 0) =
0599             inputQuantizer(clNxM.towerHits[towerIdx].l1egTowerEt + clNxM.towerHits[towerIdx].towerEm, 0.25, 10);
0600         TowerClusterImage_CB.tensor<float, 4>()(clIdx, eta, phi, 1) =
0601             inputQuantizer(clNxM.towerHits[towerIdx].towerHad, 0.25, 10);
0602       }
0603     }
0604 
0605     TowerClusterPosition_CB.tensor<float, 2>()(clIdx, 0) = clNxM.seedEta;
0606     TowerClusterPosition_CB.tensor<float, 2>()(clIdx, 1) = clNxM.seedPhi;
0607 
0608     clIdx++;  // Increase batch index
0609   }
0610 
0611   if (batchSize_CB >
0612       0)  // from CMSSW_14_0_X tensorflow does not seem to be able to deal with a tensor of dimension 0 anymore
0613   {
0614     // Apply CNN model
0615     tensorflow::NamedTensorList CNNmodel_CBinputList = {{"TowerClusterImage", TowerClusterImage_CB},
0616                                                         {"TowerClusterPosition", TowerClusterPosition_CB}};
0617     std::vector<tensorflow::Tensor> CNNmodel_CBoutputs;
0618     tensorflow::run((globalCache()->CNNmodel_CBsession),
0619                     CNNmodel_CBinputList,
0620                     {"TauMinator_CB_conv/middleMan/concat"},
0621                     &CNNmodel_CBoutputs);
0622     tensorflow::NamedTensorList DNN_CBinputsList = {{"middleMan", CNNmodel_CBoutputs[0]}};
0623 
0624     // Apply DNN for identification
0625     std::vector<tensorflow::Tensor> DNN_CBoutputsIdent;
0626     tensorflow::run((globalCache()->DNNident_CBsession),
0627                     DNN_CBinputsList,
0628                     {"TauMinator_CB_ident/sigmoid_IDout/Sigmoid"},
0629                     &DNN_CBoutputsIdent);
0630 
0631     // Apply DNN for calibration
0632     std::vector<tensorflow::Tensor> DNN_CBoutputsCalib;
0633     tensorflow::run((globalCache()->DNNcalib_CBsession),
0634                     DNN_CBinputsList,
0635                     {"TauMinator_CB_calib/LIN_DNNout/Relu"},
0636                     &DNN_CBoutputsCalib);
0637 
0638     // Fill TauMinator output variables of TowerClusters
0639     clIdx = 0;
0640     for (auto& clNxM : l1TowerClustersNxM_CB) {
0641       clNxM.IDscore = DNN_CBoutputsIdent[0].matrix<float>()(0, clIdx);
0642       clNxM.calibPt = DNN_CBoutputsCalib[0].matrix<float>()(0, clIdx);
0643       clIdx++;  // Increase batch index
0644     }
0645   }
0646 
0647   // Endcap TauMinator application
0648   int batchSize_CE = (int)(l1TowerClustersNxM_CE.size());
0649   tensorflow::TensorShape imageShape_CE({batchSize_CE, IEta_dim, IPhi_dim, 2});
0650   tensorflow::TensorShape positionShape_CE({batchSize_CE, 2});
0651   tensorflow::TensorShape cl3dfeatShape_CE({batchSize_CE, 8});
0652   tensorflow::Tensor TowerClusterImage_CE(tensorflow::DT_FLOAT, imageShape_CE);
0653   tensorflow::Tensor TowerClusterPosition_CE(tensorflow::DT_FLOAT, positionShape_CE);
0654   tensorflow::Tensor Cl3dShapeFeatures_CE(tensorflow::DT_FLOAT, cl3dfeatShape_CE);
0655 
0656   clIdx = 0;
0657   for (auto& clNxM : l1TowerClustersNxM_CE) {
0658     // Indexing of cl3ds is the same as the one of clNxMs
0659     SimpleHGCluster HGClu = HGClusters[clIdx];
0660 
0661     // Fill inputs for Tensorflow inference
0662     for (int eta = 0; eta < IEta_dim; ++eta) {
0663       for (int phi = 0; phi < IPhi_dim; ++phi) {
0664         int towerIdx = eta * IPhi_dim + phi;
0665         TowerClusterImage_CE.tensor<float, 4>()(clIdx, eta, phi, 0) =
0666             inputQuantizer(clNxM.towerHits[towerIdx].l1egTowerEt + clNxM.towerHits[towerIdx].towerEm, 0.25, 10);
0667         TowerClusterImage_CE.tensor<float, 4>()(clIdx, eta, phi, 1) =
0668             inputQuantizer(clNxM.towerHits[towerIdx].towerHad, 0.25, 10);
0669       }
0670     }
0671 
0672     TowerClusterPosition_CE.tensor<float, 2>()(clIdx, 0) = clNxM.seedEta;
0673     TowerClusterPosition_CE.tensor<float, 2>()(clIdx, 1) = clNxM.seedPhi;
0674 
0675     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 0) = inputScaler(inputQuantizer(HGClu.pt, 0.25, 14), "pt");
0676     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 1) =
0677         inputScaler(inputQuantizer(abs(HGClu.eta) - 1.321, 0.004, 9), "eta");
0678     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 2) = inputScaler(HGClu.showerlength, "showerlength");
0679     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 3) = inputScaler(HGClu.coreshowerlength, "coreshowerlength");
0680     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 4) =
0681         inputScaler(inputQuantizer(HGClu.spptot, 0.0000153, 16), "spptot");
0682     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 5) = inputScaler(inputQuantizer(HGClu.szz, 0.00153, 16), "szz");
0683     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 6) =
0684         inputScaler(inputQuantizer(HGClu.srrtot, 0.0000153, 16), "srrtot");
0685     Cl3dShapeFeatures_CE.tensor<float, 2>()(clIdx, 7) =
0686         inputScaler(inputQuantizer(10 * (abs(HGClu.meanz) - 321.05), 0.5, 12), "meanz");
0687 
0688     clIdx++;  // Increase batch index
0689   }
0690 
0691   if (batchSize_CE >
0692       0)  // from CMSSW_14_0_X tensorflow does not seem to be able to deal with a tensor of dimension 0 anymore
0693   {
0694     // Apply CNN model
0695     tensorflow::NamedTensorList CNNmodel_CEinputList = {{"TowerClusterImage", TowerClusterImage_CE},
0696                                                         {"TowerClusterPosition", TowerClusterPosition_CE},
0697                                                         {"AssociatedCl3dFeatures", Cl3dShapeFeatures_CE}};
0698     std::vector<tensorflow::Tensor> CNNmodel_CEoutputs;
0699     tensorflow::run((globalCache()->CNNmodel_CEsession),
0700                     CNNmodel_CEinputList,
0701                     {"TauMinator_CE_conv/middleMan/concat"},
0702                     &CNNmodel_CEoutputs);
0703     tensorflow::NamedTensorList DNN_CEinputsList = {{"middleMan", CNNmodel_CEoutputs[0]}};
0704 
0705     // Apply DNN for identification
0706     std::vector<tensorflow::Tensor> DNN_CEoutputsIdent;
0707     tensorflow::run((globalCache()->DNNident_CEsession),
0708                     DNN_CEinputsList,
0709                     {"TauMinator_CE_ident/sigmoid_IDout/Sigmoid"},
0710                     &DNN_CEoutputsIdent);
0711 
0712     // Apply DNN for calibration
0713     std::vector<tensorflow::Tensor> DNN_CEoutputsCalib;
0714     tensorflow::run((globalCache()->DNNcalib_CEsession),
0715                     DNN_CEinputsList,
0716                     {"TauMinator_CE_calib/LIN_DNNout/Relu"},
0717                     &DNN_CEoutputsCalib);
0718 
0719     // Fill TauMinator output variables of TowerClusters
0720     clIdx = 0;
0721     for (auto& clNxM : l1TowerClustersNxM_CE) {
0722       clNxM.IDscore = DNN_CEoutputsIdent[0].matrix<float>()(0, clIdx);
0723       clNxM.calibPt = DNN_CEoutputsCalib[0].matrix<float>()(0, clIdx);
0724       clIdx++;  // Increase batch index
0725     }
0726   }
0727 
0728   // Fill the output collection of L1 taus
0729   for (auto& clNxM : l1TowerClustersNxM_CB) {
0730     // Apply eta restriction
0731     if (abs(clNxM.seedEta) > EtaRestriction) {
0732       continue;
0733     }
0734 
0735     // Assign increasing quality to higher scoring candidates
0736     int quality = 0;
0737     // 99% WP
0738     if (clNxM.IDscore > IdWp99_CB) {
0739       quality = 1;
0740     }
0741     // 95% WP
0742     if (clNxM.IDscore > IdWp95_CB) {
0743       quality = 2;
0744     }
0745     // 90% WP
0746     if (clNxM.IDscore > IdWp90_CB) {
0747       quality = 3;
0748     }
0749 
0750     reco::Candidate::PolarLorentzVector tauP4 =
0751         reco::Candidate::PolarLorentzVector(clNxM.calibPt, clNxM.seedEta, clNxM.seedPhi, 0);
0752 
0753     // store ID score multiplied by 10E4 to have good precision even using the Phase1 tau int iso format
0754     // (this is stored just in case for possible additional offline studies)
0755     // tau initialisation =  (p4,    pt,            eta,           phi,           qual,    iso)
0756     l1t::Tau l1Tau = l1t::Tau(tauP4, clNxM.calibPt, clNxM.seedEta, clNxM.seedPhi, quality, clNxM.IDscore * 10E4);
0757     l1Tau.setTowerIEta(clNxM.seedIeta);
0758     l1Tau.setTowerIPhi(clNxM.seedIphi);
0759     l1Tau.setRawEt(clNxM.rawEt);
0760 
0761     L1NNCaloTauCollectionBXV->push_back(0, l1Tau);
0762   }
0763 
0764   for (auto& clNxM : l1TowerClustersNxM_CE) {
0765     // Apply eta restriction
0766     if (abs(clNxM.seedEta) > EtaRestriction) {
0767       continue;
0768     }
0769 
0770     // Assign increasing quality to higher scoring candidates
0771     int quality = 0;
0772     // 99% WP
0773     if (clNxM.IDscore > IdWp99_CE) {
0774       quality = 1;
0775     }
0776     // 95% WP
0777     if (clNxM.IDscore > IdWp95_CE) {
0778       quality = 2;
0779     }
0780     // 90% WP
0781     if (clNxM.IDscore > IdWp90_CE) {
0782       quality = 3;
0783     }
0784 
0785     reco::Candidate::PolarLorentzVector tauP4 =
0786         reco::Candidate::PolarLorentzVector(clNxM.calibPt, clNxM.seedEta, clNxM.seedPhi, 0);
0787 
0788     // store ID score multiplied by 10E4 to have good precision even using the Phase1 tau int iso format
0789     // (this is stored just in case for possible additional offline studies)
0790     // tau initialisation =  (p4,    pt,            eta,           phi,           qual,    iso)
0791     l1t::Tau l1Tau = l1t::Tau(tauP4, clNxM.calibPt, clNxM.seedEta, clNxM.seedPhi, quality, clNxM.IDscore * 10E4);
0792     l1Tau.setTowerIEta(clNxM.seedIeta);
0793     l1Tau.setTowerIPhi(clNxM.seedIphi);
0794     l1Tau.setRawEt(clNxM.rawEt);
0795 
0796     L1NNCaloTauCollectionBXV->push_back(0, l1Tau);
0797   }
0798 
0799   // Fill output
0800   iEvent.put(std::move(L1NNCaloTauCollectionBXV), "L1NNCaloTauCollectionBXV");
0801 
0802 }  // End of produce function
0803 
0804 int L1NNCaloTauProducer::tower_dIPhi(int& iPhi_1, int& iPhi_2) const {
0805   const int PI = 36;
0806   int result = iPhi_1 - iPhi_2;
0807   if (result > PI) {
0808     result -= 2 * PI;
0809   }
0810   if (result <= -PI) {
0811     result += 2 * PI;
0812   }
0813   return result;
0814 }
0815 
0816 int L1NNCaloTauProducer::tower_dIEta(int& iEta_1, int& iEta_2) const {
0817   if (iEta_1 * iEta_2 > 0) {
0818     return iEta_1 - iEta_2;
0819   } else {
0820     if (iEta_1 > 0) {
0821       return iEta_1 - iEta_2 - 1;
0822     } else {
0823       return iEta_1 - iEta_2 + 1;
0824     }
0825   }
0826 }
0827 
0828 int L1NNCaloTauProducer::endcap_iphi(float& phi) const {
0829   const float phi_step = 0.0872664;
0830   if (phi > 0) {
0831     return floor(phi / phi_step) + 1;
0832   } else {
0833     return floor(phi / phi_step) + 73;
0834   }
0835 }
0836 
0837 int L1NNCaloTauProducer::endcap_ieta(float& eta) const {
0838   const float eta_step = 0.0845;
0839   return floor(abs(eta) / eta_step) * std::copysign(1, eta);
0840 }
0841 
0842 float L1NNCaloTauProducer::inputQuantizer(float inputF, float LSB, int nbits) {
0843   return min(floor(inputF / LSB), float(pow(2, nbits) - 1)) * LSB;
0844 }
0845 
0846 float L1NNCaloTauProducer::inputScaler(float inputF, std::string feature) {
0847   float mean = (globalCache()->FeatScaler_CE).get_child(feature).get<float>("mean");
0848   float std = (globalCache()->FeatScaler_CE).get_child(feature).get<float>("std");
0849 
0850   return (inputF - mean) / std;
0851 }
0852 
0853 void L1NNCaloTauProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0854   edm::ParameterSetDescription desc;
0855 
0856   desc.add<edm::InputTag>("l1CaloTowers", edm::InputTag("l1tEGammaClusterEmuProducer", "L1CaloTowerCollection"));
0857   desc.add<edm::InputTag>("hgcalTowers", edm::InputTag("l1tHGCalTowerProducer", "HGCalTowerProcessor"));
0858   desc.add<edm::InputTag>("HgcalClusters",
0859                           edm::InputTag("l1tHGCalBackEndLayer2Producer", "HGCalBackendLayer2Processor3DClustering"));
0860 
0861   desc.add<std::string>("preEmId", "hOverE < 0.3 && hOverE >= 0");
0862   {
0863     edm::ParameterSetDescription psd0;
0864     psd0.add<bool>("isPUFilter", true);
0865     psd0.add<std::string>("preselection", "");
0866     psd0.add<std::string>("method", "BDT");
0867     {
0868       edm::ParameterSetDescription vpsd2;
0869       vpsd2.add<std::string>("name");
0870       vpsd2.add<std::string>("value");
0871       std::vector<edm::ParameterSet> temp2;
0872       temp2.reserve(5);
0873       {
0874         edm::ParameterSet temp3;
0875         temp3.addParameter<std::string>("name", "eMax");
0876         temp3.addParameter<std::string>("value", "eMax()");
0877         temp2.push_back(temp3);
0878       }
0879       {
0880         edm::ParameterSet temp3;
0881         temp3.addParameter<std::string>("name", "eMaxOverE");
0882         temp3.addParameter<std::string>("value", "eMax()/energy()");
0883         temp2.push_back(temp3);
0884       }
0885       {
0886         edm::ParameterSet temp3;
0887         temp3.addParameter<std::string>("name", "sigmaPhiPhiTot");
0888         temp3.addParameter<std::string>("value", "sigmaPhiPhiTot()");
0889         temp2.push_back(temp3);
0890       }
0891       {
0892         edm::ParameterSet temp3;
0893         temp3.addParameter<std::string>("name", "sigmaRRTot");
0894         temp3.addParameter<std::string>("value", "sigmaRRTot()");
0895         temp2.push_back(temp3);
0896       }
0897       {
0898         edm::ParameterSet temp3;
0899         temp3.addParameter<std::string>("name", "triggerCells90percent");
0900         temp3.addParameter<std::string>("value", "triggerCells90percent()");
0901         temp2.push_back(temp3);
0902       }
0903       psd0.addVPSet("variables", vpsd2, temp2);
0904     }
0905     psd0.add<std::string>(
0906         "weightsFile", "L1Trigger/Phase2L1ParticleFlow/data/hgcal_egID/Photon_Pion_vs_Neutrino_BDTweights_1116.xml.gz");
0907     psd0.add<std::string>("wp", "-0.10");
0908     desc.add<edm::ParameterSetDescription>("VsPuId", psd0);
0909   }
0910 
0911   desc.add<double>("EcalEtMinForClustering", 0.0);
0912   desc.add<double>("HcalEtMinForClustering", 0.0);
0913   desc.add<double>("EtMinForSeeding", 2.5);
0914   desc.add<double>("EtaRestriction", 2.4);
0915   desc.add<double>("CB_CE_split", 1.55);
0916 
0917   desc.add<std::string>("CNNmodel_CB_path", "L1Trigger/L1CaloTrigger/data/Phase2_NNCaloTaus/v22/CNNmodel_CB.pb");
0918   desc.add<std::string>("DNNident_CB_path", "L1Trigger/L1CaloTrigger/data/Phase2_NNCaloTaus/v22/DNNident_CB.pb");
0919   desc.add<std::string>("DNNcalib_CB_path", "L1Trigger/L1CaloTrigger/data/Phase2_NNCaloTaus/v22/DNNcalib_CB.pb");
0920   desc.add<std::string>("CNNmodel_CE_path", "L1Trigger/L1CaloTrigger/data/Phase2_NNCaloTaus/v22/CNNmodel_CE.pb");
0921   desc.add<std::string>("DNNident_CE_path", "L1Trigger/L1CaloTrigger/data/Phase2_NNCaloTaus/v22/DNNident_CE.pb");
0922   desc.add<std::string>("DNNcalib_CE_path", "L1Trigger/L1CaloTrigger/data/Phase2_NNCaloTaus/v22/DNNcalib_CE.pb");
0923   desc.add<std::string>("FeatScaler_CE_path", "L1Trigger/L1CaloTrigger/data/Phase2_NNCaloTaus/Cl3dFeatScaler_CE.json");
0924 
0925   desc.add<double>("IdWp90_CB", 0.706);
0926   desc.add<double>("IdWp95_CB", 0.3432);
0927   desc.add<double>("IdWp99_CB", 0.0337);
0928   desc.add<double>("IdWp90_CE", 0.5711);
0929   desc.add<double>("IdWp95_CE", 0.2742);
0930   desc.add<double>("IdWp99_CE", 0.0394);
0931 
0932   desc.add<bool>("DEBUG", false);
0933 
0934   descriptions.add("l1tNNCaloTauProducer", desc);
0935 }
0936 
0937 DEFINE_FWK_MODULE(L1NNCaloTauProducer);