Belle II Software  release-06-02-00
PIDLikelihood.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 #include <mdst/dataobjects/PIDLikelihood.h>
10 #include <framework/logging/Logger.h>
11 
12 #include <TROOT.h>
13 #include <TColor.h>
14 #include <TParticlePDG.h>
15 
16 #include <cmath>
17 #include <iostream>
18 #include <iomanip>
19 
20 using namespace std;
21 using namespace Belle2;
22 
23 PIDLikelihood::PIDLikelihood()
24 {
25  for (unsigned short i = 0; i < Const::PIDDetectors::c_size; i++) {
26  for (unsigned int k = 0; k < Const::ChargedStable::c_SetSize; k++) {
27  m_logl[i][k] = 0.0;
28  }
29  }
30 }
31 
32 
33 void PIDLikelihood::setLogLikelihood(Const::EDetector det,
34  const Const::ChargedStable& part,
35  float logl)
36 {
37  int index = Const::PIDDetectors::set().getIndex(det);
38  if (index < 0) {
39  B2ERROR("PIDLikelihood::setLogLikelihood: detector is not a PID device");
40  return;
41  }
42  if (logl != logl or logl == INFINITY) {
43  B2ERROR("PIDLikelihood::setLogLikelihood: log-likelihood for detector " << det << " is " << logl <<
44  " (i.e. +inf or NaN)! Ignoring this value. (" << Const::CDC << "=CDC, " << Const::TOP << "=TOP, " << Const::ARICH << "=ARICH, " <<
45  Const::ECL << "=ECL)");
46 
47  return;
48  }
49  m_detectors += det;
50  m_logl[index][part.getIndex()] = logl;
51 }
52 
53 
54 float PIDLikelihood::getLogL(const Const::ChargedStable& part,
55  Const::PIDDetectorSet set) const
56 {
57  float result = 0;
58  for (unsigned int index = 0; index < Const::PIDDetectorSet::set().size(); ++index) {
59  if (set.contains(Const::PIDDetectorSet::set()[index]))
60  result += m_logl[index][part.getIndex()];
61  }
62  return result;
63 }
64 
65 
66 double PIDLikelihood::getProbability(const Const::ChargedStable& p1,
67  const Const::ChargedStable& p2,
68  double ratio,
69  Const::PIDDetectorSet set) const
70 {
71  if (ratio < 0) {
72  B2ERROR("PIDLikelihood::probability argument 'ratio' is given with negative value");
73  return 0;
74  }
75  if (ratio == 0) return 0;
76 
77  double dlogl = getLogL(p2, set) - getLogL(p1, set);
78  double res;
79  if (dlogl < 0) {
80  double elogl = exp(dlogl);
81  res = ratio / (ratio + elogl);
82  } else {
83  double elogl = exp(-dlogl) * ratio; // to prevent overflow for very large dlogl
84  res = elogl / (1.0 + elogl);
85  }
86  //TODO: only necessary if one wants to use mcprod1405 MC sample. Remove when there's a good replacement.
87  if (std::isfinite(res))
88  return res;
89  return 0;
90 }
91 
92 double PIDLikelihood::getProbability(const Const::ChargedStable& part,
93  const double* fractions,
94  Const::PIDDetectorSet detSet) const
95 {
96  double prob[Const::ChargedStable::c_SetSize];
97  probability(prob, fractions, detSet);
98 
99  int k = part.getIndex();
100  if (k < 0) return 0;
101 
102  return prob[k];
103 
104 }
105 
106 Const::ChargedStable PIDLikelihood::getMostLikely(const double* fractions,
107  Const::PIDDetectorSet detSet) const
108 {
109  const unsigned int n = Const::ChargedStable::c_SetSize;
110  double prob[n];
111  probability(prob, fractions, detSet);
112 
113  int k = 0;
114  double maxProb = prob[k];
115  for (unsigned i = 0; i < n; ++i) {
116  if (prob[i] > maxProb) {maxProb = prob[i]; k = i;}
117  }
118  return Const::chargedStableSet.at(k);
119 
120 }
121 
122 
123 void PIDLikelihood::probability(double probabilities[],
124  const double* fractions,
125  Const::PIDDetectorSet detSet) const
126 {
127  const unsigned int n = Const::ChargedStable::c_SetSize;
128  double frac[n];
129  if (!fractions) {
130  for (unsigned int i = 0; i < n; ++i) frac[i] = 1.0; // normalization not needed
131  fractions = frac;
132  }
133 
134  double logL[n];
135  double logLmax = 0;
136  bool hasMax = false;
137  for (unsigned i = 0; i < n; ++i) {
138  logL[i] = 0;
139  if (fractions[i] > 0) {
140  logL[i] = getLogL(Const::chargedStableSet.at(i), detSet);
141  if (!hasMax || logL[i] > logLmax) {
142  logLmax = logL[i];
143  hasMax = true;
144  }
145  }
146  }
147 
148  double norm = 0;
149  for (unsigned i = 0; i < n; ++i) {
150  probabilities[i] = 0;
151  if (fractions[i] > 0) probabilities[i] = exp(logL[i] - logLmax) * fractions[i];
152  norm += probabilities[i];
153  }
154  if (norm == 0) return;
155 
156  for (unsigned i = 0; i < n; ++i) {
157  probabilities[i] /= norm;
158  }
159 
160 }
161 
162 
163 void PIDLikelihood::printArray() const
164 {
165 
166  string detectorName[Const::PIDDetectors::c_size] =
167  {"SVD", "CDC", "TOP", "ARICH", "ECL", "KLM"};
168  string hline("-------");
169  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++)
170  hline += "-----------";
171  string Hline("=======");
172  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++)
173  Hline += "===========";
174 
175  cout << Hline << endl;
176 
177  cout << "PDGcode";
178  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++)
179  cout << setw(10) << Const::chargedStableSet.at(i).getPDGCode() << " ";
180  cout << endl;
181 
182  cout << Hline << endl;
183 
184  float sum_logl[Const::ChargedStable::c_SetSize];
185  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++) sum_logl[i] = 0;
186 
187  for (unsigned k = 0; k < Const::PIDDetectors::c_size; k++) {
188  cout << setw(7) << detectorName[k];
189  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++) {
190  cout << setw(10) << setprecision(4) << m_logl[k][i] << " ";
191  sum_logl[i] += m_logl[k][i];
192  }
193  cout << endl;
194  }
195 
196  cout << hline << endl;
197 
198  cout << setw(7) << "sum";
199  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++)
200  cout << setw(10) << setprecision(4) << sum_logl[i] << " ";
201  cout << endl;
202 
203  if (isAvailable(Const::SVD) or isAvailable(Const::CDC) or isAvailable(Const::TOP) or
204  isAvailable(Const::ARICH) or isAvailable(Const::ECL) or isAvailable(Const::KLM)) {
205  unsigned k = 0;
206  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++)
207  if (sum_logl[i] > sum_logl[k]) k = i;
208  unsigned pos = 11 * (k + 1) + 3;
209  Hline.replace(pos, 1, "^");
210  }
211  cout << Hline << endl;
212 
213 }
214 
215 std::string PIDLikelihood::getInfoHTML() const
216 {
217  string detectorName[Const::PIDDetectors::c_size] = {"SVD", "CDC", "TOP", "ARICH", "ECL", "KLM"};
218 
219  std::stringstream stream;
220  stream << std::setprecision(4);
221 
222  //colors from event display
223  std::string colour[Const::ChargedStable::c_SetSize];
224  colour[0] = gROOT->GetColor(kAzure)->AsHexString();
225  colour[1] = gROOT->GetColor(kCyan + 1)->AsHexString();
226  colour[2] = gROOT->GetColor(kGray + 1)->AsHexString();
227  colour[3] = gROOT->GetColor(kRed + 1)->AsHexString();
228  colour[4] = gROOT->GetColor(kOrange - 2)->AsHexString();
229  colour[5] = gROOT->GetColor(kMagenta)->AsHexString();
230 
231  stream << "<b>Probabilities</b><br>";
232  stream << "<table border=0 width=100%>";
233  stream << "<tr><td>All</td><td>";
234  stream << "<table cellspacing=1 border=0 width=100%><tr>";
235  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++) {
236  double p = getProbability(Const::chargedStableSet.at(i));
237  stream << "<td bgcolor=\"" << colour[i] << "\" width=" << p * 100 << "%></td>";
238  }
239  stream << "</tr></table></td></tr>";
240 
241  for (unsigned k = 0; k < Const::PIDDetectors::c_size; k++) {
242  auto det = Const::PIDDetectorSet::set()[k];
243  stream << "<tr><td>" << detectorName[k] << "</td><td>";
244  if (!isAvailable(det)) {
245  stream << "</td></tr>";
246  continue;
247  }
248  stream << "<table cellspacing=1 border=0 width=100%><tr>";
249  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++) {
250  double p = getProbability(Const::chargedStableSet.at(i), 0, det);
251  stream << "<td bgcolor=\"" << colour[i] << "\" width=" << p * 100 << "%></td>";
252  }
253  stream << "</tr></table></td></tr>";
254  }
255  stream << "</table>\n";
256 
257  stream << "<b>Log-Likelihoods</b><br>";
258  stream << "<table>";
259  stream << "<tr><th>PID / Detector</th>";
260  for (unsigned k = 0; k < Const::PIDDetectors::c_size; k++)
261  stream << "<th>" << detectorName[k] << "</th>";
262  stream << "</tr>";
263  for (unsigned i = 0; i < Const::ChargedStable::c_SetSize; i++) {
264  stream << "<tr>";
265  stream << "<td bgcolor=\"" << colour[i] << "\">" << Const::chargedStableSet.at(i).getParticlePDG()->GetName() << "</td>";
266  for (unsigned k = 0; k < Const::PIDDetectors::c_size; k++) {
267  stream << "<td>" << m_logl[k][i] << "</td>";
268  }
269  stream << "</tr>";
270  }
271  stream << "</table>";
272  stream << "<br>";
273 
274  return stream.str();
275 }
276 
277 
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:470
bool contains(const DetectorSet &set) const
Check whether this set contains another set.
Definition: Const.h:121
A class for sets of detector IDs whose content is limited to restricted set of valid detector IDs.
Definition: Const.h:177
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition: Const.h:42
Abstract base class for different kinds of events.