-
Notifications
You must be signed in to change notification settings - Fork 68
/
TMatrixT.py
88 lines (74 loc) · 2.66 KB
/
TMatrixT.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
# BSD 3-Clause License; see https://github.com/scikit-hep/uproot5/blob/main/LICENSE
"""
This module defines versioned models for ``TLeaf`` and its subclasses.
"""
from __future__ import annotations
import struct
import numpy
import uproot
import uproot._util
import uproot.model
class Model_TMatrixTSym_3c_double_3e__v5(uproot.model.VersionedModel):
"""
A :doc:`uproot.model.VersionedModel` for ``TMatrixTSym<double>`` version 2,
which shows up as version 5 because it's reading the ``TMatrixTBase<double>``
header.
"""
def read_members(self, chunk, cursor, context, file):
if uproot._awkwardforth.get_forth_obj(context) is not None:
raise uproot.interpretation.objects.CannotBeForth()
if self.is_memberwise:
raise NotImplementedError(
f"memberwise serialization of {type(self).__name__}\nin file {self.file.file_path}"
)
self._bases.append(
file.class_named("TObject", 1).read(
chunk,
cursor,
context,
file,
self._file,
self._parent,
concrete=self.concrete,
)
)
(
self._members["fNrows"],
self._members["fNcols"],
self._members["fRowLwb"],
self._members["fColLwb"],
self._members["fNelems"],
self._members["fNrowIndex"],
self._members["fTol"],
) = cursor.fields(chunk, self._format0, context)
num_elements = self.member("fNrows") * (self.member("fNcols") + 1) // 2
self._members["fElements"] = cursor.array(
chunk, num_elements, self._dtype0, context
)
self._num_bytes += self._members["fElements"].nbytes
_format0 = struct.Struct(">iiiiiid")
_format_memberwise0 = struct.Struct(">i")
_format_memberwise1 = struct.Struct(">i")
_format_memberwise2 = struct.Struct(">i")
_format_memberwise3 = struct.Struct(">i")
_format_memberwise4 = struct.Struct(">i")
_format_memberwise5 = struct.Struct(">i")
_format_memberwise6 = struct.Struct(">d")
_dtype0 = numpy.dtype(">f8")
base_names_versions = [("TObject", 1)]
member_names = [
"fNrows",
"fNcols",
"fRowLwb",
"fColLwb",
"fNelems",
"fNrowIndex",
"fTol",
]
class_flags = {}
class Model_TMatrixTSym_3c_double_3e_(uproot.model.DispatchByVersion):
"""
A :doc:`uproot.model.DispatchByVersion` for ``TMatrixTSym<double>``.
"""
known_versions = {5: Model_TMatrixTSym_3c_double_3e__v5}
uproot.classes["TMatrixTSym<double>"] = Model_TMatrixTSym_3c_double_3e_