Belle II Software  release-08-00-10
ECLChargedPIDMVAModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 // THIS MODULE
10 #include <ecl/modules/eclChargedPIDMVA/ECLChargedPIDMVAModule.h>
11 
12 // FRAMEWORK
13 #include <framework/gearbox/Const.h>
14 #include <framework/gearbox/Unit.h>
15 #include <framework/logging/Logger.h>
16 
17 // ANALYSIS
18 #include <analysis/VariableManager/Manager.h>
19 #include <analysis/dataobjects/Particle.h>
20 
21 // MVA
22 #include <mva/interface/Interface.h>
23 #include <mva/methods/TMVA.h>
24 
25 // MDST
26 #include <mdst/dataobjects/Track.h>
27 
28 // C++
29 #include <limits>
30 #include <string>
31 #include <math.h>
32 #include "TMath.h"
33 
34 using namespace Belle2;
36 
37 REG_MODULE(ECLChargedPIDMVA);
38 
40 {
41  setDescription("This module implements charged particle identification using ECL-related observables via a multiclass MVA. For each track a set of relevant ECL variables(shower shape, PSD etc.) are fed to the MVA method which is stored in a conditions database payload. The MVA output variables are then used to construct a likelihood from pdfs also stored in the payload. The likelihood is then stored in the ECLPidLikelihood object.");
43  addParam("payloadName",
45  "The name of the database payload object with the MVA weights.",
46  std::string("ECLChargedPidMVAWeights"));
47 }
48 
50 
52 {
53  m_eventMetaData.isRequired();
54 
55  m_mvaWeights = std::make_unique<DBObjPtr<ECLChargedPIDMVAWeights>>(m_payload_name);
56 
57  m_eclPidLikelihoods.registerInDataStore();
58  m_tracks.registerRelationTo(m_eclPidLikelihoods);
59 }
60 
62 {
63  if (!m_mvaWeights) {B2FATAL("No ECLChargedPidMVAWeights payload found in database!");}
64 }
65 
67 {
68  B2DEBUG(22, "Run: " << m_eventMetaData->getRun() <<
69  ". Loading supported MVA interfaces for multi-class charged particle identification.");
70 
72  m_binningVariables = manager.getVariables((*m_mvaWeights.get())->getBinningVariables());
73 
75  auto supported_interfaces = MVA::AbstractInterface::getSupportedInterfaces();
76 
77  for (const auto& iterator : * (*m_mvaWeights.get())->getPhasespaceCategories()) {
78 
79  B2DEBUG(22, "\t\tweightfile[" << iterator.first << "]");
80 
81  // De-serialize the string into an MVA::Weightfile object.
82  std::stringstream ss(iterator.second.getSerialisedWeight());
83  auto weightfile = MVA::Weightfile::loadFromStream(ss);
84 
85  MVA::GeneralOptions general_options;
86  weightfile.getOptions(general_options);
87 
88  // Store an MVA::Expert object.
89  m_experts[iterator.first] = supported_interfaces[general_options.m_method]->getExpert();
90  m_experts.at(iterator.first)->load(weightfile);
91 
92  B2DEBUG(22, "\t\tweightfile at " << iterator.first << " successfully initialised.");
93 
94  // Get the full version of the variable names.
95  // These are stored in the xml as an additional vector.
96  const std::string identifier = std::string("de_aliased_clf_vars");
97  auto clf_vars = weightfile.getVector<std::string>(identifier);
98  m_variables[iterator.first] = manager.getVariables(clf_vars);
99 
100  B2DEBUG(22, "\t\t Weightfile at " << iterator.first << " has " << m_variables[iterator.first].size() << " features.");
101 
102  std::vector<float> features(clf_vars.size(), 0.0);
103  std::vector<float> spectators;
104 
105  m_datasets[iterator.first] = std::make_unique<MVA::SingleDataset>(general_options, features, 1.0, spectators);
106  }
107 }
108 
110 {
111  (*m_mvaWeights.get()).addCallback([this]() { checkDBPayloads();});
112  checkDBPayloads();
113  (*m_mvaWeights.get()).addCallback([this]() { initializeMVA(); });
114  initializeMVA();
115 }
116 
118 {
119  for (const auto& track : m_tracks) {
120 
121  // Load the pion fit hypothesis or closest in mass.
122  const TrackFitResult* fitRes = track.getTrackFitResultWithClosestMass(Const::pion);
123  if (fitRes == nullptr) continue;
124 
125  // If a track omega value is very small trackFitResult might return a charge of 0.
126  // The particle is then not able to be assigned a pdg code causing a crash.
127  // Skip such tracks.
128  const int charge = fitRes->getChargeSign();
129  if (charge == 0) continue;
130 
131  // Create a particle object so that we can use the variable manager.
132  const Particle particle = Particle(&track, Const::pion);
133 
134  // Get global bin index for track corresponding to N dimensional binning.
135  std::vector<float> binningVariableValues(m_binningVariables.size());
136  for (unsigned int ivar(0); ivar < binningVariableValues.size(); ivar++) {
137  auto varobj = m_binningVariables.at(ivar);
138  binningVariableValues.at(ivar) = evaluateVariable(varobj, &particle);
139  }
140  const int linearCategoryIndex = (*m_mvaWeights.get())->getLinearisedCategoryIndex(binningVariableValues);
141 
142  // Require we cover the phasespace.
143  // Alternatively could take closest covered region but behaviour will not be well understood.
144  // After this check the linearCategoryIndex is guaranteed to be positive.
145  if (!(*m_mvaWeights.get())->isPhasespaceCovered(linearCategoryIndex)) continue;
146 
147  // Get the phasespaceCategory
148  const auto phasespaceCategory = (*m_mvaWeights.get())->getPhasespaceCategory(linearCategoryIndex);
149  const auto transformMode = phasespaceCategory->getTransformMode();
150 
151  // Fill the feature vectors
152  unsigned int nvars = m_variables.at(linearCategoryIndex).size();
153  for (unsigned int ivar(0); ivar < nvars; ++ivar) {
154  auto varobj = m_variables.at(linearCategoryIndex).at(ivar);
155  m_datasets.at(linearCategoryIndex)->m_input[ivar] = evaluateVariable(varobj, &particle);
156  }
157 
158  // Get the mva response to be converted to a likelihood
159  std::vector<float> scores = m_experts.at(linearCategoryIndex)->applyMulticlass(*m_datasets.at(linearCategoryIndex))[0];
160 
161  // Log transform the scores
162  // This turns the MVA output which, if using the nominal TMVA MulticlassBDT,
163  // is heavily peaked at 0 and 1 into a smooth curve.
164  // We can then evaluate the likelihoods from this curve directly or further transform the responses.
165  for (unsigned int iResponse = 0; iResponse < scores.size(); iResponse++) {
166  if ((transformMode != transformModes::c_DirectMVAResponse) and
167  (transformMode != transformModes::c_LogMVAResponse)) {
168  scores[iResponse] = logTransformation(scores[iResponse],
169  phasespaceCategory->getLogTransformOffset(),
170  phasespaceCategory->getMaxPossibleResponseValue());
171  }
172  }
173 
174  B2DEBUG(22, "Scores: ");
175  for (unsigned int iResponse = 0; iResponse < scores.size(); iResponse++) {
176  B2DEBUG(22, "\t\t " << iResponse << " : " << scores[iResponse]);
177  }
178 
179  float logLikelihoods[Const::ChargedStable::c_SetSize];
180 
181  // Order of loop is defined in UnitConst.cc: e, mu, pi, K, p, d
182  for (const auto& hypo : Const::chargedStableSet) {
183  unsigned int hypo_idx = hypo.getIndex();
184  auto absPdgId = abs(hypo.getPDGCode());
185 
186  // Copy the scores so they aren't transformed in place
187  std::vector<float> transformed_scores;
188  transformed_scores = scores;
189 
190  // Perform extra transformations if they are booked
191  if ((transformMode == transformModes::c_GaussianTransform) or
192  (transformMode == transformModes::c_DecorrelationTransform)) {
193 
194  // Gaussian transform
195  for (unsigned int iResponse = 0; iResponse < scores.size(); iResponse++) {
196  transformed_scores[iResponse] = gaussTransformation(scores[iResponse], phasespaceCategory->getCDF(absPdgId, iResponse));
197  }
198  if (transformMode == transformModes::c_DecorrelationTransform) {
199  transformed_scores = decorrTransformation(transformed_scores, phasespaceCategory->getDecorrelationMatrix(absPdgId));
200  }
201  }
202 
203  // Get the pdf values for each response value
204  float logL = 0.0;
205  if ((transformMode == transformModes::c_DirectMVAResponse) or
206  (transformMode == transformModes::c_LogMVAResponse)) {
207 
208  logLikelihoods[hypo_idx] = scores[phasespaceCategory->getMVAIndexForHypothesis(absPdgId)];
209  if (transformMode == transformModes::c_LogMVAResponse) {
210  logLikelihoods[hypo_idx] = (std::isnormal(logLikelihoods[hypo_idx])
211  && logLikelihoods[hypo_idx] > 0) ? std::log(logLikelihoods[hypo_idx]) : c_dummyLogL;
212  }
213  logLikelihoods[hypo_idx] = logLikelihoods[hypo_idx] / phasespaceCategory->getTemperature();
214  continue;
215  }
216  B2DEBUG(22, "MVA Index for hypo " << absPdgId << " : " << phasespaceCategory->getMVAIndexForHypothesis(absPdgId));
217 
218  for (unsigned int iResponse = 0; iResponse < transformed_scores.size(); iResponse++) {
219  if ((transformMode == transformModes::c_LogTransformSingle) and
220  (phasespaceCategory->getMVAIndexForHypothesis(absPdgId) != iResponse)) {continue;}
221 
222  double xmin, xmax;
223  phasespaceCategory->getPDF(iResponse, absPdgId)->GetRange(xmin, xmax);
224  // Kick the response back into the range of the PDF (alternatively consider logL = std::numeric_limits<float>::min() if outside range)
225  float transformed_score_copy = transformed_scores[iResponse];
226  if (transformed_score_copy < xmin) transformed_score_copy = xmin + 1e-5;
227  if (transformed_score_copy > xmax) transformed_score_copy = xmax - 1e-5;
228 
229  float pdfval = phasespaceCategory->getPDF(iResponse, absPdgId)->Eval(transformed_score_copy);
230  // Don't take a log of inf or 0
231  logL += (std::isnormal(pdfval) && pdfval > 0) ? std::log(pdfval) : c_dummyLogL;
232  B2DEBUG(22, "MVA response index, MVA score, logL: " << iResponse << " " << transformed_score_copy << " " << logL << " \n");
233  }
234  logLikelihoods[hypo_idx] = logL;
235  } // hypo loop
236 
237  B2DEBUG(22, "logL: ");
238  for (unsigned int iLogL = 0; iLogL < Const::ChargedStable::c_SetSize; iLogL++) {
239  B2DEBUG(22, "\t\t " << iLogL << " : " << logLikelihoods[iLogL]);
240  }
241 
242  const auto eclPidLikelihood = m_eclPidLikelihoods.appendNew(logLikelihoods);
243 
244  track.addRelationTo(eclPidLikelihood);
245  } // tracks
246 }
247 
248 float ECLChargedPIDMVAModule::logTransformation(const float value, const float offset, const float max) const
249 {
250  return std::log((value + offset) / (max + offset - value));
251 }
252 
253 float ECLChargedPIDMVAModule::gaussTransformation(const float value, const TH1F* cdf) const
254 {
255  int binIndex = cdf->GetXaxis()->FindBin(value);
256 
257  double dx = cdf->GetBinWidth(binIndex);
258  double cumulant = cdf->GetBinContent(binIndex);
259  double dy;
260 
261  if (binIndex == cdf->GetNbinsX()) {
262  dy = 1 - cdf->GetBinContent(binIndex);
263  } else if (binIndex == 0) {
264  dy = cdf->GetBinContent(binIndex);
265  } else {
266  dy = cdf->GetBinContent(binIndex + 1) - cdf->GetBinContent(binIndex);
267  }
268 
269  cumulant = cumulant + (value - cdf->GetBinLowEdge(binIndex)) * dy / dx;
270 
271  cumulant = std::min(cumulant, 1.0 - 1e-9);
272  cumulant = std::max(cumulant, 1e-9);
273 
274  double maxErfInvArgRange = 0.99999999;
275  double arg = 2.0 * cumulant - 1.0;
276 
277  arg = std::min(maxErfInvArgRange, arg);
278  arg = std::max(-maxErfInvArgRange, arg);
279 
280  return c_sqrt2 * TMath::ErfInverse(arg);
281 }
282 
283 
284 std::vector<float> ECLChargedPIDMVAModule::decorrTransformation(const std::vector<float> scores,
285  const std::vector<float>* decorrelationMatrix) const
286 {
287  unsigned int n_scores = scores.size();
288  std::vector<float> decor_scores(n_scores, 0.0);
289 
290  for (unsigned int irow = 0; irow < n_scores; irow++) {
291  unsigned int offset = irow * n_scores;
292  for (unsigned int icol = 0; icol < n_scores; icol++) {
293  decor_scores[irow] += scores[irow] * decorrelationMatrix->at(icol + offset);
294  }
295  }
296  return decor_scores;
297 }
298 
300 {
301  float val = std::numeric_limits<float>::quiet_NaN();
302  if (std::holds_alternative<double>(varobj->function(particle))) {
303  val = (float)std::get<double>(varobj->function(particle));
304  } else if (std::holds_alternative<int>(varobj->function(particle))) {
305  val = (float)std::get<int>(varobj->function(particle));
306  } else if (std::holds_alternative<bool>(varobj->function(particle))) {
307  val = (float)std::get<bool>(varobj->function(particle));
308  } else {
309  B2ERROR("Variable '" << varobj->name << "' has wrong data type! It must be one of double, integer, or bool.");
310  }
311  if (std::isnan(val)) {
312  //NaNs cause crashed if using TMVA MulticlassBDT.
313  //In case other MVA methods are used, NaNs likely should be masked properly also.
314  B2FATAL("Variable '" << varobj->name <<
315  "' has NaN entries. Variable definitions in the MVA payload should include 'ifNaNGiveX(X, DEFAULT)'" <<
316  "with a proper default value if it is possible for the variable to evaluate as NaN.");
317  }
318 
319  return val;
320 }
static const unsigned int c_SetSize
Number of elements (for use in array bounds etc.)
Definition: Const.h:606
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:609
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
std::unordered_map< unsigned int, std::unique_ptr< MVA::Expert > > m_experts
unordered map of MVA experts.
std::vector< const Variable::Manager::Var * > m_binningVariables
Vector of variables used to define the regions in which the MVAs are trained.
void checkDBPayloads()
Check payloads for consistency.
std::vector< float > decorrTransformation(const std::vector< float > scores, const std::vector< float > *decorrelationMatrix) const
decorrelation transformation.
virtual void initialize() override
Use this to initialize resources or memory your module needs.
virtual void event() override
Called once for each event.
std::unordered_map< unsigned int, std::unique_ptr< MVA::SingleDataset > > m_datasets
MVA dataset to be passed to the expert.
static constexpr double c_sqrt2
Definition of sqrt(2)
StoreObjPtr< EventMetaData > m_eventMetaData
The event information.
void initializeMVA()
Initialise the multiclass MVAs.
virtual ~ECLChargedPIDMVAModule()
Destructor.
float evaluateVariable(const Variable::Manager::Var *varobj, const Particle *particle)
evaluate the variable for the particle and convert the return value float.
float logTransformation(const float value, const float offset, const float max) const
log transformation.
static constexpr float c_dummyLogL
Dummy value of log Likelihood for a particle hypothesis.
virtual void beginRun() override
Called once before a new run begins.
std::unordered_map< unsigned int, std::vector< const Variable::Manager::Var * > > m_variables
Map of vectors containing the variables objects to be fed to the MVA.
StoreArray< Track > m_tracks
array of track objects.
float gaussTransformation(const float value, const TH1F *cdf) const
gaussian transformation.
std::unique_ptr< DBObjPtr< ECLChargedPIDMVAWeights > > m_mvaWeights
Interface to DB payloads containing MVA weightfiles and pdfs.
StoreArray< ECLPidLikelihood > m_eclPidLikelihoods
Array of ECLPidLikelihood objects.
std::string m_payload_name
The name of the database payload object with the MVA weights.
MVAResponseTransformMode
Enum of implemented transformations which can be applied to the MVA response.
static std::map< std::string, AbstractInterface * > getSupportedInterfaces()
Returns interfaces supported by the MVA Interface.
Definition: Interface.h:53
static void initSupportedInterfaces()
Static function which initliazes all supported interfaces, has to be called once before getSupportedI...
Definition: Interface.cc:45
General options which are shared by all MVA trainings.
Definition: Options.h:62
static Weightfile loadFromStream(std::istream &stream)
Static function which deserializes a Weightfile from a stream.
Definition: Weightfile.cc:251
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Class to store reconstructed particles.
Definition: Particle.h:75
Values of the result of a track fit with a given particle hypothesis.
Global list of available variables.
Definition: Manager.h:101
static Manager & Instance()
get singleton instance.
Definition: Manager.cc:25
REG_MODULE(arichBtest)
Register the Module.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
Abstract base class for different kinds of events.
std::string name
Unique identifier of the function, used as key.
Definition: Manager.h:130
A variable returning a floating-point value for a given Particle.
Definition: Manager.h:146
FunctionPtr function
Pointer to function.
Definition: Manager.h:147