Belle II Software  release-08-01-10
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 (ROOT.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  found_obj = False
171  if not obj:
172  subdir = rootfile.GetDirectory(args[0])
173  if subdir:
174  obj = subdir.Get(name)
175  if obj:
176  found_obj = True
177  validation_metadata_set(obj, *(args[1:]), **argk)
178  subdir.WriteObject(obj, name, "Overwrite")
179  if not found_obj:
180  raise RuntimeError(
181  f"Cannot find object named {name} in {rootfile.GetName()}"
182  )
183  else:
184  validation_metadata_set(obj, *args, **argk)
185  # scope guard to avoid side effects by changing the global gDirectory
186  # in modules ...
187  # noinspection PyUnusedLocal
188  directory_guard = ROOT.TDirectory.TContext(rootfile) # noqa
189  obj.Write("", ROOT.TObject.kOverwrite)
190  if opened:
191  rootfile.Close()
192 
193 
194 class ValidationMetadataSetter(basf2.Module):
195  """
196  Simple module to set the valdiation metadata for a given list of objects
197  automatically at the end of event processing
198 
199  Just add this module **before** any
200  VariablesToNtuple/VariablesToHistogram modules and it will set the
201  correct validation metadata at the end of processing
202 
203  Warning:
204  The module needs to be before the modules creating the objects
205  as terminate() functions are executed in reverse order from last to
206  first module. If this module is after the creation modules the metadata
207  might not be set correctly
208  """
209 
210  def __init__(
211  self,
212  variables: List[Tuple[str]],
213  rootfile: Union[str, pathlib.PurePath],
214  description="",
215  ):
216  """
217  Initialize ValidationMetadataSetter
218 
219  Arguments:
220  variables (list(tuple(str))): List of objects to set the metadata
221  for. Each entry should be the name of an object followed by the
222  metadata values which will be forwarded to
223  `validation_metadata_set`:
224  ``(name, title, contact, description, check, xlabel, ylabel,
225  metaoptions)``
226  where ``xlabel``, ``ylabel`` and ``metaoptions`` are optional
227  rootfile (str or pathlib.PurePath): The name of the ROOT file where
228  the objects can be found
229  description (str): Common description/information of/about all plots
230  in this ROOT file (will be displayed above the plots)
231  """
232  super().__init__()
233 
234  self._variables = variables
235  if isinstance(rootfile, pathlib.PurePath):
236  rootfile = str(rootfile)
237 
238  self._rootfile: str = rootfile
239 
241  self._description = description
242 
244  self._tfile: ROOT.TFile = None
245 
246  def initialize(self):
247  """Make sure we keep the file open"""
248  self._tfile = Belle2.RootFileCreationManager.getInstance().getFile(
249  self._rootfile
250  )
251 
252  def terminate(self):
253  """And update the metadata at the end"""
254  for name, *metadata in self._variables:
256  validation_metadata_update(self._tfile, name, *metadata)
257  if self._description:
258  file_description_set(self._tfile, self._description)
259  del self._tfile
260 
261 
262 def create_validation_histograms(
263  path: basf2.Path,
264  rootfile: Union[str, pathlib.PurePath],
265  particlelist: str,
266  variables_1d: Optional[List[Tuple]] = None,
267  variables_2d: Optional[List[Tuple]] = None,
268  description="",
269 ) -> None:
270  """
271  Create histograms for all the variables and also label them to be useful
272  in validation plots in one go. This is similar to the
273  `modularAnalysis.variablesToHistogram` function but also sets the
274  metadata correctly to be used by the validation
275 
276  Arguments:
277  path (basf2.Path): Path where to put the modules
278  rootfile (str or pathlib.PurePath): Name of the output root file
279  particlelist (str): Name of the particle list, can be empty for event
280  dependent variables
281  variables_1d: List of 1D histogram definitions of the form
282  ``var, bins, min, max, title, contact, description, check_for
283  [, xlabel [, ylabel [, metaoptions]]]``
284  variables_2d: List of 2D histogram definitions of the form
285  ``var1, bins1, min1, max1, var2, bins2, min2, max2, title, contact,
286  description, check_for [, xlabel [, ylabel [, metaoptions]]]``
287  description: Common description/information of/about all plots in this
288  ROOT file (will be displayed above the plots)
289 
290  .. warning::
291 
292  Sadly, there are two different ways to specify latex formulas.
293 
294  1. The ROOT-style using ``#``:
295  ``"Einstein-Pythagoras a^{2} + b^{2} = #frac{E}{m}"``
296 
297  This style should be used for histogram title and labels, that is the
298  ``title``, ``xlabel`` and ``ylabel`` arguments
299 
300  2. The normal Latex style (with escaped backslashes or raw
301  string literals):
302  ``"Einstein-Pythagoras $a^2 + b^2 = \\\\frac{E}{m}$"``.
303 
304  This style should be used for all other fields like ``description``,
305  ``check_for``
306 
307  You can use the normal Latex style also for histogram title and descriptions
308  but the PDF version of the plots will still be buggy and not show all
309  formulas correctly.
310  """
311  if isinstance(rootfile, pathlib.PurePath):
312  rootfile = str(rootfile)
313  histograms_1d = []
314  histograms_2d = []
315  metadata = []
316  if variables_1d is not None:
317  for var, nbins, vmin, vmax, *data in variables_1d:
318  histograms_1d.append((var, nbins, vmin, vmax))
319  metadata.append([var] + data)
320 
321  if variables_2d is not None:
322  for row in variables_2d:
323  var1 = row[0]
324  var2 = row[4]
325  histograms_2d.append(row[:8])
326  metadata.append([var1 + var2] + list(row[8:]))
327 
328  path.add_module(
329  ValidationMetadataSetter(metadata, rootfile, description=description)
330  )
331  path.add_module(
332  "VariablesToHistogram",
333  particleList=particlelist,
334  variables=histograms_1d,
335  variables_2d=histograms_2d,
336  fileName=rootfile,
337  )
338 
339 
340 # End suppression of doxygen checks
341 # @endcond
static std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.
static RootFileCreationManager & getInstance()
Interface for the FileManager.