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