autodetect base group
This commit is contained in:
parent
35c5bfcc45
commit
cdcedd0d44
|
@ -98,7 +98,7 @@ class ConfigMaterial(Config):
|
|||
|
||||
|
||||
@staticmethod
|
||||
def load_DREAM3D(fname,base_group,data_group,ori_data,phase_id,phase_name):
|
||||
def load_DREAM3D(fname,data_group,ori_data,phase_id,phase_name,base_group=None):
|
||||
"""
|
||||
Load material data from DREAM3D file.
|
||||
|
||||
|
@ -150,10 +150,10 @@ class ConfigMaterial(Config):
|
|||
... 'EulerAngles','Phases',['Ferrite','Martensite'])
|
||||
|
||||
"""
|
||||
root_dir = 'DataContainers'
|
||||
b = util.DREAM3D_base_group(fname) if base_group is None else base_group
|
||||
hdf = h5py.File(fname,'r')
|
||||
|
||||
orientation_path = path.join(root_dir,base_group,data_group,ori_data)
|
||||
orientation_path = path.join(b,data_group,ori_data)
|
||||
if hdf[orientation_path].attrs['TupleDimensions'].shape == (3,):
|
||||
grain_orientations = np.array(hdf[orientation_path]).reshape(-1,3,order='F')
|
||||
else:
|
||||
|
@ -161,7 +161,7 @@ class ConfigMaterial(Config):
|
|||
|
||||
grain_quats = Rotation.from_Euler_angles(grain_orientations).as_quaternion()
|
||||
|
||||
phase_path = path.join(root_dir,base_group,data_group,phase_id)
|
||||
phase_path = path.join(b,data_group,phase_id)
|
||||
if hdf[phase_path].attrs['TupleDimensions'].shape == (3,):
|
||||
grain_phase = np.array(hdf[phase_path]).reshape(-1,order='F')
|
||||
else:
|
||||
|
|
|
@ -256,7 +256,7 @@ class Grid:
|
|||
|
||||
|
||||
@staticmethod
|
||||
def load_DREAM3D(fname,base_group,point_data=None,material='FeatureIds'):
|
||||
def load_DREAM3D(fname,cell_data=None,material='FeatureIds',base_group=None):
|
||||
"""
|
||||
Load from DREAM.3D file.
|
||||
|
||||
|
@ -264,26 +264,26 @@ class Grid:
|
|||
----------
|
||||
fname : str
|
||||
Filename of the DREAM.3D file
|
||||
base_group : str
|
||||
Name of the group (folder) below 'DataContainers',
|
||||
for example 'SyntheticVolumeDataContainer'.
|
||||
point_data : str, optional
|
||||
cell_data : str, optional
|
||||
Name of the group (folder) containing the pointwise material data,
|
||||
for example 'CellData'. Defaults to None, in which case points are consecutively numbered.
|
||||
material : str, optional
|
||||
Name of the dataset containing the material ID.
|
||||
Defaults to 'FeatureIds'.
|
||||
|
||||
base_group : str
|
||||
Path to the group (folder) that contains the geometry (_SIMPL_GEOMETRY),
|
||||
and, optionally, the cell data. Defaults to None, in which case
|
||||
it is set as the path that contains _SIMPL_GEOMETRY/SPACING.
|
||||
|
||||
"""
|
||||
b = util.DREAM3D_base_group(fname) if base_group is None else base_group
|
||||
f = h5py.File(fname, 'r')
|
||||
g = os.path.join('DataContainers',base_group,'_SIMPL_GEOMETRY')
|
||||
cells = f[os.path.join(g,'DIMENSIONS')][()]
|
||||
size = f[os.path.join(g,'SPACING')][()] * cells
|
||||
origin = f[os.path.join(g,'ORIGIN')][()]
|
||||
cells = f[os.path.join(b,'_SIMPL_GEOMETRY','DIMENSIONS')][()]
|
||||
size = f[os.path.join(b,'_SIMPL_GEOMETRY','SPACING')][()] * cells
|
||||
origin = f[os.path.join(b,'_SIMPL_GEOMETRY','ORIGIN')][()]
|
||||
|
||||
ma = np.arange(cells.prod(),dtype=int) \
|
||||
if point_data is None else \
|
||||
np.reshape(f[os.path.join('DataContainers',base_group,point_data,material)],cells.prod())
|
||||
ma = np.arange(cells.prod(),dtype=int) if cell_data is None else \
|
||||
np.reshape(f[os.path.join(b,cell_data,material)],cells.prod())
|
||||
|
||||
return Grid(ma.reshape(cells,order='F'),size,origin,util.execution_stamp('Grid','load_DREAM3D'))
|
||||
|
||||
|
|
|
@ -9,6 +9,7 @@ from functools import reduce
|
|||
from optparse import Option
|
||||
|
||||
import numpy as np
|
||||
import h5py
|
||||
|
||||
from . import version
|
||||
|
||||
|
@ -27,7 +28,8 @@ __all__=[
|
|||
'extendableOption',
|
||||
'execution_stamp',
|
||||
'shapeshifter', 'shapeblender',
|
||||
'extend_docstring', 'extended_docstring'
|
||||
'extend_docstring', 'extended_docstring',
|
||||
'DREAM3D_base_group'
|
||||
]
|
||||
|
||||
####################################################################################################
|
||||
|
@ -376,6 +378,15 @@ def extended_docstring(f,extra_docstring):
|
|||
return _decorator
|
||||
|
||||
|
||||
def DREAM3D_base_group(fname):
|
||||
with h5py.File(fname,'r') as f:
|
||||
base_group = f.visit(lambda path: path.rsplit('/',2)[0] if '_SIMPL_GEOMETRY/SPACING' in path else None)
|
||||
|
||||
if base_group is None:
|
||||
raise ValueError
|
||||
|
||||
return base_group
|
||||
|
||||
####################################################################################################
|
||||
# Classes
|
||||
####################################################################################################
|
||||
|
|
Loading…
Reference in New Issue