Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-31 02:19:52

0001 #include <memory>
0002 
0003 #include "RecoBTag/SecondaryVertex/interface/CandidateBoostedDoubleSecondaryVertexComputer.h"
0004 
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 
0007 #include "CondFormats/DataRecord/interface/BTauGenericMVAJetTagComputerRcd.h"
0008 #include "DataFormats/BTauReco/interface/BoostedDoubleSVTagInfo.h"
0009 #include "RecoBTau/JetTagComputer/interface/JetTagComputerRecord.h"
0010 
0011 CandidateBoostedDoubleSecondaryVertexComputer::Tokens::Tokens(const edm::ParameterSet& parameters,
0012                                                               edm::ESConsumesCollector&& cc) {
0013   if (parameters.getParameter<bool>("useCondDB")) {
0014     gbrForest_ = cc.consumes(edm::ESInputTag{"", parameters.getParameter<std::string>("gbrForestLabel")});
0015   }
0016 }
0017 
0018 CandidateBoostedDoubleSecondaryVertexComputer::CandidateBoostedDoubleSecondaryVertexComputer(
0019     const edm::ParameterSet& parameters, Tokens tokens)
0020     : weightFile_(parameters.getParameter<edm::FileInPath>("weightFile")),
0021       useGBRForest_(parameters.getParameter<bool>("useGBRForest")),
0022       useAdaBoost_(parameters.getParameter<bool>("useAdaBoost")),
0023       tokens_{tokens} {
0024   uses(0, "svTagInfos");
0025 
0026   mvaID = std::make_unique<TMVAEvaluator>();
0027 }
0028 
0029 void CandidateBoostedDoubleSecondaryVertexComputer::initialize(const JetTagComputerRecord& record) {
0030   // variable names and order need to be the same as in the training
0031   std::vector<std::string> variables({"z_ratio",
0032                                       "trackSipdSig_3",
0033                                       "trackSipdSig_2",
0034                                       "trackSipdSig_1",
0035                                       "trackSipdSig_0",
0036                                       "trackSipdSig_1_0",
0037                                       "trackSipdSig_0_0",
0038                                       "trackSipdSig_1_1",
0039                                       "trackSipdSig_0_1",
0040                                       "trackSip2dSigAboveCharm_0",
0041                                       "trackSip2dSigAboveBottom_0",
0042                                       "trackSip2dSigAboveBottom_1",
0043                                       "tau0_trackEtaRel_0",
0044                                       "tau0_trackEtaRel_1",
0045                                       "tau0_trackEtaRel_2",
0046                                       "tau1_trackEtaRel_0",
0047                                       "tau1_trackEtaRel_1",
0048                                       "tau1_trackEtaRel_2",
0049                                       "tau_vertexMass_0",
0050                                       "tau_vertexEnergyRatio_0",
0051                                       "tau_vertexDeltaR_0",
0052                                       "tau_flightDistance2dSig_0",
0053                                       "tau_vertexMass_1",
0054                                       "tau_vertexEnergyRatio_1",
0055                                       "tau_flightDistance2dSig_1",
0056                                       "jetNTracks",
0057                                       "nSV"});
0058   // book TMVA readers
0059   std::vector<std::string> spectators({"massPruned", "flavour", "nbHadrons", "ptPruned", "etaPruned"});
0060 
0061   if (tokens_.gbrForest_.isInitialized()) {
0062     mvaID->initializeGBRForest(&record.get(tokens_.gbrForest_), variables, spectators, useAdaBoost_);
0063   } else
0064     mvaID->initialize(
0065         "Color:Silent:Error", "BDT", weightFile_.fullPath(), variables, spectators, useGBRForest_, useAdaBoost_);
0066 }
0067 
0068 float CandidateBoostedDoubleSecondaryVertexComputer::discriminator(const TagInfoHelper& tagInfo) const {
0069   // get the TagInfo
0070   const reco::BoostedDoubleSVTagInfo& bdsvTagInfo = tagInfo.get<reco::BoostedDoubleSVTagInfo>(0);
0071 
0072   // get the TaggingVariables
0073   const reco::TaggingVariableList vars = bdsvTagInfo.taggingVariables();
0074 
0075   // default discriminator value
0076   float value = -10.;
0077 
0078   std::map<std::string, float> inputs;
0079   inputs["z_ratio"] = vars.get(reco::btau::z_ratio);
0080   inputs["trackSipdSig_3"] = vars.get(reco::btau::trackSip3dSig_3);
0081   inputs["trackSipdSig_2"] = vars.get(reco::btau::trackSip3dSig_2);
0082   inputs["trackSipdSig_1"] = vars.get(reco::btau::trackSip3dSig_1);
0083   inputs["trackSipdSig_0"] = vars.get(reco::btau::trackSip3dSig_0);
0084   inputs["trackSipdSig_1_0"] = vars.get(reco::btau::tau2_trackSip3dSig_0);
0085   inputs["trackSipdSig_0_0"] = vars.get(reco::btau::tau1_trackSip3dSig_0);
0086   inputs["trackSipdSig_1_1"] = vars.get(reco::btau::tau2_trackSip3dSig_1);
0087   inputs["trackSipdSig_0_1"] = vars.get(reco::btau::tau1_trackSip3dSig_1);
0088   inputs["trackSip2dSigAboveCharm_0"] = vars.get(reco::btau::trackSip2dSigAboveCharm);
0089   inputs["trackSip2dSigAboveBottom_0"] = vars.get(reco::btau::trackSip2dSigAboveBottom_0);
0090   inputs["trackSip2dSigAboveBottom_1"] = vars.get(reco::btau::trackSip2dSigAboveBottom_1);
0091   inputs["tau1_trackEtaRel_0"] = vars.get(reco::btau::tau2_trackEtaRel_0);
0092   inputs["tau1_trackEtaRel_1"] = vars.get(reco::btau::tau2_trackEtaRel_1);
0093   inputs["tau1_trackEtaRel_2"] = vars.get(reco::btau::tau2_trackEtaRel_2);
0094   inputs["tau0_trackEtaRel_0"] = vars.get(reco::btau::tau1_trackEtaRel_0);
0095   inputs["tau0_trackEtaRel_1"] = vars.get(reco::btau::tau1_trackEtaRel_1);
0096   inputs["tau0_trackEtaRel_2"] = vars.get(reco::btau::tau1_trackEtaRel_2);
0097   inputs["tau_vertexMass_0"] = vars.get(reco::btau::tau1_vertexMass);
0098   inputs["tau_vertexEnergyRatio_0"] = vars.get(reco::btau::tau1_vertexEnergyRatio);
0099   inputs["tau_vertexDeltaR_0"] = vars.get(reco::btau::tau1_vertexDeltaR);
0100   inputs["tau_flightDistance2dSig_0"] = vars.get(reco::btau::tau1_flightDistance2dSig);
0101   inputs["tau_vertexMass_1"] = vars.get(reco::btau::tau2_vertexMass);
0102   inputs["tau_vertexEnergyRatio_1"] = vars.get(reco::btau::tau2_vertexEnergyRatio);
0103   inputs["tau_flightDistance2dSig_1"] = vars.get(reco::btau::tau2_flightDistance2dSig);
0104   inputs["jetNTracks"] = vars.get(reco::btau::jetNTracks);
0105   inputs["nSV"] = vars.get(reco::btau::jetNSecondaryVertices);
0106 
0107   // evaluate the MVA
0108   value = mvaID->evaluate(inputs);
0109 
0110   // return the final discriminator value
0111   return value;
0112 }
0113 
0114 void CandidateBoostedDoubleSecondaryVertexComputer::fillPSetDescription(edm::ParameterSetDescription& desc) {
0115   desc.add<bool>("useCondDB", false);
0116   desc.add<std::string>("gbrForestLabel", "");
0117   desc.add<edm::FileInPath>("weightFile", edm::FileInPath());
0118   desc.add<bool>("useGBRForest", false);
0119   desc.add<bool>("useAdaBoost", false);
0120 }