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