Belle II Software  release-06-02-00
metadata.py
1 #!/usr/bin/env python3
2 
3 
10 
11 """Small module containing helper functions to set the metadata on objects
12 created for the validation correctly """
13 
14 # std
15 from typing import Optional, Union, List, Tuple
16 import pathlib
17 
18 import basf2
19 import ROOT
20 from ROOT import Belle2
21 
22 # circumvent BII-1264
23 ROOT.gInterpreter.Declare("#include <framework/utilities/MakeROOTCompatible.h>")
24 
25 # Unfortunately doxygen doesn't recognize our docstrings in this file :/
26 # @cond SUPPRESS_DOXYGEN
27 
28 
29 def file_description_set(
30  rootfile: Union[ROOT.TFile, str, pathlib.PurePath], description: str
31 ) -> None:
32  """
33  Add file description validation metdata to a ROOT file.
34 
35  Args:
36  rootfile (TFile, str or pathlib.PurePath): Name of the root file
37  to open or an already open TFile instance
38  description (str): Common description/information of/about all plots
39  in this ROOT file (will be displayed above the plots)
40 
41  Returns:
42  None
43  """
44  opened = False
45  if not isinstance(rootfile, ROOT.TFile):
46  if isinstance(rootfile, pathlib.PurePath):
47  rootfile = str(rootfile)
48  rootfile = ROOT.TFile(rootfile, "UPDATE")
49  opened = True
50  if not rootfile.IsOpen() or not rootfile.IsWritable():
51  raise RuntimeError(
52  f"ROOT file {rootfile.GetName()} is not open for writing"
53  )
54  # scope guard to avoid side effects by changing the global gDirectory
55  # in modules ...
56  # noinspection PyUnusedLocal
57  directory_guard = ROOT.TDirectory.TContext(rootfile) # noqa
58  desc = ROOT.TNamed("Description", description)
59  desc.Write()
60  if opened:
61  rootfile.Close()
62 
63 
64 def validation_metadata_set(
65  obj: ROOT.TObject,
66  title: str,
67  contact: str,
68  description: str,
69  check: str,
70  xlabel: Optional[str] = None,
71  ylabel: Optional[str] = None,
72  metaoptions="",
73 ) -> None:
74  """
75  Set the validation metadata for a given object by setting the necessary
76  values. This function can be used on any object supported by the
77  Validation (histograms, profiles, ntuples)
78 
79  Arguments:
80  obj: Instance of the object which should get the metadata
81  title (str): Title to use for the object
82  contact (str): Contact person, usually in the form "Name <email>"
83  description (str): Text description what can be seen in the plot
84  check (str): Text description what to look for in the validation for
85  shifters to easily see if the distribution looks ok
86  xlabel (str): If given try to set this as the label for the x axis
87  ylabel (str): If given try to set this as the label for the y axis
88  metaoptions (str): Metaoptions (additional options to influence the
89  comparison between revisions, styling of the plot, etc.)
90 
91  .. warning::
92 
93  Different ways to specify LaTeX for different arguments:
94  see `create_validation_histograms`
95  """
96  obj.SetTitle(title)
97  # For ntuples we add the metadata as aliases ...
98  try:
99  obj.SetAlias("Contact", contact)
100  obj.SetAlias("Description", description)
101  obj.SetAlias("Check", check)
102  obj.SetAlias("MetaOptions", metaoptions)
103  except AttributeError:
104  pass
105  # for TH*, TProfile we add it to the list of functions
106  try:
107  function_list = obj.GetListOfFunctions()
108  function_list.Add(ROOT.TNamed("Contact", contact))
109  function_list.Add(ROOT.TNamed("Description", description))
110  function_list.Add(ROOT.TNamed("Check", check))
111  function_list.Add(ROOT.TNamed("MetaOptions", metaoptions))
112  except AttributeError:
113  pass
114 
115  # And maybe try to set a label on the axis
116  if xlabel is not None:
117  try:
118  obj.GetXaxis().SetTitle(xlabel)
119  except AttributeError:
120  pass
121 
122  if ylabel is not None:
123  try:
124  obj.GetYaxis().SetTitle(ylabel)
125  except AttributeError:
126  pass
127 
128 
129 # noinspection PyIncorrectDocstring
130 def validation_metadata_update(
131  rootfile: Union[str, ROOT.TFile, pathlib.PurePath], name: str, *args, **argk
132 ) -> None:
133  """
134  This is a convenience helper for `validation_metadata_set` in case the
135  objects have already been saved in a ROOT file before: It will open the
136  file (or use an existing TFile), extract the object, add the metadata and
137  save the new version to the file
138 
139  Arguments:
140  rootfile (str, ROOT.TFile or pathlib.PurePath): Name of the root file
141  to open or an already open TFile instance
142  name (str): Name of the object in the file
143  title (str): Title to use for the object
144  contact (str): Contact person, usually in the form "Name <email>"
145  description (str): Text description what can be seen in the plot
146  check (str): Text description what to look for in the validation for
147  shifters to easily see if the distribution looks ok
148  xlabel (str): If given try to set this as the label for the x axis
149  ylabel (str): If given try to set this as the label for the y axis
150  metaoptions (str): Metaoptions (additional options to influence the
151  comparison between revisions, styling of the plot, etc.)
152 
153  .. warning::
154 
155  Different ways to specify LaTeX for different arguments:
156  see `create_validation_histograms`
157  """
158 
159  opened = False
160  if not isinstance(rootfile, ROOT.TFile):
161  if isinstance(rootfile, pathlib.PurePath):
162  rootfile = str(rootfile)
163  rootfile = ROOT.TFile(rootfile, "UPDATE")
164  opened = True
165  if not rootfile.IsOpen() or not rootfile.IsWritable():
166  raise RuntimeError(
167  f"ROOT file {rootfile.GetName()} is not open for writing"
168  )
169  obj = rootfile.Get(name)
170  if not obj:
171  raise RuntimeError(
172  f"Cannot find object named {name} in {rootfile.GetName()}"
173  )
174  validation_metadata_set(obj, *args, **argk)
175  # scope guard to avoid side effects by changing the global gDirectory
176  # in modules ...
177  # noinspection PyUnusedLocal
178  directory_guard = ROOT.TDirectory.TContext(rootfile) # noqa
179  obj.Write("", ROOT.TObject.kOverwrite)
180  if opened:
181  rootfile.Close()
182 
183 
184 class ValidationMetadataSetter(basf2.Module):
185  """
186  Simple module to set the valdiation metadata for a given list of objects
187  automatically at the end of event processing
188 
189  Just add this module **before** any
190  VariablesToNtuple/VariablesToHistogram modules and it will set the
191  correct validation metadata at the end of processing
192 
193  Warning:
194  The module needs to be before the modules creating the objects
195  as terminate() functions are executed in reverse order from last to
196  first module. If this module is after the creation modules the metadata
197  might not be set correctly
198  """
199 
200  def __init__(
201  self,
202  variables: List[Tuple[str]],
203  rootfile: Union[str, pathlib.PurePath],
204  description="",
205  ):
206  """
207  Initialize ValidationMetadataSetter
208 
209  Arguments:
210  variables (list(tuple(str))): List of objects to set the metadata
211  for. Each entry should be the name of an object followed by the
212  metadata values which will be forwarded to
213  `validation_metadata_set`:
214  ``(name, title, contact, description, check, xlabel, ylabel,
215  metaoptions)``
216  where ``xlabel``, ``ylabel`` and ``metaoptions`` are optional
217  rootfile (str or pathlib.PurePath): The name of the ROOT file where
218  the objects can be found
219  description (str): Common description/information of/about all plots
220  in this ROOT file (will be displayed above the plots)
221  """
222  super().__init__()
223 
224  self._variables = variables
225  if isinstance(rootfile, pathlib.PurePath):
226  rootfile = str(rootfile)
227 
228  self._rootfile: str = rootfile
229 
231  self._description = description
232 
234  self._tfile: ROOT.TFile = None
235 
236  def initialize(self):
237  """Make sure we keep the file open"""
238  self._tfile = Belle2.RootFileCreationManager.getInstance().getFile(
239  self._rootfile
240  )
241 
242  def terminate(self):
243  """And update the metadata at the end"""
244  for name, *metadata in self._variables:
245  name = Belle2.makeROOTCompatible(name)
246  validation_metadata_update(self._tfile, name, *metadata)
247  if self._description:
248  file_description_set(self._tfile, self._description)
249  del self._tfile
250 
251 
252 def create_validation_histograms(
253  path: basf2.Path,
254  rootfile: Union[str, pathlib.PurePath],
255  particlelist: str,
256  variables_1d: Optional[List[Tuple]] = None,
257  variables_2d: Optional[List[Tuple]] = None,
258  description="",
259 ) -> None:
260  """
261  Create histograms for all the variables and also label them to be useful
262  in validation plots in one go. This is similar to the
263  `modularAnalysis.variablesToHistogram` function but also sets the
264  metadata correctly to be used by the validation
265 
266  Arguments:
267  path (basf2.Path): Path where to put the modules
268  rootfile (str or pathlib.PurePath): Name of the output root file
269  particlelist (str): Name of the particle list, can be empty for event
270  dependent variables
271  variables_1d: List of 1D histogram definitions of the form
272  ``var, bins, min, max, title, contact, description, check_for
273  [, xlabel [, ylabel [, metaoptions]]]``
274  variables_2d: List of 2D histogram definitions of the form
275  ``var1, bins1, min1, max1, var2, bins2, min2, max2, title, contact,
276  description, check_for [, xlabel [, ylabel [, metaoptions]]]``
277  description: Common description/information of/about all plots in this
278  ROOT file (will be displayed above the plots)
279 
280  .. warning::
281 
282  Sadly, there are two different ways to specify latex formulas.
283 
284  1. The ROOT-style using ``#``:
285  ``"Einstein-Pythagoras a^{2} + b^{2} = #frac{E}{m}"``
286 
287  This style should be used for histogram title and labels, that is the
288  ``title``, ``xlabel`` and ``ylabel`` arguments
289 
290  2. The normal Latex style (with escaped backslashes or raw
291  string literals):
292  ``"Einstein-Pythagoras $a^2 + b^2 = \\\\frac{E}{m}$"``.
293 
294  This style should be used for all other fields like ``description``,
295  ``check_for``
296 
297  You can use the normal Latex style also for histogram title and descriptions
298  but the PDF version of the plots will still be buggy and not show all
299  formulas correctly.
300  """
301  if isinstance(rootfile, pathlib.PurePath):
302  rootfile = str(rootfile)
303  histograms_1d = []
304  histograms_2d = []
305  metadata = []
306  if variables_1d is not None:
307  for var, nbins, vmin, vmax, *data in variables_1d:
308  histograms_1d.append((var, nbins, vmin, vmax))
309  metadata.append([var] + data)
310 
311  if variables_2d is not None:
312  for row in variables_2d:
313  var1 = row[0]
314  var2 = row[4]
315  histograms_2d.append(row[:8])
316  metadata.append([var1 + var2] + list(row[8:]))
317 
318  path.add_module(
319  ValidationMetadataSetter(metadata, rootfile, description=description)
320  )
321  path.add_module(
322  "VariablesToHistogram",
323  particleList=particlelist,
324  variables=histograms_1d,
325  variables_2d=histograms_2d,
326  fileName=rootfile,
327  )
328 
329 
330 # End suppression of doxygen checks
331 # @endcond
static RootFileCreationManager & getInstance()
Interface for the FileManager.
std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.