CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CandidateBoostedDoubleSecondaryVertexComputer.cc
Go to the documentation of this file.
1 #include <memory>
2 
4 
6 
10 
13  if (parameters.getParameter<bool>("useCondDB")) {
14  gbrForest_ = cc.consumes(edm::ESInputTag{"",
15  parameters.existsAs<std::string>("gbrForestLabel")
16  ? parameters.getParameter<std::string>("gbrForestLabel")
17  : ""});
18  }
19 }
20 
22  const edm::ParameterSet& parameters, Tokens tokens)
23  : weightFile_(parameters.existsAs<edm::FileInPath>("weightFile")
24  ? parameters.getParameter<edm::FileInPath>("weightFile")
25  : edm::FileInPath()),
26  useGBRForest_(parameters.existsAs<bool>("useGBRForest") ? parameters.getParameter<bool>("useGBRForest") : false),
27  useAdaBoost_(parameters.existsAs<bool>("useAdaBoost") ? parameters.getParameter<bool>("useAdaBoost") : false),
28  tokens_{tokens} {
29  uses(0, "svTagInfos");
30 
31  mvaID = std::make_unique<TMVAEvaluator>();
32 }
33 
35  // variable names and order need to be the same as in the training
36  std::vector<std::string> variables({"z_ratio",
37  "trackSipdSig_3",
38  "trackSipdSig_2",
39  "trackSipdSig_1",
40  "trackSipdSig_0",
41  "trackSipdSig_1_0",
42  "trackSipdSig_0_0",
43  "trackSipdSig_1_1",
44  "trackSipdSig_0_1",
45  "trackSip2dSigAboveCharm_0",
46  "trackSip2dSigAboveBottom_0",
47  "trackSip2dSigAboveBottom_1",
48  "tau0_trackEtaRel_0",
49  "tau0_trackEtaRel_1",
50  "tau0_trackEtaRel_2",
51  "tau1_trackEtaRel_0",
52  "tau1_trackEtaRel_1",
53  "tau1_trackEtaRel_2",
54  "tau_vertexMass_0",
55  "tau_vertexEnergyRatio_0",
56  "tau_vertexDeltaR_0",
57  "tau_flightDistance2dSig_0",
58  "tau_vertexMass_1",
59  "tau_vertexEnergyRatio_1",
60  "tau_flightDistance2dSig_1",
61  "jetNTracks",
62  "nSV"});
63  // book TMVA readers
64  std::vector<std::string> spectators({"massPruned", "flavour", "nbHadrons", "ptPruned", "etaPruned"});
65 
67  mvaID->initializeGBRForest(&record.get(tokens_.gbrForest_), variables, spectators, useAdaBoost_);
68  } else
69  mvaID->initialize(
70  "Color:Silent:Error", "BDT", weightFile_.fullPath(), variables, spectators, useGBRForest_, useAdaBoost_);
71 }
72 
74  // get the TagInfo
75  const reco::BoostedDoubleSVTagInfo& bdsvTagInfo = tagInfo.get<reco::BoostedDoubleSVTagInfo>(0);
76 
77  // get the TaggingVariables
78  const reco::TaggingVariableList vars = bdsvTagInfo.taggingVariables();
79 
80  // default discriminator value
81  float value = -10.;
82 
83  std::map<std::string, float> inputs;
84  inputs["z_ratio"] = vars.get(reco::btau::z_ratio);
85  inputs["trackSipdSig_3"] = vars.get(reco::btau::trackSip3dSig_3);
86  inputs["trackSipdSig_2"] = vars.get(reco::btau::trackSip3dSig_2);
87  inputs["trackSipdSig_1"] = vars.get(reco::btau::trackSip3dSig_1);
88  inputs["trackSipdSig_0"] = vars.get(reco::btau::trackSip3dSig_0);
89  inputs["trackSipdSig_1_0"] = vars.get(reco::btau::tau2_trackSip3dSig_0);
90  inputs["trackSipdSig_0_0"] = vars.get(reco::btau::tau1_trackSip3dSig_0);
91  inputs["trackSipdSig_1_1"] = vars.get(reco::btau::tau2_trackSip3dSig_1);
92  inputs["trackSipdSig_0_1"] = vars.get(reco::btau::tau1_trackSip3dSig_1);
93  inputs["trackSip2dSigAboveCharm_0"] = vars.get(reco::btau::trackSip2dSigAboveCharm);
94  inputs["trackSip2dSigAboveBottom_0"] = vars.get(reco::btau::trackSip2dSigAboveBottom_0);
95  inputs["trackSip2dSigAboveBottom_1"] = vars.get(reco::btau::trackSip2dSigAboveBottom_1);
96  inputs["tau1_trackEtaRel_0"] = vars.get(reco::btau::tau2_trackEtaRel_0);
97  inputs["tau1_trackEtaRel_1"] = vars.get(reco::btau::tau2_trackEtaRel_1);
98  inputs["tau1_trackEtaRel_2"] = vars.get(reco::btau::tau2_trackEtaRel_2);
99  inputs["tau0_trackEtaRel_0"] = vars.get(reco::btau::tau1_trackEtaRel_0);
100  inputs["tau0_trackEtaRel_1"] = vars.get(reco::btau::tau1_trackEtaRel_1);
101  inputs["tau0_trackEtaRel_2"] = vars.get(reco::btau::tau1_trackEtaRel_2);
102  inputs["tau_vertexMass_0"] = vars.get(reco::btau::tau1_vertexMass);
103  inputs["tau_vertexEnergyRatio_0"] = vars.get(reco::btau::tau1_vertexEnergyRatio);
104  inputs["tau_vertexDeltaR_0"] = vars.get(reco::btau::tau1_vertexDeltaR);
105  inputs["tau_flightDistance2dSig_0"] = vars.get(reco::btau::tau1_flightDistance2dSig);
106  inputs["tau_vertexMass_1"] = vars.get(reco::btau::tau2_vertexMass);
107  inputs["tau_vertexEnergyRatio_1"] = vars.get(reco::btau::tau2_vertexEnergyRatio);
108  inputs["tau_flightDistance2dSig_1"] = vars.get(reco::btau::tau2_flightDistance2dSig);
109  inputs["jetNTracks"] = vars.get(reco::btau::jetNTracks);
110  inputs["nSV"] = vars.get(reco::btau::jetNSecondaryVertices);
111 
112  // evaluate the MVA
113  value = mvaID->evaluate(inputs);
114 
115  // return the final discriminator value
116  return value;
117 }
CandidateBoostedDoubleSecondaryVertexComputer(const edm::ParameterSet &parameters, Tokens tokens)
const T & get(unsigned int index=0) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
ProductT const & get(ESGetToken< ProductT, DepRecordT > const &iToken) const
Tokens(const edm::ParameterSet &parameters, edm::ESConsumesCollector &&cc)
constexpr bool isInitialized() const noexcept
Definition: ESGetToken.h:72
TaggingVariableList taggingVariables(void) const override
returns a description of the extended informations in a TaggingVariableList
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
float discriminator(const TagInfoHelper &tagInfos) const override
std::string fullPath() const
Definition: FileInPath.cc:161
vars
Definition: DeepTauId.cc:164