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 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258
|
"""NeXus file reading and writing."""
# -*- coding: utf-8 -*-
# Copyright 2007-2023 The HyperSpy developers
#
# This file is part of RosettaSciIO.
#
# RosettaSciIO is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# RosettaSciIO is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with RosettaSciIO. If not, see <https://www.gnu.org/licenses/#GPL>.
#
import logging
import os
import dask.array as da
import h5py
import numpy as np
from rsciio._docstrings import (
COMPRESSION_HDF5_DOC,
COMPRESSION_HDF5_NOTES_DOC,
FILENAME_DOC,
LAZY_DOC,
RETURNS_DOC,
SIGNAL_DOC,
)
from rsciio._hierarchical import get_signal_chunks
from rsciio.hspy._api import overwrite_dataset
from rsciio.utils.tools import DTBox
_logger = logging.getLogger(__name__)
def _byte_to_string(value):
"""Decode a byte string.
Parameters
----------
value : byte str
Returns
-------
str
decoded version of input value
"""
try:
text = value.decode("utf-8")
except UnicodeDecodeError:
text = value.decode("latin-1")
return text.replace("\x00", "").rstrip()
def _parse_from_file(value, lazy=False):
"""To convert values from the hdf file to compatible formats.
When reading string arrays we convert or keep string arrays as
byte_strings (some io_plugins only supports byte-strings arrays so this
ensures inter-compatibility across io_plugins)
Arrays of length 1 - return the single value stored.
Large datasets are returned as dask arrays if lazy=True.
Parameters
----------
value : input read from hdf file (array,list,tuple,string,int,float)
lazy : bool {default: False}
The lazy flag is only applied to values of size >=2
Returns
-------
str,int, float, ndarray dask Array
parsed value.
"""
toreturn = value
if isinstance(value, h5py.Dataset):
if value.size < 2:
toreturn = value[...].item()
else:
if lazy:
if value.chunks:
toreturn = da.from_array(value, value.chunks)
else:
chunks = get_signal_chunks(value.shape, value.dtype)
toreturn = da.from_array(value, chunks)
else:
toreturn = np.array(value)
if isinstance(toreturn, np.ndarray) and value.shape == (1,):
toreturn = toreturn[0]
if isinstance(toreturn, bytes):
toreturn = _byte_to_string(toreturn)
if isinstance(toreturn, (np.ndarray)) and toreturn.dtype.char == "U":
toreturn = toreturn.astype("S")
return toreturn
def _parse_to_file(value):
"""Convert to a suitable format for writing to HDF5.
For example unicode values are not compatible with hdf5 so conversion to
byte strings is required.
Parameters
----------
value - input object to write to the hdf file
Returns
-------
parsed value
"""
totest = value
toreturn = totest
if isinstance(totest, (bytes, int, float)):
toreturn = value
if isinstance(totest, (list, tuple)):
totest = np.array(value)
if isinstance(totest, np.ndarray) and totest.dtype.char == "U":
toreturn = np.array(totest).astype("S")
elif isinstance(totest, (np.ndarray, da.Array)):
toreturn = totest
if isinstance(totest, str):
toreturn = totest.encode("utf-8")
toreturn = np.bytes_(toreturn)
return toreturn
def _text_split(s, sep):
"""Split a string based of list of seperators.
Parameters
----------
s : str
sep : str - seperator or list of seperators e.g. '.' or ['_','/']
Returns
-------
list
String sections split based on the seperators
"""
stack = [s]
for char in sep:
pieces = []
for substr in stack:
pieces.extend(substr.split(char))
stack = pieces
if "" in stack:
stack.remove("")
return stack
def _getlink(h5group, rootkey, key):
"""Return the link target path.
If a hdf group is a soft link or has a target attribute
this method will return the target path. If no link is found
return None.
Returns
-------
str
Soft link path if it exists, otherwise None
"""
_target = None
if rootkey != "/":
if isinstance(h5group, h5py.Group):
_link = h5group.get(key, getlink=True)
if isinstance(_link, h5py.SoftLink):
_target = _link.path
if "target" in h5group.attrs.keys():
_target = _parse_from_file(h5group.attrs["target"])
if not _target.startswith("/"):
_target = "/" + _target
if _target == rootkey:
_target = None
return _target
def _get_nav_list(data, dataentry):
"""Get the list with information of each axes of the dataset
Parameters
----------
data : hdf dataset
the dataset to be loaded.
dataentry : hdf group
the group with corresponding attributes.
Returns
-------
nav_list : list
contains information about each axes.
"""
detector_index = 0
nav_list = []
# list indices...
axis_index_list = []
if "axes" in dataentry.attrs.keys():
axes_key = dataentry.attrs["axes"]
axes_list = ["."] * data.ndim
if isinstance(axes_key, np.ndarray):
axes_keys = axes_key[: data.ndim]
for i, num in enumerate(axes_keys):
axes_list[i] = _parse_from_file(num)
elif isinstance(axes_key, (str, bytes)):
axes_list = _parse_from_file(axes_key).split(",")[: data.ndim]
else:
axes_list[0] = _parse_from_file(axes_key)
named_axes = list(range(len(axes_list)))
for i, ax in enumerate(axes_list):
if ax != ".":
index_name = ax + "_indices"
if index_name in dataentry.attrs:
ind_in_array = dataentry.attrs[index_name]
if len(ind_in_array.shape) > 0:
ind_in_array = int(ind_in_array[0])
else:
ind_in_array = int(ind_in_array)
else:
ind_in_array = i
axis_index_list.append(ind_in_array)
if "units" in dataentry[ax].attrs:
units = _parse_from_file(dataentry[ax].attrs["units"])
else:
units = ""
navigation = True
named_axes.remove(ind_in_array)
if _is_numeric_data(dataentry[ax]):
if dataentry[ax].size > 1:
if _is_linear_axis(dataentry[ax]):
nav_list.append(
{
"size": data.shape[ind_in_array],
"index_in_array": ind_in_array,
"name": ax,
"scale": abs(dataentry[ax][1] - dataentry[ax][0]),
"offset": min(dataentry[ax][0], dataentry[ax][-1]),
"units": units,
"navigate": navigation,
}
)
else:
nav_list.append(
{
"size": data.shape[ind_in_array],
"index_in_array": ind_in_array,
"name": ax,
"scale": 1,
"offset": 0,
"navigate": navigation,
}
)
else:
nav_list.append(
{
"size": 1,
"index_in_array": ind_in_array,
"name": ax,
"scale": 1,
"offset": dataentry[ax][0],
"units": units,
"navigate": True,
}
)
else:
if len(data.shape) == len(axes_list):
nav_list.append(
{
"size": data.shape[named_axes[detector_index]],
"index_in_array": named_axes[detector_index],
"scale": 1,
"offset": 0.0,
"units": "",
"navigate": False,
}
)
detector_index = detector_index + 1
return nav_list
def _extract_hdf_dataset(group, dataset, lazy=False):
"""Import data from hdf path.
Parameters
----------
group : hdf group
group from which to load the dataset
dataset : str
path to the dataset within the group
lazy : bool {default:True}
If true use lazy opening, if false read into memory
Returns
-------
dict
A signal dictionary which can be used to instantiate a signal.
"""
data = group[dataset]
# exclude the dataset tagged by the signal attribute to avoid extracting
# duplicated dataset, which is already loaded when loading NeXus data
if "signal" in data.parent.attrs.keys():
data_key = data.parent.attrs["signal"]
if isinstance(data_key, bytes):
data_key = data_key.decode()
if dataset.split("/")[-1] == data_key:
return None
nav_list = _get_nav_list(data, data.parent)
if lazy:
if "chunks" in data.attrs.keys():
chunks = data.attrs["chunks"]
else:
signal_axes = [d["index_in_array"] for d in nav_list if not d["navigate"]]
chunks = get_signal_chunks(data.shape, data.dtype, signal_axes)
data_lazy = da.from_array(data, chunks=chunks)
else:
data_lazy = np.array(data)
dictionary = {
"data": data_lazy,
"metadata": {},
"original_metadata": {},
"axes": nav_list,
}
return dictionary
def _nexus_dataset_to_signal(group, nexus_dataset_path, lazy=False):
"""Load an NXdata set as a hyperspy signal.
Parameters
----------
group : hdf group containing the NXdata
nexus_data_path : str
Path to the NXdata set in the group
lazy : bool, default : True
lazy loading of data
Returns
-------
dict
A signal dictionary which can be used to instantiate a signal.
"""
interpretation = None
dataentry = group[nexus_dataset_path]
if "signal" in dataentry.attrs.keys():
if _is_int(dataentry.attrs["signal"]):
data_key = "data"
else:
data_key = dataentry.attrs["signal"]
else:
_logger.info(
"No signal attr associated with NXdata will\
try assume signal name is data"
)
if "data" not in dataentry.keys():
raise ValueError(
"Signal attribute not found in NXdata and "
'attempt to find a default "data" key failed'
)
else:
data_key = "data"
if "interpretation" in dataentry.attrs.keys():
interpretation = _parse_from_file(dataentry.attrs["interpretation"])
data = dataentry[data_key]
nav_list = _get_nav_list(data, dataentry)
if lazy:
if "chunks" in data.attrs.keys():
chunks = data.attrs["chunks"]
else:
signal_axes = [d["index_in_array"] for d in nav_list if not d["navigate"]]
chunks = get_signal_chunks(data.shape, data.dtype, signal_axes)
data_lazy = da.from_array(data, chunks=chunks)
else:
data_lazy = np.array(data)
if not nav_list:
for i in range(data.ndim):
nav_list.append(
{
"size": data_lazy.shape[i],
"index_in_array": i,
"scale": 1,
"offset": 0.0,
"units": "",
"navigate": True,
}
)
title = _text_split(nexus_dataset_path, "/")[-1]
metadata = {"General": {"title": title}}
#
# if interpretation - reset the nav axes
# assume the last dimensions are the signal
#
if interpretation:
for x in nav_list:
x["navigate"] = True
if interpretation == "spectrum":
nav_list[-1]["navigate"] = False
elif interpretation == "image":
nav_list[-1]["navigate"] = False
nav_list[-2]["navigate"] = False
dictionary = {"data": data_lazy, "axes": nav_list, "metadata": metadata}
return dictionary
def file_reader(
filename,
lazy=False,
dataset_key=None,
dataset_path=None,
metadata_key=None,
skip_array_metadata=False,
nxdata_only=False,
hardlinks_only=False,
use_default=False,
mapping=None,
):
"""
Read NXdata class or hdf datasets from a file and return signal(s).
Parameters
----------
%s
%s
dataset_key : None, str, list of str, default=None
If None all datasets are returned.
If a string or list of strings is provided only items
whose path contain the string(s) are returned. For example
dataset_key = ["instrument", "Fe"] will return
data entries with instrument or Fe in their hdf path.
dataset_path : None, str, list of str, default=None
If None, no absolute path is searched.
If a string or list of strings is provided items with the absolute
paths specified will be returned. For example, dataset_path =
['/data/spectrum/Mn'], it returns the exact dataset with this path.
It is not filtered by dataset_key, i.e. with dataset_key = ['Fe'],
it still returns the specific dataset at '/data/spectrum/Mn'. It is
empty if no dataset matching the absolute path provided is present.
metadata_key : None, str, list of str, default=None
Only return items from the original metadata whose path contain the
strings .e.g metadata_key = ["instrument", "Fe"] will return
all metadata entries with "instrument" or "Fe" in their hdf path.
skip_array_metadata : bool, default=False
Whether to skip loading metadata with an array entry. This is useful
as metadata may contain large array that is redundant with the data.
nxdata_only : bool, default=False
If True only NXdata will be converted into a signal
if False NXdata and any hdf datasets will be loaded as signals.
hardlinks_only : bool, default=False
If True any links (soft or External) will be ignored when loading.
use_default : bool, default=False
If True and a default NXdata is defined in the file load this as a
signal. This will ignore the other keyword options. If True and no
default is defined the file will be loaded according to
the keyword options.
mapping : None or dict
Define the mapping from the original metadata to the returned
metadata.
%s
See Also
--------
rsciio.utils.hdf5.list_datasets_in_file : Convenience function to list
datasets present in a file.
rsciio.utils.hdf5.read_metadata_from_file : Convenience function to read
metadata present in a file.
Notes
-----
Loading all datasets can result in a large number of signals
Please review your datasets and use the dataset_key to target
the datasets of interest.
"keys" is a special keywords and prepended with "fix" in the metadata
structure to avoid any issues.
Datasets are all arrays with size>2 (arrays, lists)
"""
# search for NXdata sets...
if mapping is None:
mapping = {}
original_metadata = {}
learning = {}
fin = h5py.File(filename, "r")
signal_dict_list = []
dataset_key = _check_search_keys(dataset_key)
dataset_path = _check_search_keys(dataset_path)
metadata_key = _check_search_keys(metadata_key)
original_metadata = _load_metadata(
fin, lazy=lazy, skip_array_metadata=skip_array_metadata
)
# some default values...
nexus_data_paths = []
hdf_data_paths = []
# check if a default dataset is defined
if use_default:
nexus_data_paths, hdf_data_paths = _find_data(
fin, search_keys=None, hardlinks_only=False
)
nxentry = None
nxdata = None
if "attrs" in original_metadata:
if "default" in original_metadata["attrs"]:
nxentry = original_metadata["attrs"]["default"]
else:
rootlist = list(original_metadata.keys())
rootlist.remove("attrs")
if rootlist and len(rootlist) == 1:
nxentry == rootlist[0]
if nxentry:
if "default" in original_metadata[nxentry]["attrs"]:
nxdata = original_metadata[nxentry]["attrs"]["default"]
if nxentry and nxdata:
nxdata = "/" + nxentry + "/" + nxdata
if nxdata:
hdf_data_paths = []
nexus_data_paths = [
nxpath for nxpath in nexus_data_paths if nxdata in nxpath
]
# if no default found then search for the data as normal
if not nexus_data_paths and not hdf_data_paths:
nexus_data_paths, hdf_data_paths = _find_data(
fin,
search_keys=dataset_key,
hardlinks_only=hardlinks_only,
absolute_path=dataset_path,
)
for data_path in nexus_data_paths:
dictionary = _nexus_dataset_to_signal(fin, data_path, lazy=lazy)
entryname = _text_split(data_path, "/")[0]
dictionary["mapping"] = mapping
title = dictionary["metadata"]["General"]["title"]
if entryname in original_metadata:
if metadata_key is None:
dictionary["original_metadata"] = original_metadata[entryname]
else:
dictionary["original_metadata"] = _find_search_keys_in_dict(
original_metadata, search_keys=metadata_key
)
# test if it's a hyperspy_nexus format and update metadata
# as appropriate.
if (
"attrs" in original_metadata
and "file_writer" in original_metadata["attrs"]
):
if original_metadata["attrs"]["file_writer"] == "hyperspy_nexus_v3":
orig_metadata = original_metadata[entryname]
if "auxiliary" in orig_metadata:
oma = orig_metadata["auxiliary"]
if "learning_results" in oma:
learning = oma["learning_results"]
dictionary["attributes"] = {}
dictionary["attributes"]["learning_results"] = learning
if "original_metadata" in oma:
if metadata_key is None:
dictionary["original_metadata"] = oma[
"original_metadata"
]
else:
dictionary["original_metadata"] = (
_find_search_keys_in_dict(
(oma["original_metadata"]),
search_keys=metadata_key,
)
)
# reconstruct the axes_list for axes_manager
for k, v in oma["original_metadata"].items():
if k.startswith("_sig_"):
hyper_ax = [
ax_v
for ax_k, ax_v in v.items()
if ax_k.startswith("_axes")
]
oma["original_metadata"][k]["axes"] = hyper_ax
if "hyperspy_metadata" in oma:
hyper_metadata = oma["hyperspy_metadata"]
hyper_metadata.update(dictionary["metadata"])
dictionary["metadata"] = hyper_metadata
else:
dictionary["original_metadata"] = {}
signal_dict_list.append(dictionary)
if not nxdata_only:
for data_path in hdf_data_paths:
datadict = _extract_hdf_dataset(fin, data_path, lazy=lazy)
if datadict:
title = data_path[1:].replace("/", "_")
basic_metadata = {
"General": {
"original_filename": os.path.split(filename)[1],
"title": title,
}
}
datadict["metadata"].update(basic_metadata)
signal_dict_list.append(datadict)
return signal_dict_list
file_reader.__doc__ %= (FILENAME_DOC, LAZY_DOC, RETURNS_DOC)
def _is_linear_axis(data):
"""Check if the data is linearly incrementing.
Parameters
----------
data : dask or numpy array
Returns
-------
bool
True or False
"""
steps = np.diff(data)
est_steps = np.array([steps[0]] * len(steps))
return np.allclose(est_steps, steps, rtol=1.0e-5)
def _is_numeric_data(data):
"""Check that data contains numeric data.
Parameters
----------
data : dask or numpy array
Returns
-------
bool
True or False
"""
try:
data.astype(float)
return True
except ValueError:
return False
def _is_int(s):
"""Check that s in an integer.
Parameters
----------
s : python object to test
Returns
-------
bool
True or False
"""
try:
int(s)
return True
except ValueError:
return False
def _check_search_keys(search_keys):
if isinstance(search_keys, str):
return [search_keys]
elif isinstance(search_keys, list):
if all(isinstance(key, str) for key in search_keys):
return search_keys
else:
raise ValueError("key list provided is not a list of strings")
elif search_keys is None:
return search_keys
else:
raise ValueError("search keys must be None, a string, " "or a list of strings")
def _find_data(group, search_keys=None, hardlinks_only=False, absolute_path=None):
"""Read from a nexus or hdf file and return a list of the dataset entries.
The method iterates through group attributes and returns NXdata or
hdf datasets of size >=2 if they're not already NXdata blocks
and returns a list of the entries
This is a convenience method to inspect a file to see which datasets
are present rather than loading all the sets in the file as signals
h5py.visit or visititems does not visit soft
links or external links so an implementation of a recursive
search is required. See https://github.com/h5py/h5py/issues/671
Parameters
----------
group : hdf group or File
search_keys : string, list of strings or None, default: None
Only return items which contain the strings
.e.g search_list = ["instrument","Fe"] will return
hdf entries with instrument or Fe in their hdf path.
hardlinks_only : bool , default : False
Option to ignore links (soft or External) within the file.
absolute_path : string, list of strings or None, default: None
Return items with the exact specified absolute path
Returns
-------
nx_dataset_list, hdf_dataset_list
nx_dataset_list is a list of all NXdata paths
hdf_dataset_list is a list of all hdf_datasets not linked to an
NXdata set.
"""
_check_search_keys(search_keys)
_check_search_keys(absolute_path)
all_hdf_datasets = []
unique_hdf_datasets = []
all_nx_datasets = []
unique_nx_datasets = []
rootname = ""
def find_data_in_tree(group, rootname):
for key, value in group.items():
if rootname != "":
rootkey = rootname + "/" + key
else:
rootkey = "/" + key
if isinstance(value, h5py.Group):
target = _getlink(group, rootkey, key)
if "NX_class" in value.attrs:
if (
value.attrs["NX_class"] in [b"NXdata", "NXdata"]
and "signal" in value.attrs.keys()
):
all_nx_datasets.append(rootkey)
if target is None:
unique_nx_datasets.append(rootkey)
if hardlinks_only:
if target is None:
find_data_in_tree(value, rootkey)
else:
find_data_in_tree(value, rootkey)
else:
if isinstance(value, h5py.Dataset):
if value.size >= 2:
target = _getlink(group, rootkey, key)
if not (value.dtype.type is str or value.dtype.type is object):
all_hdf_datasets.append(rootkey)
if target is None:
unique_hdf_datasets.append(rootkey)
# need to use custom recursive function as visititems in h5py
# does not visit links
find_data_in_tree(group, rootname)
if search_keys is None and absolute_path is None:
# return all datasets
if hardlinks_only:
# return only the stored data, no linked data
return unique_nx_datasets, unique_hdf_datasets
else:
return all_nx_datasets, all_hdf_datasets
elif isinstance(search_keys, list) or isinstance(absolute_path, list):
if hardlinks_only:
# return only the stored data, no linked data
nx_datasets = unique_nx_datasets
hdf_datasets = unique_hdf_datasets
else:
nx_datasets = all_nx_datasets
hdf_datasets = all_hdf_datasets
matched_hdf = set()
matched_nexus = set()
# return data having the specified absolute paths
if absolute_path is not None:
matched_hdf.update(
[j for j in hdf_datasets if any(s == j for s in absolute_path)]
)
matched_nexus.update(
[j for j in nx_datasets if any(s == j for s in absolute_path)]
)
# return data which contains a search string
if search_keys is not None:
matched_hdf.update(
[j for j in hdf_datasets if any(s in j for s in search_keys)]
)
matched_nexus.update(
[j for j in nx_datasets if any(s in j for s in search_keys)]
)
matched_nexus = list(matched_nexus)
matched_hdf = list(matched_hdf)
return matched_nexus, matched_hdf
def _load_metadata(group, lazy=False, skip_array_metadata=False):
"""Search through a hdf group and return the group structure.
h5py.visit or visititems does not visit soft
links or external links so an implementation of a recursive
search is required. See https://github.com/h5py/h5py/issues/671
Parameters
----------
group : hdf group
location to load the metadata from
lazy : bool , default : False
Option for lazy loading
skip_array_metadata : bool, default : False
whether to skip loading array metadata
Returns
-------
dict
dictionary of group contents
"""
rootname = ""
def find_meta_in_tree(group, rootname, lazy=False, skip_array_metadata=False):
tree = {}
for key, item in group.attrs.items():
new_key = _fix_exclusion_keys(key)
if "attrs" not in tree.keys():
tree["attrs"] = {}
tree["attrs"][new_key] = _parse_from_file(item, lazy=lazy)
for key, item in group.items():
if rootname != "":
rootkey = rootname + "/" + key
else:
rootkey = "/" + key
new_key = _fix_exclusion_keys(key)
if isinstance(item, h5py.Dataset):
if item.attrs:
if new_key not in tree.keys():
tree[new_key] = {}
# avoid loading array as metadata
if skip_array_metadata:
if item.size < 2 and item.ndim == 0:
tree[new_key]["value"] = _parse_from_file(item, lazy=lazy)
else:
tree[new_key]["value"] = _parse_from_file(item, lazy=lazy)
for k, v in item.attrs.items():
if "attrs" not in tree[new_key].keys():
tree[new_key]["attrs"] = {}
tree[new_key]["attrs"][k] = _parse_from_file(v, lazy=lazy)
else:
# this is to support hyperspy where datasets are not saved
# with attributes
if skip_array_metadata:
if item.size < 2 and item.ndim == 0:
tree[new_key] = _parse_from_file(item, lazy=lazy)
else:
tree[new_key] = _parse_from_file(item, lazy=lazy)
elif isinstance(item, h5py.Group):
if "NX_class" in item.attrs:
if item.attrs["NX_class"] not in [b"NXdata", "NXdata"]:
tree[new_key] = find_meta_in_tree(
item,
rootkey,
lazy=lazy,
skip_array_metadata=skip_array_metadata,
)
else:
tree[new_key] = find_meta_in_tree(
item,
rootkey,
lazy=lazy,
skip_array_metadata=skip_array_metadata,
)
return tree
extracted_tree = find_meta_in_tree(
group, rootname, lazy=lazy, skip_array_metadata=skip_array_metadata
)
return extracted_tree
def _fix_exclusion_keys(key):
"""Exclude hyperspy specific keys.
Signal and DictionaryBrowser break if a
a key is a dict method - e.g. {"keys":2.0}.
This method prepends the key with ``fix_`` so the information is
still present to work around this issue
Parameters
----------
key : str
Returns
-------
str
"""
if key.startswith("keys"):
return "fix_" + key
else:
return key
def _find_search_keys_in_dict(tree, search_keys=None):
"""Search through a dict for search keys.
This is a convenience method to inspect a file for a value
rather than loading the file as a signal
Parameters
----------
tree : h5py File object
search_keys : string or list of strings
Only return items which contain the strings
.e.g search_keys = ["instrument","Fe"] will return
hdf entries with instrument or Fe in their hdf path.
Returns
-------
dict
When search_list is specified only full paths
containing one or more search_keys will be returned
"""
_check_search_keys(search_keys)
metadata_dict = {}
rootname = ""
# recursive function
def find_searchkeys_in_tree(myDict, rootname):
for key, value in myDict.items():
if rootname != "":
rootkey = rootname + "/" + key
else:
rootkey = key
if isinstance(search_keys, list) and any(
[s1 in rootkey for s1 in search_keys]
):
mod_keys = _text_split(rootkey, (".", "/"))
# create the key, values in the dict
p = metadata_dict
for d in mod_keys[:-1]:
p = p.setdefault(d, {})
p[mod_keys[-1]] = value
if isinstance(value, dict):
find_searchkeys_in_tree(value, rootkey)
if search_keys is None:
return tree
else:
find_searchkeys_in_tree(tree, rootname)
return metadata_dict
def _write_nexus_groups(dictionary, group, skip_keys=None, **kwds):
"""Recursively iterate throuh dictionary and write groups to nexus.
Parameters
----------
dictionary : dict
dictionary contents to store to hdf group
group : hdf group
location to store dictionary
skip_keys : str or list of str
the key(s) to skip when writing into the group
**kwds : additional keywords
additional keywords to pass to h5py.create_dataset method
"""
if skip_keys is None:
skip_keys = []
elif isinstance(skip_keys, str):
skip_keys = [skip_keys]
for key, value in dictionary.items():
if key == "attrs" or key in skip_keys:
# we will handle attrs later... and skip unwanted keys
continue
if isinstance(value, dict):
if "attrs" in value:
if (
"NX_class" in value["attrs"]
and value["attrs"]["NX_class"] == "NXdata"
):
continue
if (
"value" in value.keys()
and not isinstance(value["value"], dict)
and len(set(list(value.keys()) + ["attrs", "value"])) == 2
):
value = value["value"]
else:
_write_nexus_groups(
value, group.require_group(key), skip_keys=skip_keys, **kwds
)
if isinstance(value, (list, tuple, np.ndarray, da.Array)):
if all(isinstance(v, dict) for v in value):
# a list of dictionary is from the axes of HyperSpy signal
for i, ax_dict in enumerate(value):
ax_key = "_axes_" + str(i)
_write_nexus_groups(
ax_dict,
group.require_group(ax_key),
skip_keys=skip_keys,
**kwds,
)
else:
data = _parse_to_file(value)
overwrite_dataset(group, data, key, chunks=None, **kwds)
elif isinstance(value, (int, float, str, bytes)):
group.create_dataset(key, data=_parse_to_file(value))
else:
if value is not None and key not in group:
_write_nexus_groups(
value, group.require_group(key), skip_keys=skip_keys, **kwds
)
def _write_nexus_attr(dictionary, group, skip_keys=None):
"""Recursively iterate through dictionary and write "attrs" dictionaries.
This step is called after the groups and datasets have been created
Parameters
----------
dictionary : dict
Input dictionary to be written to the hdf group
group : hdf group
location to store the attrs sections of the dictionary
"""
if skip_keys is None:
skip_keys = []
elif isinstance(skip_keys, str):
skip_keys = [skip_keys]
for key, value in dictionary.items():
if key == "attrs":
for k, v in value.items():
group.attrs[k] = _parse_to_file(v)
else:
if isinstance(value, dict):
if "attrs" in value:
if (
"NX_class" in value["attrs"]
and value["attrs"]["NX_class"] == "NXdata"
):
continue
if key not in skip_keys:
_write_nexus_attr(dictionary[key], group[key], skip_keys=skip_keys)
def _write_signal(signal, nxgroup, signal_name, **kwds):
"""Store the signal data as an NXdata dataset.
Parameters
----------
signal : HyperSpy signal dictionary
nxgroup : HDF group
Entry at which to save signal data
signal_name : str
Name under which to store the signal entry in the file
"""
signal_axes = [ax for ax in signal["axes"] if not ax["navigate"]]
data = signal["data"]
if len(signal_axes) == 1:
record_by = "spectrum"
elif len(signal_axes) == 2:
record_by = "image"
else:
record_by = ""
nxdata = nxgroup.require_group(signal_name)
nxdata.attrs["NX_class"] = _parse_to_file("NXdata")
nxdata.attrs["signal"] = _parse_to_file("data")
if record_by:
nxdata.attrs["interpretation"] = _parse_to_file(record_by)
overwrite_dataset(
nxdata,
data,
"data",
chunks=None,
signal_axes=[signal["axes"].index(ax) for ax in signal_axes],
**kwds,
)
axis_names = [_parse_to_file(".")] * len(data.shape)
for i, ax in enumerate(signal["axes"]):
if ax["name"] is not None:
index_in_array = signal["axes"].index(ax)
indices = _parse_to_file(ax["name"] + "_indices")
nxdata.attrs[indices] = _parse_to_file(index_in_array)
if ax["_type"] == "DataAxis":
data = ax["axis"]
elif ax["_type"] == "UniformDataAxis":
data = ax["offset"] + ax["scale"] * np.arange(ax["size"])
elif ax["_type"] == "FunctionalDataAxis":
raise ValueError("Not supported")
else:
raise RuntimeError(f"Axis of type ({ax['_type']}) are not supported.")
dset = nxdata.require_dataset(
ax["name"], data=data, shape=data.shape, dtype=data.dtype
)
if ax["units"] is not None:
dset.attrs["units"] = ax["units"]
axis_names[index_in_array] = ax["name"]
nxdata.attrs["axes"] = _parse_to_file(axis_names)
return nxdata
def file_writer(
filename,
signal,
save_original_metadata=True,
skip_metadata_key=None,
use_default=False,
compression="gzip",
):
"""
Write the signal and metadata as a NeXus file.
This will save the signal in NXdata format in the file.
As the form of the metadata can vary and is not validated it will
be stored as an NXcollection (an unvalidated collection)
Parameters
----------
%s
%s
save_original_metadata : bool , default : True
Option to save hyperspy.original_metadata with the signal.
A loaded NeXus file may have a large amount of data
when loaded which you may wish to omit on saving.
skip_metadata_key : str or list of str, default : None
The key(s) to skip when saving original metadata. This is useful
when some metadata keys should be ignored.
use_default : bool , default : False
Option to define the default dataset in the file.
If set to True the signal or first signal in the list of signals
will be defined as the default (following NeXus v3 data rules).
%s
See Also
--------
rsciio.utils.hdf5.list_datasets_in_file : Convenience function to list
datasets present in a file.
rsciio.utils.hdf5.read_metadata_from_file : Convenience function to read
metadata present in a file.
Notes
-----
%s
"""
if not isinstance(signal, list):
signal = [signal]
with h5py.File(filename, mode="w") as f:
f.attrs["file_format"] = "nexus"
f.attrs["file_writer"] = "hyperspy_nexus_v3"
if use_default:
f.attrs["default"] = "entry1"
#
# write the signal
#
for i, sig in enumerate(signal):
md = DTBox(sig["metadata"], box_dots=True)
nxentry = f.create_group(f"entry{i + 1}")
nxentry.attrs["NX_class"] = _parse_to_file("NXentry")
signal_name = md.get("General.title")
if not signal_name:
signal_name = f"unnamed__{i}"
if "/" in signal_name:
signal_name = signal_name.replace("/", "_")
if signal_name.startswith("__"):
signal_name = signal_name[2:]
if i == 0 and use_default:
nxentry.attrs["default"] = signal_name
nxaux = nxentry.create_group("auxiliary")
nxaux.attrs["NX_class"] = _parse_to_file("NXentry")
_write_signal(sig, nxentry, signal_name, compression=compression)
learn = sig.get("learning_results")
if learn:
nxlearn = nxaux.create_group("learning_results")
nxlearn.attrs["NX_class"] = _parse_to_file("NXcollection")
_write_nexus_groups(learn, nxlearn, compression=compression)
_write_nexus_attr(learn, nxlearn)
#
# write metadata
#
if save_original_metadata:
om = sig.get("original_metadata")
if om:
nxom = nxaux.create_group("original_metadata")
nxom.attrs["NX_class"] = _parse_to_file("NXcollection")
# write the groups and structure
_write_nexus_groups(
om, nxom, skip_keys=skip_metadata_key, compression=compression
)
_write_nexus_attr(om, nxom, skip_keys=skip_metadata_key)
md = sig.get("metadata")
if md:
nxmd = nxaux.create_group("hyperspy_metadata")
nxmd.attrs["NX_class"] = _parse_to_file("NXcollection")
# write the groups and structure
_write_nexus_groups(md, nxmd, compression=compression)
_write_nexus_attr(md, nxmd)
file_writer.__doc__ %= (
FILENAME_DOC.replace("read", "write to"),
SIGNAL_DOC,
COMPRESSION_HDF5_DOC,
COMPRESSION_HDF5_NOTES_DOC,
)
|