10 #include <Math/VectorUtil.h> 14 #include "tensorflow/core/util/memmapped_file_system.h" 38 #include "oneapi/tbb/concurrent_unordered_set.h" 52 using GraphPtr = std::shared_ptr<tensorflow::GraphDef>;
55 for (
const auto& graph_entry : graph_names) {
64 if (!mmap_status.ok()) {
65 throw cms::Exception(
"DeepTauCache: unable to initalize memmapped environment for ")
67 << mmap_status.ToString();
70 graphs_[entry_name] = std::make_unique<tensorflow::GraphDef>();
73 tensorflow::MemmappedFileSystem::kMemmappedPackageDefaultGraphDef,
75 if (!load_graph_status.ok())
77 << load_graph_status.ToString();
79 options.getSessionOptions().config.mutable_graph_options()->mutable_optimizer_options()->set_opt_level(
80 ::tensorflow::OptimizerOptions::L0);
104 std::map<std::string, std::unique_ptr<tensorflow::MemmappedEnv>>
memmappedEnv_;
113 constexpr
int number_of_inner_cell = 11;
114 constexpr
int number_of_outer_cell = 21;
115 constexpr
int number_of_conv_features = 64;
116 namespace TauBlockInputs {
125 tau_n_charged_prongs,
126 tau_n_neutral_prongs,
128 chargedIsoPtSumdR03_over_dR05,
131 neutralIsoPtSumWeight_over_neutralIsoPtSum,
132 neutralIsoPtSumWeightdR03_over_neutralIsoPtSum,
133 neutralIsoPtSumdR03_over_dR05,
134 photonPtSumOutsideSignalCone,
151 tau_flightLength_sig,
152 tau_pt_weighted_deta_strip,
153 tau_pt_weighted_dphi_strip,
154 tau_pt_weighted_dr_signal,
155 tau_pt_weighted_dr_iso,
156 tau_leadingTrackNormChi2,
159 tau_gj_angle_diff_valid,
163 tau_inside_ecal_crack,
164 leadChargedCand_etaAtEcalEntrance_minus_tau_eta,
167 std::vector<int> varsToDrop = {
168 tau_phi, tau_dxy_pca_x, tau_dxy_pca_y, tau_dxy_pca_z};
171 namespace EgammaBlockInputs {
176 tau_inside_ecal_crack,
181 pfCand_ele_pvAssociationQuality,
182 pfCand_ele_puppiWeight,
184 pfCand_ele_lostInnerHits,
185 pfCand_ele_numberOfPixelHits,
186 pfCand_ele_vertex_dx,
187 pfCand_ele_vertex_dy,
188 pfCand_ele_vertex_dz,
189 pfCand_ele_vertex_dx_tauFL,
190 pfCand_ele_vertex_dy_tauFL,
191 pfCand_ele_vertex_dz_tauFL,
192 pfCand_ele_hasTrackDetails,
197 pfCand_ele_track_chi2_ndof,
198 pfCand_ele_track_ndof,
204 ele_cc_ele_rel_energy,
205 ele_cc_gamma_rel_energy,
207 ele_rel_trackMomentumAtVtx,
208 ele_rel_trackMomentumAtCalo,
209 ele_rel_trackMomentumOut,
210 ele_rel_trackMomentumAtEleClus,
211 ele_rel_trackMomentumAtVtxWithConstraint,
214 ele_eSuperClusterOverP,
215 ele_eSeedClusterOverP,
216 ele_eSeedClusterOverPout,
217 ele_eEleClusterOverPout,
218 ele_deltaEtaSuperClusterTrackAtVtx,
219 ele_deltaEtaSeedClusterTrackAtCalo,
220 ele_deltaEtaEleClusterTrackAtCalo,
221 ele_deltaPhiEleClusterTrackAtCalo,
222 ele_deltaPhiSuperClusterTrackAtVtx,
223 ele_deltaPhiSeedClusterTrackAtCalo,
224 ele_mvaInput_earlyBrem,
225 ele_mvaInput_lateBrem,
226 ele_mvaInput_sigmaEtaEta,
227 ele_mvaInput_hadEnergy,
228 ele_mvaInput_deltaEta,
229 ele_gsfTrack_normalizedChi2,
230 ele_gsfTrack_numberOfValidHits,
233 ele_has_closestCtfTrack,
234 ele_closestCtfTrack_normalizedChi2,
235 ele_closestCtfTrack_numberOfValidHits,
240 pfCand_gamma_pvAssociationQuality,
242 pfCand_gamma_puppiWeight,
243 pfCand_gamma_puppiWeightNoLep,
244 pfCand_gamma_lostInnerHits,
245 pfCand_gamma_numberOfPixelHits,
246 pfCand_gamma_vertex_dx,
247 pfCand_gamma_vertex_dy,
248 pfCand_gamma_vertex_dz,
249 pfCand_gamma_vertex_dx_tauFL,
250 pfCand_gamma_vertex_dy_tauFL,
251 pfCand_gamma_vertex_dz_tauFL,
252 pfCand_gamma_hasTrackDetails,
254 pfCand_gamma_dxy_sig,
257 pfCand_gamma_track_chi2_ndof,
258 pfCand_gamma_track_ndof,
263 namespace MuonBlockInputs {
268 tau_inside_ecal_crack,
273 pfCand_muon_pvAssociationQuality,
275 pfCand_muon_puppiWeight,
277 pfCand_muon_lostInnerHits,
278 pfCand_muon_numberOfPixelHits,
279 pfCand_muon_vertex_dx,
280 pfCand_muon_vertex_dy,
281 pfCand_muon_vertex_dz,
282 pfCand_muon_vertex_dx_tauFL,
283 pfCand_muon_vertex_dy_tauFL,
284 pfCand_muon_vertex_dz_tauFL,
285 pfCand_muon_hasTrackDetails,
290 pfCand_muon_track_chi2_ndof,
291 pfCand_muon_track_ndof,
298 muon_normalizedChi2_valid,
300 muon_numberOfValidHits,
301 muon_segmentCompatibility,
302 muon_caloCompatibility,
303 muon_pfEcalEnergy_valid,
304 muon_rel_pfEcalEnergy,
309 muon_n_matches_CSC_1,
310 muon_n_matches_CSC_2,
311 muon_n_matches_CSC_3,
312 muon_n_matches_CSC_4,
313 muon_n_matches_RPC_1,
314 muon_n_matches_RPC_2,
315 muon_n_matches_RPC_3,
316 muon_n_matches_RPC_4,
333 namespace HadronBlockInputs {
338 tau_inside_ecal_crack,
343 pfCand_chHad_leadChargedHadrCand,
344 pfCand_chHad_pvAssociationQuality,
346 pfCand_chHad_puppiWeight,
347 pfCand_chHad_puppiWeightNoLep,
349 pfCand_chHad_lostInnerHits,
350 pfCand_chHad_numberOfPixelHits,
351 pfCand_chHad_vertex_dx,
352 pfCand_chHad_vertex_dy,
353 pfCand_chHad_vertex_dz,
354 pfCand_chHad_vertex_dx_tauFL,
355 pfCand_chHad_vertex_dy_tauFL,
356 pfCand_chHad_vertex_dz_tauFL,
357 pfCand_chHad_hasTrackDetails,
359 pfCand_chHad_dxy_sig,
362 pfCand_chHad_track_chi2_ndof,
363 pfCand_chHad_track_ndof,
364 pfCand_chHad_hcalFraction,
365 pfCand_chHad_rawCaloFraction,
370 pfCand_nHad_puppiWeight,
371 pfCand_nHad_puppiWeightNoLep,
372 pfCand_nHad_hcalFraction,
379 static tbb::concurrent_unordered_set<std::string> isFirstWarning;
385 <<
"Exception in <getTauID>: No tauID '" <<
tauID <<
"' available in pat::Tau given as function argument.";
387 if (isFirstWarning.insert(
tauID).second) {
389 <<
"' available in pat::Tau given as function argument." 390 <<
" Using default_value = " << default_value <<
" instead." << std::endl;
392 return default_value;
403 std::map<BasicDiscr, size_t> indexMap;
404 std::map<BasicDiscr, size_t> indexMapdR03;
410 return getTauID(
tau,
"chargedIsoPtSum");
416 return getTauID(
tau,
"chargedIsoPtSumdR03");
419 return (*basicTauDiscriminatordR03Collection)[tau_ref].rawValues.at(
423 return getTauID(
tau,
"footprintCorrectiondR03");
426 return getTauID(
tau,
"footprintCorrection");
432 return getTauID(
tau,
"neutralIsoPtSum");
438 return getTauID(
tau,
"neutralIsoPtSumdR03");
444 return getTauID(
tau,
"neutralIsoPtSumWeight");
448 return (*basicTauDiscriminatordR03Collection)[tau_ref].rawValues.at(
452 return getTauID(
tau,
"neutralIsoPtSumWeightdR03");
454 const float getPhotonPtSumOutsideSignalCone(
const reco::PFTau&
tau,
456 return (*basicTauDiscriminatorCollection)[tau_ref].rawValues.at(
459 const float getPhotonPtSumOutsideSignalCone(
const pat::Tau&
tau,
461 return getTauID(
tau,
"photonPtSumOutsideSignalCone");
463 const float getPhotonPtSumOutsideSignalConedR03(
const reco::PFTau&
tau,
465 return (*basicTauDiscriminatordR03Collection)[tau_ref].rawValues.at(
468 const float getPhotonPtSumOutsideSignalConedR03(
const pat::Tau&
tau,
470 return getTauID(
tau,
"photonPtSumOutsideSignalConedR03");
476 return getTauID(
tau,
"puCorrPtSum");
479 auto getdxyPCA(
const reco::PFTau&
tau,
const size_t tau_index)
const {
482 auto getdxyPCA(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.dxy_PCA(); }
483 auto getdxy(
const reco::PFTau&
tau,
const size_t tau_index)
const {
486 auto getdxy(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.dxy(); }
487 auto getdxyError(
const reco::PFTau&
tau,
const size_t tau_index)
const {
490 auto getdxyError(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.dxy_error(); }
491 auto getdxySig(
const reco::PFTau&
tau,
const size_t tau_index)
const {
494 auto getdxySig(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.dxy_Sig(); }
495 auto getip3d(
const reco::PFTau&
tau,
const size_t tau_index)
const {
498 auto getip3d(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.ip3d(); }
499 auto getip3dError(
const reco::PFTau&
tau,
const size_t tau_index)
const {
502 auto getip3dError(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.ip3d_error(); }
503 auto getip3dSig(
const reco::PFTau&
tau,
const size_t tau_index)
const {
506 auto getip3dSig(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.ip3d_Sig(); }
507 auto getHasSecondaryVertex(
const reco::PFTau&
tau,
const size_t tau_index)
const {
510 auto getHasSecondaryVertex(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.hasSecondaryVertex(); }
511 auto getFlightLength(
const reco::PFTau&
tau,
const size_t tau_index)
const {
514 auto getFlightLength(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.flightLength(); }
515 auto getFlightLengthSig(
const reco::PFTau&
tau,
const size_t tau_index)
const {
518 auto getFlightLengthSig(
const pat::Tau&
tau,
const size_t tau_index)
const {
return tau.flightLengthSig(); }
521 auto getLeadingTrackNormChi2(
const pat::Tau&
tau) {
return tau.leadingTrackNormChi2(); }
522 auto getEmFraction(
const pat::Tau&
tau) {
return tau.emFraction_MVA(); }
524 auto getEtaAtEcalEntrance(
const pat::Tau&
tau) {
return tau.etaAtEcalEntranceLeadChargedCand(); }
526 return tau.leadPFChargedHadrCand()->positionAtECALEntrance().eta();
528 auto getEcalEnergyLeadingChargedHadr(
const reco::PFTau&
tau) {
return tau.leadPFChargedHadrCand()->ecalEnergy(); }
529 auto getEcalEnergyLeadingChargedHadr(
const pat::Tau&
tau) {
return tau.ecalEnergyLeadChargedHadrCand(); }
530 auto getHcalEnergyLeadingChargedHadr(
const reco::PFTau&
tau) {
return tau.leadPFChargedHadrCand()->hcalEnergy(); }
531 auto getHcalEnergyLeadingChargedHadr(
const pat::Tau&
tau) {
return tau.hcalEnergyLeadChargedHadrCand(); }
533 template <
typename PreDiscrType>
534 bool passPrediscriminants(
const PreDiscrType prediscriminants,
535 const size_t andPrediscriminants,
537 bool passesPrediscriminants = (andPrediscriminants ? 1 : 0);
539 size_t nPrediscriminants = prediscriminants.size();
540 for (
size_t iDisc = 0; iDisc < nPrediscriminants; ++iDisc) {
542 double discResult = (*prediscriminants[iDisc].handle)[tau_ref];
543 uint8_t thisPasses = (discResult > prediscriminants[iDisc].cut) ? 1 : 0;
560 if (thisPasses ^ andPrediscriminants)
562 passesPrediscriminants = (andPrediscriminants ? 0 : 1);
566 return passesPrediscriminants;
574 return cand.bestTrack() !=
nullptr && std::isnormal(
cand.bestTrack()->dz()) && std::isnormal(
cand.dzError()) &&
578 return cand.hasTrackDetails() && std::isnormal(
cand.dz()) && std::isnormal(
cand.dzError()) &&
cand.dzError() > 0;
588 return cand.puppiWeightNoLep();
591 return cand.bestTrack() !=
nullptr 597 return cand.bestTrack() !=
nullptr 602 return cand.numberOfPixelHits();
611 return cand.rawHcalEnergy() / (
cand.rawHcalEnergy() +
cand.rawEcalEnergy());
614 float hcal_fraction = 0.;
618 hcal_fraction =
cand.hcalFraction();
621 if (
cand.pdgId() == 1 ||
cand.pdgId() == 130) {
622 hcal_fraction =
cand.hcalFraction();
623 }
else if (
cand.isIsolatedChargedHadron()) {
624 hcal_fraction =
cand.rawHcalFraction();
627 return hcal_fraction;
630 return (
cand.rawEcalEnergy() +
cand.rawHcalEnergy()) /
cand.energy();
635 template <
typename LVector1,
typename LVector2>
636 float dEta(
const LVector1& p4,
const LVector2& tau_p4) {
637 return static_cast<float>(p4.eta() - tau_p4.eta());
640 template <
typename LVector1,
typename LVector2>
641 float dPhi(
const LVector1& p4_1,
const LVector2& p4_2) {
642 return static_cast<float>(
reco::deltaPhi(p4_2.phi(), p4_1.phi()));
645 struct MuonHitMatchV2 {
646 static constexpr
size_t n_muon_stations = 4;
647 static constexpr
int first_station_id = 1;
648 static constexpr
int last_station_id = first_station_id + n_muon_stations - 1;
649 using CountArray = std::array<unsigned, n_muon_stations>;
650 using CountMap = std::map<int, CountArray>;
652 const std::vector<int>& consideredSubdets() {
658 static const std::map<int, std::string> subdet_names = {
660 if (!subdet_names.count(subdet))
661 throw cms::Exception(
"MuonHitMatch") <<
"Subdet name for subdet id " << subdet <<
" not found.";
662 return subdet_names.at(subdet);
665 size_t getStationIndex(
int station,
bool throw_exception)
const {
666 if (station < first_station_id || station > last_station_id) {
668 throw cms::Exception(
"MuonHitMatch") <<
"Station id is out of range";
671 return static_cast<size_t>(
station - 1);
675 for (
int subdet : consideredSubdets()) {
676 n_matches[subdet].fill(0);
677 n_hits[subdet].fill(0);
685 for (
const auto& segment :
muon.matches()) {
686 if (segment.segmentMatches.empty() && segment.rpcMatches.empty())
688 if (n_matches.count(segment.detector())) {
689 const size_t station_index = getStationIndex(segment.station(),
true);
690 ++n_matches.at(segment.detector()).at(station_index);
696 if (
muon.outerTrack().isNonnull()) {
697 const auto& hit_pattern =
muon.outerTrack()->hitPattern();
704 const size_t station_index = getStationIndex(hit_pattern.getMuonStation(hit_id),
false);
705 if (station_index < n_muon_stations) {
706 CountArray* muon_n_hits =
nullptr;
707 if (hit_pattern.muonDTHitFilter(hit_id))
709 else if (hit_pattern.muonCSCHitFilter(hit_id))
711 else if (hit_pattern.muonRPCHitFilter(hit_id))
715 ++muon_n_hits->at(station_index);
722 unsigned nMatches(
int subdet,
int station)
const {
723 if (!n_matches.count(subdet))
724 throw cms::Exception(
"MuonHitMatch") <<
"Subdet " << subdet <<
" not found.";
725 const size_t station_index = getStationIndex(
station,
true);
726 return n_matches.at(subdet).at(station_index);
730 if (!n_hits.count(subdet))
731 throw cms::Exception(
"MuonHitMatch") <<
"Subdet " << subdet <<
" not found.";
732 const size_t station_index = getStationIndex(
station,
true);
733 return n_hits.at(subdet).at(station_index);
736 unsigned countMuonStationsWithMatches(
int first_station,
int last_station)
const {
737 static const std::map<int, std::vector<bool>> masks = {
742 const size_t first_station_index = getStationIndex(first_station,
true);
743 const size_t last_station_index = getStationIndex(last_station,
true);
745 for (
size_t n = first_station_index;
n <= last_station_index; ++
n) {
746 for (
const auto&
match : n_matches) {
747 if (!masks.at(
match.first).at(
n) &&
match.second.at(
n) > 0)
754 unsigned countMuonStationsWithHits(
int first_station,
int last_station)
const {
755 static const std::map<int, std::vector<bool>> masks = {
761 const size_t first_station_index = getStationIndex(first_station,
true);
762 const size_t last_station_index = getStationIndex(last_station,
true);
764 for (
size_t n = first_station_index;
n <= last_station_index; ++
n) {
765 for (
const auto&
hit : n_hits) {
766 if (!masks.at(
hit.first).at(
n) &&
hit.second.at(
n) > 0)
774 CountMap n_matches, n_hits;
780 PfCand_chargedHadron,
781 PfCand_neutralHadron,
788 template <
typename Object>
801 static const std::map<int, CellObjectType> obj_types = {{11, CellObjectType::PfCand_electron},
802 {13, CellObjectType::PfCand_muon},
803 {22, CellObjectType::PfCand_gamma},
804 {130, CellObjectType::PfCand_neutralHadron},
805 {211, CellObjectType::PfCand_chargedHadron}};
808 if (iter == obj_types.end())
809 return CellObjectType::Other;
813 using Cell = std::map<CellObjectType, size_t>;
826 using Map = std::map<CellIndex, Cell>;
827 using const_iterator = Map::const_iterator;
829 CellGrid(
unsigned n_cells_eta,
830 unsigned n_cells_phi,
831 double cell_size_eta,
832 double cell_size_phi,
834 : nCellsEta(n_cells_eta),
835 nCellsPhi(n_cells_phi),
836 nTotal(nCellsEta * nCellsPhi),
837 cellSizeEta(cell_size_eta),
838 cellSizePhi(cell_size_phi),
840 if (nCellsEta % 2 != 1 || nCellsEta < 1)
841 throw cms::Exception(
"DeepTauId") <<
"Invalid number of eta cells.";
842 if (nCellsPhi % 2 != 1 || nCellsPhi < 1)
843 throw cms::Exception(
"DeepTauId") <<
"Invalid number of phi cells.";
844 if (cellSizeEta <= 0 || cellSizePhi <= 0)
848 int maxEtaIndex()
const {
return static_cast<int>((nCellsEta - 1) / 2); }
849 int maxPhiIndex()
const {
return static_cast<int>((nCellsPhi - 1) / 2); }
850 double maxDeltaEta()
const {
return cellSizeEta * (0.5 + maxEtaIndex()); }
851 double maxDeltaPhi()
const {
return cellSizePhi * (0.5 + maxPhiIndex()); }
852 int getEtaTensorIndex(
const CellIndex& cellIndex)
const {
return cellIndex.eta + maxEtaIndex(); }
853 int getPhiTensorIndex(
const CellIndex& cellIndex)
const {
return cellIndex.phi + maxPhiIndex(); }
855 bool tryGetCellIndex(
double deltaEta,
double deltaPhi, CellIndex& cellIndex)
const {
856 const auto getCellIndex = [
this](
double x,
double maxX,
double size,
int&
index) {
861 if (disable_CellIndex_workaround_) {
864 absIndex = std::floor(absX /
size + 0.5);
869 index =
static_cast<int>(std::copysign(absIndex,
x));
877 size_t num_valid_cells()
const {
return cells.size(); }
879 const Cell& at(
const CellIndex& cellIndex)
const {
return cells.at(cellIndex); }
880 size_t count(
const CellIndex& cellIndex)
const {
return cells.count(cellIndex); }
881 const_iterator
find(
const CellIndex& cellIndex)
const {
return cells.find(cellIndex); }
882 const_iterator
begin()
const {
return cells.begin(); }
883 const_iterator
end()
const {
return cells.end(); }
886 const unsigned nCellsEta, nCellsPhi, nTotal;
887 const double cellSizeEta, cellSizePhi;
890 std::map<CellIndex, Cell>
cells;
891 const bool disable_CellIndex_workaround_;
915 const tensorflow::Tensor& pred,
918 std::vector<reco::SingleTauDiscriminatorContainer> outputbuffer(
taus->size());
920 for (
size_t tau_index = 0; tau_index <
taus->size(); ++tau_index) {
922 for (
size_t num_elem :
num_)
923 x += pred.matrix<
float>()(tau_index, num_elem);
924 if (
x != 0 && !
den_.empty()) {
926 for (
size_t den_elem :
den_)
927 den_val += pred.matrix<
float>()(tau_index, den_elem);
930 outputbuffer[tau_index].rawValues.push_back(
x);
934 outputbuffer[tau_index].workingPoints.push_back(pass);
938 std::unique_ptr<TauDiscriminator>
output = std::make_unique<TauDiscriminator>();
940 filler.insert(
taus, outputbuffer.begin(), outputbuffer.end());
951 static constexpr
size_t e_index = 0, mu_index = 1, tau_index = 2, jet_index = 3;
953 {
"VSe",
Output({tau_index}, {e_index, tau_index})},
954 {
"VSmu",
Output({tau_index}, {mu_index, tau_index})},
955 {
"VSjet",
Output({tau_index}, {jet_index, tau_index})},
965 std::vector<BasicDiscriminator> requiredDiscr) {
966 std::map<std::string, size_t> discrIndexMapStr;
967 auto const aHandle =
event.getHandle(discriminatorContainerToken);
968 auto const aProv = aHandle.provenance();
969 if (aProv ==
nullptr)
970 aHandle.whyFailed()->raise();
971 const auto& psetsFromProvenance =
edm::parameterSet(aProv->stable(),
event.processHistory());
972 auto const idlist = psetsFromProvenance.
getParameter<std::vector<edm::ParameterSet>>(
"IDdefinitions");
973 for (
size_t j = 0;
j < idlist.size(); ++
j) {
975 if (discrIndexMapStr.count(idname)) {
977 <<
"basic discriminator " << idname <<
" appears more than once in the input.";
979 discrIndexMapStr[idname] =
j;
983 std::map<BasicDiscriminator, size_t> discrIndexMap;
984 for (
size_t i = 0;
i < requiredDiscr.size();
i++) {
987 <<
" was not provided in the config file.";
991 return discrIndexMap;
1002 desc.add<std::vector<std::string>>(
"graph_file",
1003 {
"RecoTauTag/TrainingFiles/data/DeepTauId/deepTau_2017v2p6_e6.pb"});
1004 desc.add<
bool>(
"mem_mapped",
false);
1005 desc.add<
unsigned>(
"year", 2017);
1006 desc.add<
unsigned>(
"version", 2);
1007 desc.add<
unsigned>(
"sub_version", 1);
1008 desc.add<
int>(
"debug_level", 0);
1009 desc.add<
bool>(
"disable_dxy_pca",
false);
1010 desc.add<
bool>(
"disable_hcalFraction_workaround",
false);
1011 desc.add<
bool>(
"disable_CellIndex_workaround",
false);
1012 desc.add<
bool>(
"save_inputs",
false);
1013 desc.add<
bool>(
"is_online",
false);
1015 desc.add<std::vector<std::string>>(
"VSeWP", {
"-1."});
1016 desc.add<std::vector<std::string>>(
"VSmuWP", {
"-1."});
1017 desc.add<std::vector<std::string>>(
"VSjetWP", {
"-1."});
1025 pset_Prediscriminants.
add<
std::string>(
"BooleanOperator",
"and");
1028 psd1.
add<
double>(
"cut");
1035 descriptions.
add(
"DeepTau",
desc);
1050 cfg.getUntrackedParameter<
edm::
InputTag>(
"basicTauDiscriminators"))),
1052 cfg.getUntrackedParameter<
edm::
InputTag>(
"basicTauDiscriminatorsdR03"))),
1055 cfg.getParameter<
edm::
InputTag>(
"pfTauTransverseImpactParameters"))),
1056 year_(
cfg.getParameter<unsigned>(
"year")),
1066 for (
const auto& output_desc :
outputs_) {
1067 produces<TauDiscriminator>(output_desc.first);
1068 const auto& cut_list =
cfg.getParameter<std::vector<std::string>>(output_desc.first +
"WP");
1070 workingPoints_[output_desc.first].push_back(std::make_unique<Cutter>(cut_str));
1081 transform(pdBoolOperator.begin(), pdBoolOperator.end(), pdBoolOperator.begin(), ::tolower);
1083 if (pdBoolOperator ==
"and") {
1085 }
else if (pdBoolOperator ==
"or") {
1089 <<
"PrediscriminantBooleanOperator defined incorrectly, options are: AND,OR";
1093 std::vector<std::string> prediscriminantsNames =
1096 for (
auto const& iDisc : prediscriminantsNames) {
1104 thisDiscriminator.
cut =
cut;
1105 thisDiscriminator.
disc_token = consumes<reco::PFTauDiscriminator>(
label);
1110 thisDiscriminator.
cut =
cut;
1111 thisDiscriminator.
disc_token = consumes<pat::PATTauDiscriminator>(
label);
1123 tensorflow::DT_FLOAT, tensorflow::TensorShape{1, TauBlockInputs::NumberOfInputs});
1126 std::sort(TauBlockInputs::varsToDrop.begin(), TauBlockInputs::varsToDrop.
end());
1127 for (
auto v : TauBlockInputs::varsToDrop) {
1129 for (std::size_t
i =
v + 1;
i < TauBlockInputs::NumberOfInputs; ++
i)
1133 tensorflow::DT_FLOAT,
1134 tensorflow::TensorShape{1,
1135 static_cast<int>(TauBlockInputs::NumberOfInputs) -
1136 static_cast<int>(TauBlockInputs::varsToDrop.
size())});
1137 if (
year_ == 2026) {
1145 std::map<std::vector<bool>, std::vector<sc::FeatureT>> GridFeatureTypes_map = {
1146 {{
false}, {sc::FeatureT::TauFlat, sc::FeatureT::GridGlobal}},
1148 {sc::FeatureT::PfCand_electron,
1149 sc::FeatureT::PfCand_muon,
1150 sc::FeatureT::PfCand_chHad,
1151 sc::FeatureT::PfCand_nHad,
1152 sc::FeatureT::PfCand_gamma,
1157 for (
const auto&
p : GridFeatureTypes_map) {
1158 for (
auto is_inner :
p.first) {
1159 for (
auto featureType :
p.second) {
1160 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(featureType, is_inner));
1161 if (!(sp.mean_.size() == sp.std_.size() && sp.mean_.size() == sp.lim_min_.size() &&
1162 sp.mean_.size() == sp.lim_max_.size()))
1163 throw cms::Exception(
"DeepTauId") <<
"sizes of scaling parameter vectors do not match between each other";
1168 for (
size_t n = 0;
n < 2; ++
n) {
1169 const bool is_inner =
n == 0;
1170 const auto n_cells = is_inner ? number_of_inner_cell : number_of_outer_cell;
1171 eGammaTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1172 tensorflow::DT_FLOAT, tensorflow::TensorShape{1, 1, 1, EgammaBlockInputs::NumberOfInputs});
1173 muonTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1174 tensorflow::DT_FLOAT, tensorflow::TensorShape{1, 1, 1, MuonBlockInputs::NumberOfInputs});
1176 tensorflow::DT_FLOAT, tensorflow::TensorShape{1, 1, 1, HadronBlockInputs::NumberOfInputs});
1177 convTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1178 tensorflow::DT_FLOAT, tensorflow::TensorShape{1, n_cells, n_cells, number_of_conv_features});
1180 tensorflow::DT_FLOAT, tensorflow::TensorShape{1, 1, 1, number_of_conv_features});
1198 if (
taus->empty()) {
1207 size_t nPrediscriminants =
1209 for (
size_t iDisc = 0; iDisc < nPrediscriminants; ++iDisc) {
1221 if (tauProductID != discKeyId) {
1223 <<
"The tau collection has product ID: " << tauProductID
1224 <<
" but the pre-discriminator is keyed with product ID: " << discKeyId << std::endl;
1233 const auto graph_name_vector =
cfg.getParameter<std::vector<std::string>>(
"graph_file");
1234 std::map<std::string, std::string> graph_names;
1235 for (
const auto&
entry : graph_name_vector) {
1236 const size_t sep_pos =
entry.find(
':');
1238 if (sep_pos != std::string::npos) {
1239 entry_name =
entry.substr(0, sep_pos);
1246 if (graph_names.count(entry_name))
1247 throw cms::Exception(
"DeepTauCache") <<
"Duplicated graph entries";
1251 return std::make_unique<deep_tau::DeepTauCache>(graph_names,
mem_mapped);
1256 template <
typename ConsumeType>
1274 for (
const auto& output_desc :
outputs_) {
1284 template <
typename T>
1286 return std::isnormal(
value) ?
static_cast<float>(
value) : 0.
f;
1289 template <
typename T>
1295 transformed_value = transformed_value * 2 - 1;
1296 return transformed_value;
1299 template <
typename T>
1302 const float norm_value = (fixed_value -
mean) / sigma;
1303 return std::clamp(norm_value, -n_sigmas_max, n_sigmas_max);
1309 float& cc_ele_energy,
1310 float& cc_gamma_energy,
1312 cc_ele_energy = cc_gamma_energy = 0;
1315 if (superCluster.isNonnull() && superCluster.isAvailable() && superCluster->clusters().isNonnull() &&
1316 superCluster->clusters().isAvailable()) {
1317 for (
auto iter = superCluster->clustersBegin(); iter != superCluster->clustersEnd(); ++iter) {
1318 const float energy =
static_cast<float>((*iter)->energy());
1319 if (iter == superCluster->clustersBegin())
1322 cc_gamma_energy +=
energy;
1334 const CellGrid*
grid =
nullptr)
const {
1336 std::cout <<
"<checkInputs>: block_name = " << block_name << std::endl;
1337 if (block_name ==
"input_tau") {
1338 for (
int input_index = 0; input_index < n_inputs; ++input_index) {
1339 float input =
inputs.matrix<
float>()(0, input_index);
1342 <<
"in the " << block_name
1343 <<
", input is not finite, i.e. infinite or NaN, for input_index = " << input_index;
1353 if (block_name.find(
"input_inner") != std::string::npos) {
1356 }
else if (block_name.find(
"input_outer") != std::string::npos) {
1361 int eta_phi_index = 0;
1362 for (
int eta = -n_eta;
eta <= n_eta; ++
eta) {
1363 for (
int phi = -n_phi;
phi <= n_phi; ++
phi) {
1364 const CellIndex cell_index{
eta,
phi};
1365 const auto cell_iter =
grid->find(cell_index);
1366 if (cell_iter !=
grid->end()) {
1367 for (
int input_index = 0; input_index < n_inputs; ++input_index) {
1371 <<
"in the " << block_name <<
", input is not finite, i.e. infinite or NaN, for eta = " <<
eta 1372 <<
", phi = " <<
phi <<
", input_index = " << input_index;
1375 std::cout << block_name <<
"[eta = " <<
eta <<
"][phi = " <<
phi <<
"][var = " << input_index
1376 <<
"] = " << std::setprecision(5) <<
std::fixed <<
input << std::endl;
1390 const CellGrid*
grid =
nullptr) {
1392 std::cout <<
"<saveInputs>: block_name = " << block_name << std::endl;
1395 (*json_file_) <<
", ";
1396 (*json_file_) <<
"\"" << block_name <<
"\": [";
1397 if (block_name ==
"input_tau") {
1398 for (
int input_index = 0; input_index < n_inputs; ++input_index) {
1399 float input =
inputs.matrix<
float>()(0, input_index);
1400 if (input_index != 0)
1401 (*json_file_) <<
", ";
1402 (*json_file_) <<
input;
1407 if (block_name.find(
"input_inner") != std::string::npos) {
1410 }
else if (block_name.find(
"input_outer") != std::string::npos) {
1415 int eta_phi_index = 0;
1416 for (
int eta = -n_eta;
eta <= n_eta; ++
eta) {
1418 (*json_file_) <<
", ";
1419 (*json_file_) <<
"[";
1420 for (
int phi = -n_phi;
phi <= n_phi; ++
phi) {
1422 (*json_file_) <<
", ";
1423 (*json_file_) <<
"[";
1424 const CellIndex cell_index{
eta,
phi};
1425 const auto cell_iter =
grid->find(cell_index);
1426 for (
int input_index = 0; input_index < n_inputs; ++input_index) {
1428 if (cell_iter !=
grid->end()) {
1431 if (input_index != 0)
1432 (*json_file_) <<
", ";
1433 (*json_file_) <<
input;
1435 if (cell_iter !=
grid->end()) {
1438 (*json_file_) <<
"]";
1440 (*json_file_) <<
"]";
1443 (*json_file_) <<
"]";
1450 const std::vector<pat::Electron> electron_collection_default;
1451 const std::vector<pat::Muon> muon_collection_default;
1455 pfTauTransverseImpactParameters_default;
1457 const std::vector<pat::Electron>* electron_collection;
1458 const std::vector<pat::Muon>* muon_collection;
1471 electron_collection = &electron_collection_default;
1472 muon_collection = &muon_collection_default;
1502 auto const& eventnr =
event.id().event();
1506 for (
size_t tau_index = 0; tau_index <
taus->size(); ++tau_index) {
1509 std::vector<tensorflow::Tensor> pred_vector;
1511 bool passesPrediscriminants;
1513 passesPrediscriminants = tauIDs.passPrediscriminants<std::vector<TauDiscInfo<reco::PFTauDiscriminator>>>(
1516 passesPrediscriminants = tauIDs.passPrediscriminants<std::vector<TauDiscInfo<pat::PATTauDiscriminator>>>(
1520 if (passesPrediscriminants) {
1523 getPredictionsV2<reco::PFCandidate, reco::PFTau>(
taus->at(tau_index),
1526 electron_collection,
1535 getPredictionsV2<pat::PackedCandidate, pat::Tau>(
taus->at(tau_index),
1538 electron_collection,
1551 const float pred = pred_vector[0].flat<
float>()(
k);
1552 if (!(pred >= 0 && pred <= 1))
1554 <<
"invalid prediction = " << pred <<
" for tau_index = " << tau_index <<
", pred_index = " <<
k;
1555 predictions.matrix<
float>()(tau_index,
k) = pred;
1564 predictions.matrix<
float>()(tau_index,
k) = (
k == 2) ? -1.
f : 2.
f;
1571 template <
typename Cand
idateCastType,
typename TauCastType>
1573 const size_t tau_index,
1575 const std::vector<pat::Electron>*
electrons,
1576 const std::vector<pat::Muon>*
muons,
1581 std::vector<tensorflow::Tensor>& pred_vector,
1582 TauFunc tau_funcs) {
1585 std::cout <<
"<DeepTauId::getPredictionsV2 (moduleLabel = " << moduleDescription().moduleLabel()
1586 <<
")>:" << std::endl;
1587 std::cout <<
" tau: pT = " <<
tau.pt() <<
", eta = " <<
tau.eta() <<
", phi = " <<
tau.phi()
1588 <<
", eventnr = " << eventnr << std::endl;
1593 fillGrids(dynamic_cast<const TauCastType&>(
tau), *
muons, inner_grid, outer_grid);
1594 fillGrids(dynamic_cast<const TauCastType&>(
tau), pfCands, inner_grid, outer_grid);
1596 createTauBlockInputs<CandidateCastType>(
1597 dynamic_cast<const TauCastType&
>(
tau), tau_index, tau_ref,
pv,
rho, tau_funcs);
1599 createConvFeatures<CandidateCastType>(
dynamic_cast<const TauCastType&
>(
tau),
1613 createConvFeatures<CandidateCastType>(
dynamic_cast<const TauCastType&
>(
tau),
1630 json_file_ =
new std::ofstream(json_file_name.data());
1632 (*json_file_) <<
"{";
1635 *
eGammaTensor_[
true],
"input_inner_egamma", dnn_inputs_v2::EgammaBlockInputs::NumberOfInputs, &inner_grid);
1636 saveInputs(*
muonTensor_[
true],
"input_inner_muon", dnn_inputs_v2::MuonBlockInputs::NumberOfInputs, &inner_grid);
1638 *
hadronsTensor_[
true],
"input_inner_hadrons", dnn_inputs_v2::HadronBlockInputs::NumberOfInputs, &inner_grid);
1640 *
eGammaTensor_[
false],
"input_outer_egamma", dnn_inputs_v2::EgammaBlockInputs::NumberOfInputs, &outer_grid);
1641 saveInputs(*
muonTensor_[
false],
"input_outer_muon", dnn_inputs_v2::MuonBlockInputs::NumberOfInputs, &outer_grid);
1643 *
hadronsTensor_[
false],
"input_outer_hadrons", dnn_inputs_v2::HadronBlockInputs::NumberOfInputs, &outer_grid);
1644 (*json_file_) <<
"}";
1653 {
"main_output/Softmax"},
1677 template <
typename Collection,
typename TauCastType>
1678 void fillGrids(
const TauCastType&
tau,
const Collection&
objects, CellGrid& inner_grid, CellGrid& outer_grid) {
1679 static constexpr
double outer_dR2 = 0.25;
1681 const double inner_dR2 =
std::pow(inner_radius, 2);
1683 const auto addObject = [&](
size_t n,
double deta,
double dphi, CellGrid&
grid) {
1686 if (obj_type == CellObjectType::Other)
1688 CellIndex cell_index;
1689 if (
grid.tryGetCellIndex(deta, dphi, cell_index)) {
1691 auto iter = cell.find(obj_type);
1692 if (iter != cell.end()) {
1693 const auto& prev_obj =
objects.at(iter->second);
1694 if (
obj.polarP4().pt() > prev_obj.polarP4().pt())
1702 for (
size_t n = 0;
n <
objects.size(); ++
n) {
1704 const double deta =
obj.polarP4().eta() -
tau.polarP4().eta();
1707 if (dR2 < inner_dR2)
1708 addObject(
n, deta, dphi, inner_grid);
1709 if (dR2 < outer_dR2)
1710 addObject(
n, deta, dphi, outer_grid);
1715 std::vector<tensorflow::Tensor> pred_vector;
1723 {
"inner_all_dropout_4/Identity"},
1732 {
"outer_all_dropout_4/Identity"},
1735 return pred_vector.at(0);
1738 template <
typename Cand
idateCastType,
typename TauCastType>
1740 const size_t tau_index,
1744 const std::vector<pat::Electron>*
electrons,
1745 const std::vector<pat::Muon>*
muons,
1747 const CellGrid&
grid,
1751 std::cout <<
"<DeepTauId::createConvFeatures (is_inner = " << is_inner <<
")>:" << std::endl;
1753 tensorflow::Tensor& convTensor = *
convTensor_.at(is_inner);
1754 eGammaTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1755 tensorflow::DT_FLOAT,
1756 tensorflow::TensorShape{
1757 (
long long int)
grid.num_valid_cells(), 1, 1, dnn_inputs_v2::EgammaBlockInputs::NumberOfInputs});
1758 muonTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1759 tensorflow::DT_FLOAT,
1760 tensorflow::TensorShape{
1761 (
long long int)
grid.num_valid_cells(), 1, 1, dnn_inputs_v2::MuonBlockInputs::NumberOfInputs});
1763 tensorflow::DT_FLOAT,
1764 tensorflow::TensorShape{
1765 (
long long int)
grid.num_valid_cells(), 1, 1, dnn_inputs_v2::HadronBlockInputs::NumberOfInputs});
1775 std::cout <<
"processing ( eta = " <<
eta <<
", phi = " <<
phi <<
" )" << std::endl;
1777 const CellIndex cell_index{
eta,
phi};
1778 const auto cell_iter =
grid.find(cell_index);
1779 if (cell_iter !=
grid.end()) {
1781 std::cout <<
" creating inputs for ( eta = " <<
eta <<
", phi = " <<
phi <<
" ): idx = " <<
idx 1784 const Cell& cell = cell_iter->second;
1785 createEgammaBlockInputs<CandidateCastType>(
1786 idx,
tau, tau_index, tau_ref,
pv,
rho,
electrons, pfCands, cell, tau_funcs, is_inner);
1787 createMuonBlockInputs<CandidateCastType>(
1788 idx,
tau, tau_index, tau_ref,
pv,
rho,
muons, pfCands, cell, tau_funcs, is_inner);
1789 createHadronsBlockInputs<CandidateCastType>(
1790 idx,
tau, tau_index, tau_ref,
pv,
rho, pfCands, cell, tau_funcs, is_inner);
1794 std::cout <<
" skipping creation of inputs, because ( eta = " <<
eta <<
", phi = " <<
phi 1795 <<
" ) is not in the grid !!" << std::endl;
1805 const CellIndex cell_index{
eta,
phi};
1806 const int eta_index =
grid.getEtaTensorIndex(cell_index);
1807 const int phi_index =
grid.getPhiTensorIndex(cell_index);
1809 const auto cell_iter =
grid.find(cell_index);
1810 if (cell_iter !=
grid.end()) {
1821 const tensorflow::Tensor&
features,
1825 for (
int n = 0;
n < dnn_inputs_v2::number_of_conv_features; ++
n) {
1826 convTensor.tensor<
float, 4>()(0, eta_index, phi_index,
n) =
features.tensor<
float, 4>()(batch_idx, 0, 0,
n);
1830 template <
typename Cand
idateCastType,
typename TauCastType>
1832 const size_t& tau_index,
1836 TauFunc tau_funcs) {
1843 inputs.flat<
float>().setZero();
1845 const auto&
get = [&](
int var_index) ->
float& {
1849 auto leadChargedHadrCand =
dynamic_cast<const CandidateCastType*
>(
tau.leadChargedHadrCand().get());
1860 get(dnn::tau_n_charged_prongs) = sp.scale(
tau.decayMode() / 5 + 1,
tauInputs_indices_[dnn::tau_n_charged_prongs]);
1861 get(dnn::tau_n_neutral_prongs) = sp.scale(
tau.decayMode() % 5,
tauInputs_indices_[dnn::tau_n_neutral_prongs]);
1862 get(dnn::chargedIsoPtSum) =
1864 get(dnn::chargedIsoPtSumdR03_over_dR05) =
1865 sp.scale(tau_funcs.getChargedIsoPtSumdR03(
tau, tau_ref) / tau_funcs.getChargedIsoPtSum(
tau, tau_ref),
1868 get(dnn::footprintCorrection) =
1869 sp.scale(tau_funcs.getFootprintCorrectiondR03(
tau, tau_ref),
tauInputs_indices_[dnn::footprintCorrection]);
1871 get(dnn::footprintCorrection) =
1872 sp.scale(tau_funcs.getFootprintCorrection(
tau, tau_ref),
tauInputs_indices_[dnn::footprintCorrection]);
1874 get(dnn::neutralIsoPtSum) =
1876 get(dnn::neutralIsoPtSumWeight_over_neutralIsoPtSum) =
1877 sp.scale(tau_funcs.getNeutralIsoPtSumWeight(
tau, tau_ref) / tau_funcs.getNeutralIsoPtSum(
tau, tau_ref),
1879 get(dnn::neutralIsoPtSumWeightdR03_over_neutralIsoPtSum) =
1880 sp.scale(tau_funcs.getNeutralIsoPtSumdR03Weight(
tau, tau_ref) / tau_funcs.getNeutralIsoPtSum(
tau, tau_ref),
1882 get(dnn::neutralIsoPtSumdR03_over_dR05) =
1883 sp.scale(tau_funcs.getNeutralIsoPtSumdR03(
tau, tau_ref) / tau_funcs.getNeutralIsoPtSum(
tau, tau_ref),
1885 get(dnn::photonPtSumOutsideSignalCone) = sp.scale(tau_funcs.getPhotonPtSumOutsideSignalCone(
tau, tau_ref),
1887 get(dnn::puCorrPtSum) = sp.scale(tau_funcs.getPuCorrPtSum(
tau, tau_ref),
tauInputs_indices_[dnn::puCorrPtSum]);
1894 auto const pca = tau_funcs.getdxyPCA(
tau, tau_index);
1899 get(dnn::tau_dxy_pca_x) = 0;
1900 get(dnn::tau_dxy_pca_y) = 0;
1901 get(dnn::tau_dxy_pca_z) = 0;
1905 const bool tau_dxy_valid =
1906 isAbove(tau_funcs.getdxy(
tau, tau_index), -10) &&
isAbove(tau_funcs.getdxyError(
tau, tau_index), 0);
1907 if (tau_dxy_valid) {
1908 get(dnn::tau_dxy_valid) = sp.scale(tau_dxy_valid,
tauInputs_indices_[dnn::tau_dxy_valid]);
1910 get(dnn::tau_dxy_sig) =
1911 sp.scale(
std::abs(tau_funcs.getdxy(
tau, tau_index)) / tau_funcs.getdxyError(
tau, tau_index),
1914 const bool tau_ip3d_valid =
1915 isAbove(tau_funcs.getip3d(
tau, tau_index), -10) &&
isAbove(tau_funcs.getip3dError(
tau, tau_index), 0);
1916 if (tau_ip3d_valid) {
1917 get(dnn::tau_ip3d_valid) = sp.scale(tau_ip3d_valid,
tauInputs_indices_[dnn::tau_ip3d_valid]);
1918 get(dnn::tau_ip3d) = sp.scale(tau_funcs.getip3d(
tau, tau_index),
tauInputs_indices_[dnn::tau_ip3d]);
1919 get(dnn::tau_ip3d_sig) =
1920 sp.scale(
std::abs(tau_funcs.getip3d(
tau, tau_index)) / tau_funcs.getip3dError(
tau, tau_index),
1923 if (leadChargedHadrCand) {
1924 const bool hasTrackDetails = candFunc::getHasTrackDetails(*leadChargedHadrCand);
1925 const float tau_dz = (
is_online_ && !hasTrackDetails) ? 0 : candFunc::getTauDz(*leadChargedHadrCand);
1927 get(dnn::tau_dz_sig_valid) =
1928 sp.scale(candFunc::getTauDZSigValid(*leadChargedHadrCand),
tauInputs_indices_[dnn::tau_dz_sig_valid]);
1929 const double dzError = hasTrackDetails ? leadChargedHadrCand->dzError() : -999.;
1932 get(dnn::tau_flightLength_x) =
1933 sp.scale(tau_funcs.getFlightLength(
tau, tau_index).x(),
tauInputs_indices_[dnn::tau_flightLength_x]);
1934 get(dnn::tau_flightLength_y) =
1935 sp.scale(tau_funcs.getFlightLength(
tau, tau_index).y(),
tauInputs_indices_[dnn::tau_flightLength_y]);
1936 get(dnn::tau_flightLength_z) =
1937 sp.scale(tau_funcs.getFlightLength(
tau, tau_index).z(),
tauInputs_indices_[dnn::tau_flightLength_z]);
1939 get(dnn::tau_flightLength_sig) = 0.55756444;
1941 get(dnn::tau_flightLength_sig) =
1942 sp.scale(tau_funcs.getFlightLengthSig(
tau, tau_index),
tauInputs_indices_[dnn::tau_flightLength_sig]);
1951 get(dnn::tau_pt_weighted_dr_iso) =
1953 get(dnn::tau_leadingTrackNormChi2) =
1954 sp.scale(tau_funcs.getLeadingTrackNormChi2(
tau),
tauInputs_indices_[dnn::tau_leadingTrackNormChi2]);
1956 const bool tau_e_ratio_valid = std::isnormal(
eratio) &&
eratio > 0.f;
1957 get(dnn::tau_e_ratio_valid) = sp.scale(tau_e_ratio_valid,
tauInputs_indices_[dnn::tau_e_ratio_valid]);
1960 const bool tau_gj_angle_diff_valid = (std::isnormal(gj_angle_diff) || gj_angle_diff == 0) && gj_angle_diff >= 0;
1961 get(dnn::tau_gj_angle_diff_valid) =
1963 get(dnn::tau_gj_angle_diff) =
1964 tau_gj_angle_diff_valid ? sp.scale(gj_angle_diff,
tauInputs_indices_[dnn::tau_gj_angle_diff]) : 0;
1966 get(dnn::tau_emFraction) = sp.scale(tau_funcs.getEmFraction(
tau),
tauInputs_indices_[dnn::tau_emFraction]);
1968 get(dnn::tau_inside_ecal_crack) =
1970 get(dnn::leadChargedCand_etaAtEcalEntrance_minus_tau_eta) =
1971 sp.scale(tau_funcs.getEtaAtEcalEntrance(
tau) -
tau.p4().eta(),
1975 template <
typename Cand
idateCastType,
typename TauCastType>
1977 const TauCastType&
tau,
1978 const size_t tau_index,
1982 const std::vector<pat::Electron>*
electrons,
1984 const Cell& cell_map,
1995 int PFe_index_offset =
scalingParamsMap_->at(std::make_pair(ft_global,
false)).mean_.size();
1996 int e_index_offset = PFe_index_offset +
scalingParamsMap_->at(std::make_pair(ft_PFe,
false)).mean_.size();
1997 int PFg_index_offset = e_index_offset +
scalingParamsMap_->at(std::make_pair(ft_e,
false)).mean_.size();
2000 int fill_index_offset_e = 0;
2001 int fill_index_offset_PFg = 0;
2003 fill_index_offset_e =
2005 fill_index_offset_PFg =
2011 const auto&
get = [&](
int var_index) ->
float& {
return inputs.tensor<
float, 4>()(
idx, 0, 0, var_index); };
2013 const bool valid_index_pf_ele = cell_map.count(CellObjectType::PfCand_electron);
2014 const bool valid_index_pf_gamma = cell_map.count(CellObjectType::PfCand_gamma);
2017 if (!cell_map.empty()) {
2018 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_global,
false));
2022 get(dnn::tau_inside_ecal_crack) = sp.scale(
isInEcalCrack(
tau.polarP4().eta()), dnn::tau_inside_ecal_crack);
2024 if (valid_index_pf_ele) {
2025 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_PFe, is_inner));
2026 size_t index_pf_ele = cell_map.at(CellObjectType::PfCand_electron);
2027 const auto& ele_cand =
dynamic_cast<const CandidateCastType&
>(pfCands.
at(index_pf_ele));
2029 get(dnn::pfCand_ele_valid) = sp.scale(valid_index_pf_ele, dnn::pfCand_ele_valid - PFe_index_offset);
2030 get(dnn::pfCand_ele_rel_pt) =
2031 sp.scale(ele_cand.polarP4().pt() /
tau.polarP4().pt(), dnn::pfCand_ele_rel_pt - PFe_index_offset);
2032 get(dnn::pfCand_ele_deta) =
2033 sp.scale(ele_cand.polarP4().eta() -
tau.polarP4().eta(), dnn::pfCand_ele_deta - PFe_index_offset);
2034 get(dnn::pfCand_ele_dphi) =
2035 sp.scale(
dPhi(
tau.polarP4(), ele_cand.polarP4()), dnn::pfCand_ele_dphi - PFe_index_offset);
2036 get(dnn::pfCand_ele_pvAssociationQuality) = sp.scale<
int>(
2037 candFunc::getPvAssocationQuality(ele_cand), dnn::pfCand_ele_pvAssociationQuality - PFe_index_offset);
2038 get(dnn::pfCand_ele_puppiWeight) = is_inner ? sp.scale(candFunc::getPuppiWeight(ele_cand, 0.9906834f),
2039 dnn::pfCand_ele_puppiWeight - PFe_index_offset)
2040 : sp.scale(candFunc::getPuppiWeight(ele_cand, 0.9669586
f),
2041 dnn::pfCand_ele_puppiWeight - PFe_index_offset);
2042 get(dnn::pfCand_ele_charge) = sp.scale(ele_cand.charge(), dnn::pfCand_ele_charge - PFe_index_offset);
2043 get(dnn::pfCand_ele_lostInnerHits) =
2044 sp.scale<
int>(candFunc::getLostInnerHits(ele_cand, 0), dnn::pfCand_ele_lostInnerHits - PFe_index_offset);
2045 get(dnn::pfCand_ele_numberOfPixelHits) =
2046 sp.scale(candFunc::getNumberOfPixelHits(ele_cand, 0), dnn::pfCand_ele_numberOfPixelHits - PFe_index_offset);
2047 get(dnn::pfCand_ele_vertex_dx) =
2048 sp.scale(ele_cand.vertex().x() -
pv.position().x(), dnn::pfCand_ele_vertex_dx - PFe_index_offset);
2049 get(dnn::pfCand_ele_vertex_dy) =
2050 sp.scale(ele_cand.vertex().y() -
pv.position().y(), dnn::pfCand_ele_vertex_dy - PFe_index_offset);
2051 get(dnn::pfCand_ele_vertex_dz) =
2052 sp.scale(ele_cand.vertex().z() -
pv.position().z(), dnn::pfCand_ele_vertex_dz - PFe_index_offset);
2053 get(dnn::pfCand_ele_vertex_dx_tauFL) =
2054 sp.scale(ele_cand.vertex().x() -
pv.position().x() - tau_funcs.getFlightLength(
tau, tau_index).x(),
2055 dnn::pfCand_ele_vertex_dx_tauFL - PFe_index_offset);
2056 get(dnn::pfCand_ele_vertex_dy_tauFL) =
2057 sp.scale(ele_cand.vertex().y() -
pv.position().y() - tau_funcs.getFlightLength(
tau, tau_index).y(),
2058 dnn::pfCand_ele_vertex_dy_tauFL - PFe_index_offset);
2059 get(dnn::pfCand_ele_vertex_dz_tauFL) =
2060 sp.scale(ele_cand.vertex().z() -
pv.position().z() - tau_funcs.getFlightLength(
tau, tau_index).z(),
2061 dnn::pfCand_ele_vertex_dz_tauFL - PFe_index_offset);
2063 const bool hasTrackDetails = candFunc::getHasTrackDetails(ele_cand);
2064 if (hasTrackDetails) {
2065 get(dnn::pfCand_ele_hasTrackDetails) =
2066 sp.scale(hasTrackDetails, dnn::pfCand_ele_hasTrackDetails - PFe_index_offset);
2067 get(dnn::pfCand_ele_dxy) = sp.scale(candFunc::getTauDxy(ele_cand), dnn::pfCand_ele_dxy - PFe_index_offset);
2068 get(dnn::pfCand_ele_dxy_sig) = sp.scale(
std::abs(candFunc::getTauDxy(ele_cand)) / ele_cand.dxyError(),
2069 dnn::pfCand_ele_dxy_sig - PFe_index_offset);
2070 get(dnn::pfCand_ele_dz) = sp.scale(candFunc::getTauDz(ele_cand), dnn::pfCand_ele_dz - PFe_index_offset);
2071 get(dnn::pfCand_ele_dz_sig) = sp.scale(
std::abs(candFunc::getTauDz(ele_cand)) / ele_cand.dzError(),
2072 dnn::pfCand_ele_dz_sig - PFe_index_offset);
2073 get(dnn::pfCand_ele_track_chi2_ndof) =
2074 candFunc::getPseudoTrack(ele_cand).ndof() > 0
2075 ? sp.scale(candFunc::getPseudoTrack(ele_cand).
chi2() / candFunc::getPseudoTrack(ele_cand).
ndof(),
2076 dnn::pfCand_ele_track_chi2_ndof - PFe_index_offset)
2078 get(dnn::pfCand_ele_track_ndof) =
2079 candFunc::getPseudoTrack(ele_cand).ndof() > 0
2080 ? sp.scale(candFunc::getPseudoTrack(ele_cand).
ndof(), dnn::pfCand_ele_track_ndof - PFe_index_offset)
2084 if (valid_index_pf_gamma) {
2085 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_PFg, is_inner));
2086 size_t index_pf_gamma = cell_map.at(CellObjectType::PfCand_gamma);
2087 const auto& gamma_cand =
dynamic_cast<const CandidateCastType&
>(pfCands.
at(index_pf_gamma));
2089 get(dnn::pfCand_gamma_valid + fill_index_offset_PFg) =
2090 sp.scale(valid_index_pf_gamma, dnn::pfCand_gamma_valid - PFg_index_offset);
2091 get(dnn::pfCand_gamma_rel_pt + fill_index_offset_PFg) =
2092 sp.scale(gamma_cand.polarP4().pt() /
tau.polarP4().pt(), dnn::pfCand_gamma_rel_pt - PFg_index_offset);
2093 get(dnn::pfCand_gamma_deta + fill_index_offset_PFg) =
2094 sp.scale(gamma_cand.polarP4().eta() -
tau.polarP4().eta(), dnn::pfCand_gamma_deta - PFg_index_offset);
2095 get(dnn::pfCand_gamma_dphi + fill_index_offset_PFg) =
2096 sp.scale(
dPhi(
tau.polarP4(), gamma_cand.polarP4()), dnn::pfCand_gamma_dphi - PFg_index_offset);
2097 get(dnn::pfCand_gamma_pvAssociationQuality + fill_index_offset_PFg) = sp.scale<
int>(
2098 candFunc::getPvAssocationQuality(gamma_cand), dnn::pfCand_gamma_pvAssociationQuality - PFg_index_offset);
2099 get(dnn::pfCand_gamma_fromPV + fill_index_offset_PFg) =
2100 sp.scale<
int>(candFunc::getFromPV(gamma_cand), dnn::pfCand_gamma_fromPV - PFg_index_offset);
2101 get(dnn::pfCand_gamma_puppiWeight + fill_index_offset_PFg) =
2102 is_inner ? sp.scale(candFunc::getPuppiWeight(gamma_cand, 0.9084110f),
2103 dnn::pfCand_gamma_puppiWeight - PFg_index_offset)
2104 : sp.scale(candFunc::getPuppiWeight(gamma_cand, 0.4211567
f),
2105 dnn::pfCand_gamma_puppiWeight - PFg_index_offset);
2106 get(dnn::pfCand_gamma_puppiWeightNoLep + fill_index_offset_PFg) =
2107 is_inner ? sp.scale(candFunc::getPuppiWeightNoLep(gamma_cand, 0.8857716f),
2108 dnn::pfCand_gamma_puppiWeightNoLep - PFg_index_offset)
2109 : sp.scale(candFunc::getPuppiWeightNoLep(gamma_cand, 0.3822604
f),
2110 dnn::pfCand_gamma_puppiWeightNoLep - PFg_index_offset);
2111 get(dnn::pfCand_gamma_lostInnerHits + fill_index_offset_PFg) =
2112 sp.scale<
int>(candFunc::getLostInnerHits(gamma_cand, 0), dnn::pfCand_gamma_lostInnerHits - PFg_index_offset);
2113 get(dnn::pfCand_gamma_numberOfPixelHits + fill_index_offset_PFg) = sp.scale(
2114 candFunc::getNumberOfPixelHits(gamma_cand, 0), dnn::pfCand_gamma_numberOfPixelHits - PFg_index_offset);
2115 get(dnn::pfCand_gamma_vertex_dx + fill_index_offset_PFg) =
2116 sp.scale(gamma_cand.vertex().x() -
pv.position().x(), dnn::pfCand_gamma_vertex_dx - PFg_index_offset);
2117 get(dnn::pfCand_gamma_vertex_dy + fill_index_offset_PFg) =
2118 sp.scale(gamma_cand.vertex().y() -
pv.position().y(), dnn::pfCand_gamma_vertex_dy - PFg_index_offset);
2119 get(dnn::pfCand_gamma_vertex_dz + fill_index_offset_PFg) =
2120 sp.scale(gamma_cand.vertex().z() -
pv.position().z(), dnn::pfCand_gamma_vertex_dz - PFg_index_offset);
2121 get(dnn::pfCand_gamma_vertex_dx_tauFL + fill_index_offset_PFg) =
2122 sp.scale(gamma_cand.vertex().x() -
pv.position().x() - tau_funcs.getFlightLength(
tau, tau_index).x(),
2123 dnn::pfCand_gamma_vertex_dx_tauFL - PFg_index_offset);
2124 get(dnn::pfCand_gamma_vertex_dy_tauFL + fill_index_offset_PFg) =
2125 sp.scale(gamma_cand.vertex().y() -
pv.position().y() - tau_funcs.getFlightLength(
tau, tau_index).y(),
2126 dnn::pfCand_gamma_vertex_dy_tauFL - PFg_index_offset);
2127 get(dnn::pfCand_gamma_vertex_dz_tauFL + fill_index_offset_PFg) =
2128 sp.scale(gamma_cand.vertex().z() -
pv.position().z() - tau_funcs.getFlightLength(
tau, tau_index).z(),
2129 dnn::pfCand_gamma_vertex_dz_tauFL - PFg_index_offset);
2130 const bool hasTrackDetails = candFunc::getHasTrackDetails(gamma_cand);
2131 if (hasTrackDetails) {
2132 get(dnn::pfCand_gamma_hasTrackDetails + fill_index_offset_PFg) =
2133 sp.scale(hasTrackDetails, dnn::pfCand_gamma_hasTrackDetails - PFg_index_offset);
2134 get(dnn::pfCand_gamma_dxy + fill_index_offset_PFg) =
2135 sp.scale(candFunc::getTauDxy(gamma_cand), dnn::pfCand_gamma_dxy - PFg_index_offset);
2136 get(dnn::pfCand_gamma_dxy_sig + fill_index_offset_PFg) =
2137 sp.scale(
std::abs(candFunc::getTauDxy(gamma_cand)) / gamma_cand.dxyError(),
2138 dnn::pfCand_gamma_dxy_sig - PFg_index_offset);
2139 get(dnn::pfCand_gamma_dz + fill_index_offset_PFg) =
2140 sp.scale(candFunc::getTauDz(gamma_cand), dnn::pfCand_gamma_dz - PFg_index_offset);
2141 get(dnn::pfCand_gamma_dz_sig + fill_index_offset_PFg) =
2142 sp.scale(
std::abs(candFunc::getTauDz(gamma_cand)) / gamma_cand.dzError(),
2143 dnn::pfCand_gamma_dz_sig - PFg_index_offset);
2144 get(dnn::pfCand_gamma_track_chi2_ndof + fill_index_offset_PFg) =
2145 candFunc::getPseudoTrack(gamma_cand).ndof() > 0
2146 ? sp.scale(candFunc::getPseudoTrack(gamma_cand).
chi2() / candFunc::getPseudoTrack(gamma_cand).
ndof(),
2147 dnn::pfCand_gamma_track_chi2_ndof - PFg_index_offset)
2149 get(dnn::pfCand_gamma_track_ndof + fill_index_offset_PFg) =
2150 candFunc::getPseudoTrack(gamma_cand).ndof() > 0
2151 ? sp.scale(candFunc::getPseudoTrack(gamma_cand).
ndof(), dnn::pfCand_gamma_track_ndof - PFg_index_offset)
2155 if (valid_index_ele) {
2156 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_e, is_inner));
2158 const auto& ele =
electrons->at(index_ele);
2160 get(dnn::ele_valid + fill_index_offset_e) = sp.scale(valid_index_ele, dnn::ele_valid - e_index_offset);
2161 get(dnn::ele_rel_pt + fill_index_offset_e) =
2162 sp.scale(ele.polarP4().pt() /
tau.polarP4().pt(), dnn::ele_rel_pt - e_index_offset);
2163 get(dnn::ele_deta + fill_index_offset_e) =
2164 sp.scale(ele.polarP4().eta() -
tau.polarP4().eta(), dnn::ele_deta - e_index_offset);
2165 get(dnn::ele_dphi + fill_index_offset_e) =
2166 sp.scale(
dPhi(
tau.polarP4(), ele.polarP4()), dnn::ele_dphi - e_index_offset);
2168 float cc_ele_energy, cc_gamma_energy;
2172 get(dnn::ele_cc_valid + fill_index_offset_e) = sp.scale(cc_valid, dnn::ele_cc_valid - e_index_offset);
2173 get(dnn::ele_cc_ele_rel_energy + fill_index_offset_e) =
2174 sp.scale(cc_ele_energy / ele.polarP4().pt(), dnn::ele_cc_ele_rel_energy - e_index_offset);
2175 get(dnn::ele_cc_gamma_rel_energy + fill_index_offset_e) =
2176 sp.scale(cc_gamma_energy / cc_ele_energy, dnn::ele_cc_gamma_rel_energy - e_index_offset);
2177 get(dnn::ele_cc_n_gamma + fill_index_offset_e) = sp.scale(cc_n_gamma, dnn::ele_cc_n_gamma - e_index_offset);
2179 get(dnn::ele_rel_trackMomentumAtVtx + fill_index_offset_e) =
2180 sp.scale(ele.trackMomentumAtVtx().R() / ele.polarP4().pt(), dnn::ele_rel_trackMomentumAtVtx - e_index_offset);
2181 get(dnn::ele_rel_trackMomentumAtCalo + fill_index_offset_e) = sp.scale(
2182 ele.trackMomentumAtCalo().R() / ele.polarP4().pt(), dnn::ele_rel_trackMomentumAtCalo - e_index_offset);
2183 get(dnn::ele_rel_trackMomentumOut + fill_index_offset_e) =
2184 sp.scale(ele.trackMomentumOut().R() / ele.polarP4().pt(), dnn::ele_rel_trackMomentumOut - e_index_offset);
2185 get(dnn::ele_rel_trackMomentumAtEleClus + fill_index_offset_e) = sp.scale(
2186 ele.trackMomentumAtEleClus().R() / ele.polarP4().pt(), dnn::ele_rel_trackMomentumAtEleClus - e_index_offset);
2187 get(dnn::ele_rel_trackMomentumAtVtxWithConstraint + fill_index_offset_e) =
2188 sp.scale(ele.trackMomentumAtVtxWithConstraint().R() / ele.polarP4().pt(),
2189 dnn::ele_rel_trackMomentumAtVtxWithConstraint - e_index_offset);
2190 get(dnn::ele_rel_ecalEnergy + fill_index_offset_e) =
2191 sp.scale(ele.ecalEnergy() / ele.polarP4().pt(), dnn::ele_rel_ecalEnergy - e_index_offset);
2192 get(dnn::ele_ecalEnergy_sig + fill_index_offset_e) =
2193 sp.scale(ele.ecalEnergy() / ele.ecalEnergyError(), dnn::ele_ecalEnergy_sig - e_index_offset);
2194 get(dnn::ele_eSuperClusterOverP + fill_index_offset_e) =
2195 sp.scale(ele.eSuperClusterOverP(), dnn::ele_eSuperClusterOverP - e_index_offset);
2196 get(dnn::ele_eSeedClusterOverP + fill_index_offset_e) =
2197 sp.scale(ele.eSeedClusterOverP(), dnn::ele_eSeedClusterOverP - e_index_offset);
2198 get(dnn::ele_eSeedClusterOverPout + fill_index_offset_e) =
2199 sp.scale(ele.eSeedClusterOverPout(), dnn::ele_eSeedClusterOverPout - e_index_offset);
2200 get(dnn::ele_eEleClusterOverPout + fill_index_offset_e) =
2201 sp.scale(ele.eEleClusterOverPout(), dnn::ele_eEleClusterOverPout - e_index_offset);
2202 get(dnn::ele_deltaEtaSuperClusterTrackAtVtx + fill_index_offset_e) =
2203 sp.scale(ele.deltaEtaSuperClusterTrackAtVtx(), dnn::ele_deltaEtaSuperClusterTrackAtVtx - e_index_offset);
2204 get(dnn::ele_deltaEtaSeedClusterTrackAtCalo + fill_index_offset_e) =
2205 sp.scale(ele.deltaEtaSeedClusterTrackAtCalo(), dnn::ele_deltaEtaSeedClusterTrackAtCalo - e_index_offset);
2206 get(dnn::ele_deltaEtaEleClusterTrackAtCalo + fill_index_offset_e) =
2207 sp.scale(ele.deltaEtaEleClusterTrackAtCalo(), dnn::ele_deltaEtaEleClusterTrackAtCalo - e_index_offset);
2208 get(dnn::ele_deltaPhiEleClusterTrackAtCalo + fill_index_offset_e) =
2209 sp.scale(ele.deltaPhiEleClusterTrackAtCalo(), dnn::ele_deltaPhiEleClusterTrackAtCalo - e_index_offset);
2210 get(dnn::ele_deltaPhiSuperClusterTrackAtVtx + fill_index_offset_e) =
2211 sp.scale(ele.deltaPhiSuperClusterTrackAtVtx(), dnn::ele_deltaPhiSuperClusterTrackAtVtx - e_index_offset);
2212 get(dnn::ele_deltaPhiSeedClusterTrackAtCalo + fill_index_offset_e) =
2213 sp.scale(ele.deltaPhiSeedClusterTrackAtCalo(), dnn::ele_deltaPhiSeedClusterTrackAtCalo - e_index_offset);
2214 const bool mva_valid =
2215 (ele.mvaInput().earlyBrem > -2) ||
2219 get(dnn::ele_mvaInput_earlyBrem + fill_index_offset_e) =
2220 sp.scale(ele.mvaInput().earlyBrem, dnn::ele_mvaInput_earlyBrem - e_index_offset);
2221 get(dnn::ele_mvaInput_lateBrem + fill_index_offset_e) =
2222 sp.scale(ele.mvaInput().lateBrem, dnn::ele_mvaInput_lateBrem - e_index_offset);
2223 get(dnn::ele_mvaInput_sigmaEtaEta + fill_index_offset_e) =
2224 sp.scale(ele.mvaInput().sigmaEtaEta, dnn::ele_mvaInput_sigmaEtaEta - e_index_offset);
2225 get(dnn::ele_mvaInput_hadEnergy + fill_index_offset_e) =
2226 sp.scale(ele.mvaInput().hadEnergy, dnn::ele_mvaInput_hadEnergy - e_index_offset);
2227 get(dnn::ele_mvaInput_deltaEta + fill_index_offset_e) =
2228 sp.scale(ele.mvaInput().deltaEta, dnn::ele_mvaInput_deltaEta - e_index_offset);
2230 const auto& gsfTrack = ele.gsfTrack();
2231 if (gsfTrack.isNonnull()) {
2232 get(dnn::ele_gsfTrack_normalizedChi2 + fill_index_offset_e) =
2233 sp.scale(gsfTrack->normalizedChi2(), dnn::ele_gsfTrack_normalizedChi2 - e_index_offset);
2234 get(dnn::ele_gsfTrack_numberOfValidHits + fill_index_offset_e) =
2235 sp.scale(gsfTrack->numberOfValidHits(), dnn::ele_gsfTrack_numberOfValidHits - e_index_offset);
2236 get(dnn::ele_rel_gsfTrack_pt + fill_index_offset_e) =
2237 sp.scale(gsfTrack->pt() / ele.polarP4().pt(), dnn::ele_rel_gsfTrack_pt - e_index_offset);
2238 get(dnn::ele_gsfTrack_pt_sig + fill_index_offset_e) =
2239 sp.scale(gsfTrack->pt() / gsfTrack->ptError(), dnn::ele_gsfTrack_pt_sig - e_index_offset);
2241 const auto& closestCtfTrack = ele.closestCtfTrackRef();
2242 const bool has_closestCtfTrack = closestCtfTrack.isNonnull();
2243 if (has_closestCtfTrack) {
2244 get(dnn::ele_has_closestCtfTrack + fill_index_offset_e) =
2245 sp.scale(has_closestCtfTrack, dnn::ele_has_closestCtfTrack - e_index_offset);
2246 get(dnn::ele_closestCtfTrack_normalizedChi2 + fill_index_offset_e) =
2247 sp.scale(closestCtfTrack->normalizedChi2(), dnn::ele_closestCtfTrack_normalizedChi2 - e_index_offset);
2248 get(dnn::ele_closestCtfTrack_numberOfValidHits + fill_index_offset_e) =
2249 sp.scale(closestCtfTrack->numberOfValidHits(), dnn::ele_closestCtfTrack_numberOfValidHits - e_index_offset);
2254 template <
typename Cand
idateCastType,
typename TauCastType>
2256 const TauCastType&
tau,
2257 const size_t tau_index,
2261 const std::vector<pat::Muon>*
muons,
2263 const Cell& cell_map,
2273 int PFmu_index_offset =
scalingParamsMap_->at(std::make_pair(ft_global,
false)).mean_.size();
2274 int mu_index_offset = PFmu_index_offset +
scalingParamsMap_->at(std::make_pair(ft_PFmu,
false)).mean_.size();
2278 const auto&
get = [&](
int var_index) ->
float& {
return inputs.tensor<
float, 4>()(
idx, 0, 0, var_index); };
2280 const bool valid_index_pf_muon = cell_map.count(CellObjectType::PfCand_muon);
2283 if (!cell_map.empty()) {
2284 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_global,
false));
2288 get(dnn::tau_inside_ecal_crack) = sp.scale(
isInEcalCrack(
tau.polarP4().eta()), dnn::tau_inside_ecal_crack);
2290 if (valid_index_pf_muon) {
2291 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_PFmu, is_inner));
2292 size_t index_pf_muon = cell_map.at(CellObjectType::PfCand_muon);
2293 const auto& muon_cand =
dynamic_cast<const CandidateCastType&
>(pfCands.
at(index_pf_muon));
2295 get(dnn::pfCand_muon_valid) = sp.scale(valid_index_pf_muon, dnn::pfCand_muon_valid - PFmu_index_offset);
2296 get(dnn::pfCand_muon_rel_pt) =
2297 sp.scale(muon_cand.polarP4().pt() /
tau.polarP4().pt(), dnn::pfCand_muon_rel_pt - PFmu_index_offset);
2298 get(dnn::pfCand_muon_deta) =
2299 sp.scale(muon_cand.polarP4().eta() -
tau.polarP4().eta(), dnn::pfCand_muon_deta - PFmu_index_offset);
2300 get(dnn::pfCand_muon_dphi) =
2301 sp.scale(
dPhi(
tau.polarP4(), muon_cand.polarP4()), dnn::pfCand_muon_dphi - PFmu_index_offset);
2302 get(dnn::pfCand_muon_pvAssociationQuality) = sp.scale<
int>(
2303 candFunc::getPvAssocationQuality(muon_cand), dnn::pfCand_muon_pvAssociationQuality - PFmu_index_offset);
2304 get(dnn::pfCand_muon_fromPV) =
2305 sp.scale<
int>(candFunc::getFromPV(muon_cand), dnn::pfCand_muon_fromPV - PFmu_index_offset);
2306 get(dnn::pfCand_muon_puppiWeight) = is_inner ? sp.scale(candFunc::getPuppiWeight(muon_cand, 0.9786588f),
2307 dnn::pfCand_muon_puppiWeight - PFmu_index_offset)
2308 : sp.scale(candFunc::getPuppiWeight(muon_cand, 0.8132477
f),
2309 dnn::pfCand_muon_puppiWeight - PFmu_index_offset);
2310 get(dnn::pfCand_muon_charge) = sp.scale(muon_cand.charge(), dnn::pfCand_muon_charge - PFmu_index_offset);
2311 get(dnn::pfCand_muon_lostInnerHits) =
2312 sp.scale<
int>(candFunc::getLostInnerHits(muon_cand, 0), dnn::pfCand_muon_lostInnerHits - PFmu_index_offset);
2313 get(dnn::pfCand_muon_numberOfPixelHits) = sp.scale(candFunc::getNumberOfPixelHits(muon_cand, 0),
2314 dnn::pfCand_muon_numberOfPixelHits - PFmu_index_offset);
2315 get(dnn::pfCand_muon_vertex_dx) =
2316 sp.scale(muon_cand.vertex().x() -
pv.position().x(), dnn::pfCand_muon_vertex_dx - PFmu_index_offset);
2317 get(dnn::pfCand_muon_vertex_dy) =
2318 sp.scale(muon_cand.vertex().y() -
pv.position().y(), dnn::pfCand_muon_vertex_dy - PFmu_index_offset);
2319 get(dnn::pfCand_muon_vertex_dz) =
2320 sp.scale(muon_cand.vertex().z() -
pv.position().z(), dnn::pfCand_muon_vertex_dz - PFmu_index_offset);
2321 get(dnn::pfCand_muon_vertex_dx_tauFL) =
2322 sp.scale(muon_cand.vertex().x() -
pv.position().x() - tau_funcs.getFlightLength(
tau, tau_index).x(),
2323 dnn::pfCand_muon_vertex_dx_tauFL - PFmu_index_offset);
2324 get(dnn::pfCand_muon_vertex_dy_tauFL) =
2325 sp.scale(muon_cand.vertex().y() -
pv.position().y() - tau_funcs.getFlightLength(
tau, tau_index).y(),
2326 dnn::pfCand_muon_vertex_dy_tauFL - PFmu_index_offset);
2327 get(dnn::pfCand_muon_vertex_dz_tauFL) =
2328 sp.scale(muon_cand.vertex().z() -
pv.position().z() - tau_funcs.getFlightLength(
tau, tau_index).z(),
2329 dnn::pfCand_muon_vertex_dz_tauFL - PFmu_index_offset);
2331 const bool hasTrackDetails = candFunc::getHasTrackDetails(muon_cand);
2332 if (hasTrackDetails) {
2333 get(dnn::pfCand_muon_hasTrackDetails) =
2334 sp.scale(hasTrackDetails, dnn::pfCand_muon_hasTrackDetails - PFmu_index_offset);
2335 get(dnn::pfCand_muon_dxy) = sp.scale(candFunc::getTauDxy(muon_cand), dnn::pfCand_muon_dxy - PFmu_index_offset);
2336 get(dnn::pfCand_muon_dxy_sig) = sp.scale(
std::abs(candFunc::getTauDxy(muon_cand)) / muon_cand.dxyError(),
2337 dnn::pfCand_muon_dxy_sig - PFmu_index_offset);
2338 get(dnn::pfCand_muon_dz) = sp.scale(candFunc::getTauDz(muon_cand), dnn::pfCand_muon_dz - PFmu_index_offset);
2339 get(dnn::pfCand_muon_dz_sig) = sp.scale(
std::abs(candFunc::getTauDz(muon_cand)) / muon_cand.dzError(),
2340 dnn::pfCand_muon_dz_sig - PFmu_index_offset);
2341 get(dnn::pfCand_muon_track_chi2_ndof) =
2342 candFunc::getPseudoTrack(muon_cand).ndof() > 0
2343 ? sp.scale(candFunc::getPseudoTrack(muon_cand).
chi2() / candFunc::getPseudoTrack(muon_cand).
ndof(),
2344 dnn::pfCand_muon_track_chi2_ndof - PFmu_index_offset)
2346 get(dnn::pfCand_muon_track_ndof) =
2347 candFunc::getPseudoTrack(muon_cand).ndof() > 0
2348 ? sp.scale(candFunc::getPseudoTrack(muon_cand).
ndof(), dnn::pfCand_muon_track_ndof - PFmu_index_offset)
2352 if (valid_index_muon) {
2353 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_mu, is_inner));
2355 const auto&
muon =
muons->at(index_muon);
2357 get(dnn::muon_valid) = sp.scale(valid_index_muon, dnn::muon_valid - mu_index_offset);
2358 get(dnn::muon_rel_pt) = sp.scale(
muon.polarP4().pt() /
tau.polarP4().pt(), dnn::muon_rel_pt - mu_index_offset);
2359 get(dnn::muon_deta) = sp.scale(
muon.polarP4().eta() -
tau.polarP4().eta(), dnn::muon_deta - mu_index_offset);
2360 get(dnn::muon_dphi) = sp.scale(
dPhi(
tau.polarP4(),
muon.polarP4()), dnn::muon_dphi - mu_index_offset);
2362 get(dnn::muon_dxy_sig) =
2365 const bool normalizedChi2_valid =
muon.globalTrack().isNonnull() &&
muon.normChi2() >= 0;
2366 if (normalizedChi2_valid) {
2367 get(dnn::muon_normalizedChi2_valid) =
2368 sp.scale(normalizedChi2_valid, dnn::muon_normalizedChi2_valid - mu_index_offset);
2369 get(dnn::muon_normalizedChi2) = sp.scale(
muon.normChi2(), dnn::muon_normalizedChi2 - mu_index_offset);
2370 if (
muon.innerTrack().isNonnull())
2371 get(dnn::muon_numberOfValidHits) =
2372 sp.scale(
muon.numberOfValidHits(), dnn::muon_numberOfValidHits - mu_index_offset);
2374 get(dnn::muon_segmentCompatibility) =
2375 sp.scale(
muon.segmentCompatibility(), dnn::muon_segmentCompatibility - mu_index_offset);
2376 get(dnn::muon_caloCompatibility) =
2377 sp.scale(
muon.caloCompatibility(), dnn::muon_caloCompatibility - mu_index_offset);
2379 const bool pfEcalEnergy_valid =
muon.pfEcalEnergy() >= 0;
2380 if (pfEcalEnergy_valid) {
2381 get(dnn::muon_pfEcalEnergy_valid) =
2382 sp.scale(pfEcalEnergy_valid, dnn::muon_pfEcalEnergy_valid - mu_index_offset);
2383 get(dnn::muon_rel_pfEcalEnergy) =
2384 sp.scale(
muon.pfEcalEnergy() /
muon.polarP4().pt(), dnn::muon_rel_pfEcalEnergy - mu_index_offset);
2387 MuonHitMatchV2 hit_match(
muon);
2388 static const std::map<int, std::pair<int, int>> muonMatchHitVars = {
2393 for (
int subdet : hit_match.MuonHitMatchV2::consideredSubdets()) {
2394 const auto& matchHitVar = muonMatchHitVars.at(subdet);
2395 for (
int station = MuonHitMatchV2::first_station_id;
station <= MuonHitMatchV2::last_station_id; ++
station) {
2396 const unsigned n_matches = hit_match.nMatches(subdet,
station);
2397 const unsigned n_hits = hit_match.nHits(subdet,
station);
2398 get(matchHitVar.first +
station - 1) = sp.scale(n_matches, matchHitVar.first +
station - 1 - mu_index_offset);
2399 get(matchHitVar.second +
station - 1) = sp.scale(n_hits, matchHitVar.second +
station - 1 - mu_index_offset);
2405 template <
typename Cand
idateCastType,
typename TauCastType>
2407 const TauCastType&
tau,
2408 const size_t tau_index,
2413 const Cell& cell_map,
2423 int PFchH_index_offset =
scalingParamsMap_->at(std::make_pair(ft_global,
false)).mean_.size();
2424 int PFnH_index_offset = PFchH_index_offset +
scalingParamsMap_->at(std::make_pair(ft_PFchH,
false)).mean_.size();
2428 const auto&
get = [&](
int var_index) ->
float& {
return inputs.tensor<
float, 4>()(
idx, 0, 0, var_index); };
2430 const bool valid_chH = cell_map.count(CellObjectType::PfCand_chargedHadron);
2431 const bool valid_nH = cell_map.count(CellObjectType::PfCand_neutralHadron);
2433 if (!cell_map.empty()) {
2434 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_global,
false));
2438 get(dnn::tau_inside_ecal_crack) = sp.scale(
isInEcalCrack(
tau.polarP4().eta()), dnn::tau_inside_ecal_crack);
2441 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_PFchH, is_inner));
2442 size_t index_chH = cell_map.at(CellObjectType::PfCand_chargedHadron);
2443 const auto& chH_cand =
dynamic_cast<const CandidateCastType&
>(pfCands.
at(index_chH));
2445 get(dnn::pfCand_chHad_valid) = sp.scale(valid_chH, dnn::pfCand_chHad_valid - PFchH_index_offset);
2446 get(dnn::pfCand_chHad_rel_pt) =
2447 sp.scale(chH_cand.polarP4().pt() /
tau.polarP4().pt(), dnn::pfCand_chHad_rel_pt - PFchH_index_offset);
2448 get(dnn::pfCand_chHad_deta) =
2449 sp.scale(chH_cand.polarP4().eta() -
tau.polarP4().eta(), dnn::pfCand_chHad_deta - PFchH_index_offset);
2450 get(dnn::pfCand_chHad_dphi) =
2451 sp.scale(
dPhi(
tau.polarP4(), chH_cand.polarP4()), dnn::pfCand_chHad_dphi - PFchH_index_offset);
2452 get(dnn::pfCand_chHad_leadChargedHadrCand) =
2453 sp.scale(&chH_cand == dynamic_cast<const CandidateCastType*>(
tau.leadChargedHadrCand().get()),
2454 dnn::pfCand_chHad_leadChargedHadrCand - PFchH_index_offset);
2455 get(dnn::pfCand_chHad_pvAssociationQuality) = sp.scale<
int>(
2456 candFunc::getPvAssocationQuality(chH_cand), dnn::pfCand_chHad_pvAssociationQuality - PFchH_index_offset);
2457 get(dnn::pfCand_chHad_fromPV) =
2458 sp.scale<
int>(candFunc::getFromPV(chH_cand), dnn::pfCand_chHad_fromPV - PFchH_index_offset);
2459 const float default_chH_pw_inner = 0.7614090f;
2460 const float default_chH_pw_outer = 0.1974930f;
2461 get(dnn::pfCand_chHad_puppiWeight) = is_inner ? sp.scale(candFunc::getPuppiWeight(chH_cand, default_chH_pw_inner),
2462 dnn::pfCand_chHad_puppiWeight - PFchH_index_offset)
2463 : sp.scale(candFunc::getPuppiWeight(chH_cand, default_chH_pw_outer),
2464 dnn::pfCand_chHad_puppiWeight - PFchH_index_offset);
2465 get(dnn::pfCand_chHad_puppiWeightNoLep) =
2466 is_inner ? sp.scale(candFunc::getPuppiWeightNoLep(chH_cand, default_chH_pw_inner),
2467 dnn::pfCand_chHad_puppiWeightNoLep - PFchH_index_offset)
2468 : sp.scale(candFunc::getPuppiWeightNoLep(chH_cand, default_chH_pw_outer),
2469 dnn::pfCand_chHad_puppiWeightNoLep - PFchH_index_offset);
2470 get(dnn::pfCand_chHad_charge) = sp.scale(chH_cand.charge(), dnn::pfCand_chHad_charge - PFchH_index_offset);
2471 get(dnn::pfCand_chHad_lostInnerHits) =
2472 sp.scale<
int>(candFunc::getLostInnerHits(chH_cand, 0), dnn::pfCand_chHad_lostInnerHits - PFchH_index_offset);
2473 get(dnn::pfCand_chHad_numberOfPixelHits) = sp.scale(candFunc::getNumberOfPixelHits(chH_cand, 0),
2474 dnn::pfCand_chHad_numberOfPixelHits - PFchH_index_offset);
2475 get(dnn::pfCand_chHad_vertex_dx) =
2476 sp.scale(chH_cand.vertex().x() -
pv.position().x(), dnn::pfCand_chHad_vertex_dx - PFchH_index_offset);
2477 get(dnn::pfCand_chHad_vertex_dy) =
2478 sp.scale(chH_cand.vertex().y() -
pv.position().y(), dnn::pfCand_chHad_vertex_dy - PFchH_index_offset);
2479 get(dnn::pfCand_chHad_vertex_dz) =
2480 sp.scale(chH_cand.vertex().z() -
pv.position().z(), dnn::pfCand_chHad_vertex_dz - PFchH_index_offset);
2481 get(dnn::pfCand_chHad_vertex_dx_tauFL) =
2482 sp.scale(chH_cand.vertex().x() -
pv.position().x() - tau_funcs.getFlightLength(
tau, tau_index).x(),
2483 dnn::pfCand_chHad_vertex_dx_tauFL - PFchH_index_offset);
2484 get(dnn::pfCand_chHad_vertex_dy_tauFL) =
2485 sp.scale(chH_cand.vertex().y() -
pv.position().y() - tau_funcs.getFlightLength(
tau, tau_index).y(),
2486 dnn::pfCand_chHad_vertex_dy_tauFL - PFchH_index_offset);
2487 get(dnn::pfCand_chHad_vertex_dz_tauFL) =
2488 sp.scale(chH_cand.vertex().z() -
pv.position().z() - tau_funcs.getFlightLength(
tau, tau_index).z(),
2489 dnn::pfCand_chHad_vertex_dz_tauFL - PFchH_index_offset);
2491 const bool hasTrackDetails = candFunc::getHasTrackDetails(chH_cand);
2492 if (hasTrackDetails) {
2493 get(dnn::pfCand_chHad_hasTrackDetails) =
2494 sp.scale(hasTrackDetails, dnn::pfCand_chHad_hasTrackDetails - PFchH_index_offset);
2495 get(dnn::pfCand_chHad_dxy) =
2496 sp.scale(candFunc::getTauDxy(chH_cand), dnn::pfCand_chHad_dxy - PFchH_index_offset);
2497 get(dnn::pfCand_chHad_dxy_sig) = sp.scale(
std::abs(candFunc::getTauDxy(chH_cand)) / chH_cand.dxyError(),
2498 dnn::pfCand_chHad_dxy_sig - PFchH_index_offset);
2499 get(dnn::pfCand_chHad_dz) = sp.scale(candFunc::getTauDz(chH_cand), dnn::pfCand_chHad_dz - PFchH_index_offset);
2500 get(dnn::pfCand_chHad_dz_sig) = sp.scale(
std::abs(candFunc::getTauDz(chH_cand)) / chH_cand.dzError(),
2501 dnn::pfCand_chHad_dz_sig - PFchH_index_offset);
2502 get(dnn::pfCand_chHad_track_chi2_ndof) =
2503 candFunc::getPseudoTrack(chH_cand).ndof() > 0
2504 ? sp.scale(candFunc::getPseudoTrack(chH_cand).
chi2() / candFunc::getPseudoTrack(chH_cand).
ndof(),
2505 dnn::pfCand_chHad_track_chi2_ndof - PFchH_index_offset)
2507 get(dnn::pfCand_chHad_track_ndof) =
2508 candFunc::getPseudoTrack(chH_cand).ndof() > 0
2509 ? sp.scale(candFunc::getPseudoTrack(chH_cand).
ndof(), dnn::pfCand_chHad_track_ndof - PFchH_index_offset)
2513 get(dnn::pfCand_chHad_hcalFraction) =
2514 sp.scale(hcal_fraction, dnn::pfCand_chHad_hcalFraction - PFchH_index_offset);
2515 get(dnn::pfCand_chHad_rawCaloFraction) =
2516 sp.scale(candFunc::getRawCaloFraction(chH_cand), dnn::pfCand_chHad_rawCaloFraction - PFchH_index_offset);
2519 const sc::ScalingParams& sp =
scalingParamsMap_->at(std::make_pair(ft_PFnH, is_inner));
2520 size_t index_nH = cell_map.at(CellObjectType::PfCand_neutralHadron);
2521 const auto& nH_cand =
dynamic_cast<const CandidateCastType&
>(pfCands.
at(index_nH));
2523 get(dnn::pfCand_nHad_valid) = sp.scale(valid_nH, dnn::pfCand_nHad_valid - PFnH_index_offset);
2524 get(dnn::pfCand_nHad_rel_pt) =
2525 sp.scale(nH_cand.polarP4().pt() /
tau.polarP4().pt(), dnn::pfCand_nHad_rel_pt - PFnH_index_offset);
2526 get(dnn::pfCand_nHad_deta) =
2527 sp.scale(nH_cand.polarP4().eta() -
tau.polarP4().eta(), dnn::pfCand_nHad_deta - PFnH_index_offset);
2528 get(dnn::pfCand_nHad_dphi) =
2529 sp.scale(
dPhi(
tau.polarP4(), nH_cand.polarP4()), dnn::pfCand_nHad_dphi - PFnH_index_offset);
2530 get(dnn::pfCand_nHad_puppiWeight) = is_inner ? sp.scale(candFunc::getPuppiWeight(nH_cand, 0.9798355f),
2531 dnn::pfCand_nHad_puppiWeight - PFnH_index_offset)
2532 : sp.scale(candFunc::getPuppiWeight(nH_cand, 0.7813260
f),
2533 dnn::pfCand_nHad_puppiWeight - PFnH_index_offset);
2534 get(dnn::pfCand_nHad_puppiWeightNoLep) = is_inner
2535 ? sp.scale(candFunc::getPuppiWeightNoLep(nH_cand, 0.9046796f),
2536 dnn::pfCand_nHad_puppiWeightNoLep - PFnH_index_offset)
2537 : sp.scale(candFunc::getPuppiWeightNoLep(nH_cand, 0.6554860
f),
2538 dnn::pfCand_nHad_puppiWeightNoLep - PFnH_index_offset);
2540 get(dnn::pfCand_nHad_hcalFraction) = sp.scale(hcal_fraction, dnn::pfCand_nHad_hcalFraction - PFnH_index_offset);
2546 elecEe = elecEgamma = 0;
2548 if (superCluster.isNonnull() && superCluster.isAvailable() && superCluster->clusters().isNonnull() &&
2549 superCluster->clusters().isAvailable()) {
2550 for (
auto iter = superCluster->clustersBegin(); iter != superCluster->clustersEnd(); ++iter) {
2551 const double energy = (*iter)->energy();
2552 if (iter == superCluster->clustersBegin())
2563 template <
typename Cand
idateCollection,
typename TauCastType>
2586 const bool isInside_innerSigCone =
dR < innerSigCone_radius;
2587 if (isInside_innerSigCone) {
2588 p4_inner +=
cand->p4();
2591 p4_outer +=
cand->p4();
2607 template <
typename Cand
idateCollection,
typename TauCastType>
2631 static constexpr
double min_pt = 30., min_radius = 0.05, cone_opening_coef = 3.;
2637 template <
typename TauCastType>
2639 const size_t tau_index,
2641 TauFunc tau_funcs) {
2642 if (tau_funcs.getHasSecondaryVertex(
tau, tau_index)) {
2643 static constexpr
double mTau = 1.77682;
2644 const double mAOne =
tau.p4().M();
2645 const double pAOneMag =
tau.p();
2646 const double argumentThetaGJmax = (
std::pow(mTau, 2) -
std::pow(mAOne, 2)) / (2 * mTau * pAOneMag);
2647 const double argumentThetaGJmeasured =
tau.p4().Vect().Dot(tau_funcs.getFlightLength(
tau, tau_index)) /
2648 (pAOneMag * tau_funcs.getFlightLength(
tau, tau_index).R());
2649 if (
std::abs(argumentThetaGJmax) <= 1. &&
std::abs(argumentThetaGJmeasured) <= 1.) {
2650 double thetaGJmax = std::asin(argumentThetaGJmax);
2651 double thetaGJmeasured = std::acos(argumentThetaGJmeasured);
2652 gj_diff = thetaGJmeasured - thetaGJmax;
2659 template <
typename TauCastType>
2661 const size_t tau_index,
2662 TauFunc tau_funcs) {
2665 return static_cast<float>(gj_diff);
2671 return abs_eta > 1.46 && abs_eta < 1.558;
2674 template <
typename TauCastType>
2676 const std::vector<pat::Electron>*
electrons,
2681 if (
reco::deltaR2(
tau.p4(), ele.p4()) < dR2 && (!matched_ele || matched_ele->
pt() < ele.pt())) {
void createConvFeatures(const TauCastType &tau, const size_t tau_index, const edm::RefToBase< reco::BaseTau > tau_ref, const reco::Vertex &pv, double rho, const std::vector< pat::Electron > *electrons, const std::vector< pat::Muon > *muons, const edm::View< reco::Candidate > &pfCands, const CellGrid &grid, TauFunc tau_funcs, bool is_inner)
static constexpr float default_value
constexpr double deltaPhi(double phi1, double phi2)
std::unique_ptr< Cutter > CutterPtr
void createMuonBlockInputs(unsigned idx, const TauCastType &tau, const size_t tau_index, const edm::RefToBase< reco::BaseTau > tau_ref, const reco::Vertex &pv, double rho, const std::vector< pat::Muon > *muons, const edm::View< reco::Candidate > &pfCands, const Cell &cell_map, TauFunc tau_funcs, bool is_inner)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
tensorflow::Tensor getPredictions(edm::Event &event, edm::Handle< TauCollection > taus)
unsigned int n_photons_total(const reco::PFTau &tau)
return total number of pf photon candidates with pT>500 MeV, which are associated to signal ...
T getParameter(std::string const &) const
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > hadronsTensor_
std::map< BasicDiscriminator, size_t > basicDiscrdR03IndexMap_
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
std::map< std::string, WPList > workingPoints_
static void processIsolationPFComponents(const TauCastType &tau, const CandidateCollection &candidates, LorentzVectorXYZ &p4, float &pt, float &d_eta, float &d_phi, float &m, float &n)
static std::unique_ptr< deep_tau::DeepTauCache > initializeGlobalCache(const edm::ParameterSet &cfg)
double pt() const final
transverse momentum
static void processSignalPFComponents(const TauCastType &tau, const CandidateCollection &candidates, LorentzVectorXYZ &p4_inner, LorentzVectorXYZ &p4_outer, float &pt_inner, float &dEta_inner, float &dPhi_inner, float &m_inner, float &pt_outer, float &dEta_outer, float &dPhi_outer, float &m_outer, float &n_inner, float &n_outer)
DeepTauCache(const std::map< std::string, std::string > &graph_names, bool mem_mapped)
float pt_weighted_dr_signal(const reco::PFTau &tau, int dm)
std::string fullPath() const
bool operator<(IOVSyncValue const &iLHS, IOVSyncValue const &iRHS)
const bool disable_dxy_pca_
static float getValueLinear(T value, float min_value, float max_value, bool positive)
static const pat::Electron * findMatchedElectron(const TauCastType &tau, const std::vector< pat::Electron > *electrons, double deltaR)
constexpr bool isNotFinite(T x)
void createHadronsBlockInputs(unsigned idx, const TauCastType &tau, const size_t tau_index, const edm::RefToBase< reco::BaseTau > tau_ref, const reco::Vertex &pv, double rho, const edm::View< reco::Candidate > &pfCands, const Cell &cell_map, TauFunc tau_funcs, bool is_inner)
void fillGrids(const TauCastType &tau, const Collection &objects, CellGrid &inner_grid, CellGrid &outer_grid)
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > zeroOutputTensor_
std::map< BasicDiscriminator, size_t > basicDiscrIndexMap_
std::map< std::string, tensorflow::Session * > sessions_
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
GraphDef * loadGraphDef(const std::string &pbFile)
static double getInnerSignalConeRadius(double pt)
std::vector< TauDiscInfo< reco::PFTauDiscriminator > > recoPrediscriminants_
uint8_t andPrediscriminants_
std::ofstream * json_file_
unsigned long long EventNumber_t
OutputCollection outputs_
static bool isAbove(double value, double min)
static bool isInEcalCrack(double eta)
void setCellConvFeatures(tensorflow::Tensor &convTensor, const tensorflow::Tensor &features, unsigned batch_idx, int eta_index, int phi_index)
std::vector< TauDiscInfo< pat::PATTauDiscriminator > > patPrediscriminants_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
pat::ElectronCollection ElectronCollection
const unsigned sub_version_
pfTauTransverseImpactParameters
reco::SuperClusterRef superCluster() const override
override the reco::GsfElectron::superCluster method, to access the internal storage of the superclust...
void countMatches(const reco::Muon &muon, std::vector< int > &numMatchesDT, std::vector< int > &numMatchesCSC, std::vector< int > &numMatchesRPC)
void getPredictionsV2(TauCollection::const_reference &tau, const size_t tau_index, const edm::RefToBase< reco::BaseTau > tau_ref, const std::vector< pat::Electron > *electrons, const std::vector< pat::Muon > *muons, const edm::View< reco::Candidate > &pfCands, const reco::Vertex &pv, double rho, const edm::EventNumber_t &eventnr, std::vector< tensorflow::Tensor > &pred_vector, TauFunc tau_funcs)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
std::vector< Vertex > VertexCollection
edm::RefProd< PFTauCollection > PFTauRefProd
references to PFTau collection
const std::map< ValueQuantityType, double > min_value
ParameterSet const & parameterSet(StableProvenance const &provenance, ProcessHistory const &history)
edm::EDGetTokenT< TauCollection > tausToken_
static const double deltaEta
static std::string to_string(const XMLCh *ch)
const std::map< std::pair< deep_tau::Scaling::FeatureT, bool >, deep_tau::Scaling::ScalingParams > * scalingParamsMap_
static const OutputCollection & GetOutputs()
static std::string const input
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > convTensor_
DeepTauId(const edm::ParameterSet &cfg, const deep_tau::DeepTauCache *cache)
const std::map< std::pair< FeatureT, bool >, ScalingParams > scalingParamsMap_v2p5
std::unique_ptr< TauDiscriminator > get_value(const edm::Handle< TauCollection > &taus, const tensorflow::Tensor &pred, const WPList *working_points, bool is_online) const
void fill(const edm::Event &evt)
const deep_tau::DeepTauCache * cache_
edm::EDGetTokenT< CandidateCollection > pfcandToken_
std::string output_layer_
std::vector< Electron > ElectronCollection
std::shared_ptr< tensorflow::GraphDef > GraphPtr
static constexpr float pi
static float getValue(T value)
const tensorflow::GraphDef & getGraph(const std::string &name="") const
float pt_weighted_deta_strip(const reco::PFTau &tau, int dm)
std::vector< float > features(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
const std::map< BasicDiscriminator, size_t > matchDiscriminatorIndices(edm::Event &event, edm::EDGetTokenT< reco::TauDiscriminatorContainer > discriminatorContainerToken, std::vector< BasicDiscriminator > requiredDiscr)
edm::Handle< ConsumeType > handle
std::vector< std::string > getParameterNamesForType(bool trackiness=true) const
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
void saveInputs(const tensorflow::Tensor &inputs, const std::string &block_name, int n_inputs, const CellGrid *grid=nullptr)
const bool disable_CellIndex_workaround_
bool closeSession(Session *&session)
static float getValueNorm(T value, float mean, float sigma, float n_sigmas_max=5)
Abs< T >::type abs(const T &t)
float pt_weighted_dr_iso(const reco::PFTau &tau, int dm)
basicTauDiscriminatorsdR03
void createEgammaBlockInputs(unsigned idx, const TauCastType &tau, const size_t tau_index, const edm::RefToBase< reco::BaseTau > tau_ref, const reco::Vertex &pv, double rho, const std::vector< pat::Electron > *electrons, const edm::View< reco::Candidate > &pfCands, const Cell &cell_map, TauFunc tau_funcs, bool is_inner)
#define DEFINE_FWK_MODULE(type)
static void globalEndJob(const deep_tau::DeepTauCache *cache_)
std::vector< size_t > den_
std::vector< size_t > num_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void checkInputs(const tensorflow::Tensor &inputs, const std::string &block_name, int n_inputs, const CellGrid *grid=nullptr) const
tensorflow::Session & getSession(const std::string &name="") const
const_reference at(size_type pos) const
Analysis-level tau class.
const std::map< std::pair< FeatureT, bool >, ScalingParams > scalingParamsMap_v2p1
Session * createSession()
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
const std::vector< BasicDiscriminator > requiredBasicDiscriminators_
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > LorentzVectorXYZ
tensorflow::Tensor getPartialPredictions(bool is_inner)
const std::vector< BasicDiscriminator > requiredBasicDiscriminatorsdR03_
static bool calculateGottfriedJacksonAngleDifference(const TauCastType &tau, const size_t tau_index, double &gj_diff, TauFunc tau_funcs)
edm::EDGetTokenT< edm::AssociationVector< reco::PFTauRefProd, std::vector< reco::PFTauTransverseImpactParameterRef > > > pfTauTransverseImpactParameters_token_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
edm::ValueMap< SingleTauDiscriminatorContainer > TauDiscriminatorContainer
void createTauBlockInputs(const TauCastType &tau, const size_t &tau_index, const edm::RefToBase< reco::BaseTau > tau_ref, const reco::Vertex &pv, double rho, TauFunc tau_funcs)
float pt_weighted_dphi_strip(const reco::PFTau &tau, int dm)
float eratio(const reco::PFTau &tau)
return ratio of energy in ECAL over sum of energy in ECAL and HCAL
void createOutputs(edm::Event &event, const tensorflow::Tensor &pred, edm::Handle< TauCollection > taus)
Analysis-level electron class.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
disable_CellIndex_workaround
pat::MuonCollection MuonCollection
edm::Ref< PFTauTransverseImpactParameterCollection > PFTauTransverseImpactParameterRef
presistent reference to a PFTauTransverseImpactParameter
edm::EDGetTokenT< reco::TauDiscriminatorContainer > basicTauDiscriminators_inputToken_
std::vector< Muon > MuonCollection
static float calculateGottfriedJacksonAngleDifference(const TauCastType &tau, const size_t tau_index, TauFunc tau_funcs)
static bool calculateElectronClusterVarsV2(const pat::Electron &ele, float &cc_ele_energy, float &cc_gamma_energy, int &cc_n_gamma)
Particle reconstructed by the particle flow algorithm.
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > eGammaTensor_
const bool disable_hcalFraction_workaround_
edm::EDGetTokenT< reco::TauDiscriminatorContainer > basicTauDiscriminatorsdR03_inputToken_
std::map< std::string, GraphPtr > graphs_
const std::map< ValueQuantityType, double > max_value
disable_hcalFraction_workaround
const std::map< BasicDiscriminator, std::string > stringFromDiscriminator_
std::vector< int > tauInputs_indices_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
float lead_track_chi2(const reco::PFTau &tau)
return chi2 of the leading track ==> deprecated? <==
std::map< std::string, std::unique_ptr< tensorflow::MemmappedEnv > > memmappedEnv_
Log< level::Warning, false > LogWarning
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > muonTensor_
edm::EDGetTokenT< std::vector< pat::Electron > > electrons_token_
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
Output(const std::vector< size_t > &num, const std::vector< size_t > &den)
void produce(edm::Event &event, const edm::EventSetup &es) override
edm::EDGetTokenT< std::vector< pat::Muon > > muons_token_
T operator[](int i) const
edm::EDGetTokenT< double > rho_token_
static void calculateElectronClusterVars(const pat::Electron *ele, float &elecEe, float &elecEgamma)
T const & const_reference
static const std::string subdets[7]
const std::map< std::pair< FeatureT, bool >, ScalingParams > scalingParamsMap_PhaseIIv2p5
Analysis-level muon class.
std::map< std::string, Output > OutputCollection
void countHits(const reco::Muon &muon, std::vector< int > &numHitsDT, std::vector< int > &numHitsCSC, std::vector< int > &numHitsRPC)
std::vector< CutterPtr > WPList
edm::EDGetTokenT< ConsumeType > disc_token
constexpr int NumberOfOutputs
std::unique_ptr< tensorflow::Tensor > tauBlockTensor_