Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:10

0001 #include "RecoEgamma/PhotonIdentification/interface/PhotonDNNEstimator.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "FWCore/Utilities/interface/FileInPath.h"
0004 
0005 #include <iostream>
0006 #include <fstream>
0007 #include <memory>
0008 
0009 using namespace std::placeholders;
0010 
0011 inline uint photonModelSelector(const std::map<std::string, float>& vars, float etaThr) {
0012   /* 
0013   Selection of the model to be applied on the photon based on eta limit
0014   */
0015   const auto absEta = std::abs(vars.at("eta"));
0016   if (absEta <= etaThr) {
0017     return 0;
0018   } else {
0019     return 1;
0020   }
0021 }
0022 
0023 PhotonDNNEstimator::PhotonDNNEstimator(const egammaTools::DNNConfiguration& cfg, const bool useEBModelInGap)
0024     : dnnHelper_(cfg,
0025                  std::bind(photonModelSelector,
0026                            _1,
0027                            (useEBModelInGap) ? PhotonDNNEstimator::ecalBarrelMaxEtaWithGap
0028                                              : PhotonDNNEstimator::ecalBarrelMaxEtaNoGap),
0029                  PhotonDNNEstimator::dnnAvaibleInputs),
0030       useEBModelInGap_(useEBModelInGap) {}
0031 
0032 std::vector<tensorflow::Session*> PhotonDNNEstimator::getSessions() const { return dnnHelper_.getSessions(); };
0033 
0034 const std::vector<std::string> PhotonDNNEstimator::dnnAvaibleInputs = {{"pt",
0035                                                                         "eta",
0036                                                                         "hadTowOverEm",
0037                                                                         "trkSumPtHollowConeDR03",
0038                                                                         "EcalRecHit",
0039                                                                         "SigmaIetaIeta",
0040                                                                         "SigmaIetaIetaFull5x5",
0041                                                                         "SigmaIEtaIPhiFull5x5",
0042                                                                         "EcalPFClusterIso",
0043                                                                         "HcalPFClusterIso",
0044                                                                         "HasPixelSeed",
0045                                                                         "R9Full5x5",
0046                                                                         "hcalTower"}};
0047 
0048 std::map<std::string, float> PhotonDNNEstimator::getInputsVars(const reco::Photon& photon) const {
0049   // Prepare a map with all the defined variables
0050   std::map<std::string, float> variables;
0051   variables["pt"] = photon.pt();
0052   variables["eta"] = photon.eta();
0053   variables["hadTowOverEm"] = photon.hadTowOverEmValid() ? photon.hadTowOverEm() : 0;
0054   variables["trkSumPtHollowConeDR03"] = photon.trkSumPtHollowConeDR03();
0055   variables["EcalRecHit"] = photon.ecalRecHitSumEtConeDR03();
0056   variables["SigmaIetaIeta"] = photon.sigmaIetaIeta();
0057   variables["SigmaIetaIetaFull5x5"] = photon.full5x5_sigmaIetaIeta();
0058   variables["SigmaIEtaIPhiFull5x5"] = photon.full5x5_showerShapeVariables().sigmaIetaIphi;
0059   variables["EcalPFClusterIso"] = photon.ecalPFClusterIso();
0060   variables["HcalPFClusterIso"] = photon.hcalPFClusterIso();
0061   variables["HasPixelSeed"] = (Int_t)photon.hasPixelSeed();
0062   variables["R9Full5x5"] = photon.full5x5_r9();
0063   variables["hcalTower"] = photon.hcalTowerSumEtConeDR03();
0064   variables["R9Full5x5"] = photon.full5x5_r9();
0065   // Define more variables here and use them directly in the model config!
0066   return variables;
0067 }
0068 
0069 std::vector<std::pair<uint, std::vector<float>>> PhotonDNNEstimator::evaluate(
0070     const reco::PhotonCollection& photons, const std::vector<tensorflow::Session*>& sessions) const {
0071   // Collect the map of variables for each candidate and call the dnnHelper
0072   // Scaling, model selection and running is performed in the helper
0073   std::vector<std::map<std::string, float>> inputs;
0074   for (const auto& photon : photons) {
0075     inputs.push_back(getInputsVars(photon));
0076   }
0077   return dnnHelper_.evaluate(inputs, sessions);
0078 }