let the computer do the work

This commit is contained in:
Martin Diehl 2021-03-23 14:28:56 +01:00
parent e0e088eaa8
commit 572c3204d0
5 changed files with 58 additions and 15 deletions

View File

@ -99,7 +99,7 @@ class ConfigMaterial(Config):
@staticmethod
def load_DREAM3D(fname,
grain_data=None,cell_data='CellData',cell_ensemble_data='CellEnsembleData',
grain_data=None,cell_data=None,cell_ensemble_data='CellEnsembleData',
phases='Phases',Euler_angles='EulerAngles',phase_names='PhaseName',
base_group=None):
"""
@ -120,7 +120,8 @@ class ConfigMaterial(Config):
Name of the group (folder) containing grain-wise data. Defaults
to None, in which case cell-wise data is used.
cell_data : str
Name of the group (folder) containing cell-wise data. Defaults to 'CellData'.
Name of the group (folder) containing cell-wise data. Defaults to
None in wich case it is automatically detected.
cell_ensemble_data : str
Name of the group (folder) containing data of cell ensembles. This
group is used to inquire the name of the phases. Phases will get
@ -141,11 +142,12 @@ class ConfigMaterial(Config):
"""
b = util.DREAM3D_base_group(fname) if base_group is None else base_group
c = util.DREAM3D_cell_data_group(fname) if cell_data is None else cell_data
f = h5py.File(fname,'r')
if grain_data is None:
phase = f[os.path.join(b,cell_data,phases)][()].flatten()
O = Rotation.from_Euler_angles(f[os.path.join(b,cell_data,Euler_angles)]).as_quaternion().reshape(-1,4) # noqa
phase = f[os.path.join(b,c,phases)][()].flatten()
O = Rotation.from_Euler_angles(f[os.path.join(b,c,Euler_angles)]).as_quaternion().reshape(-1,4) # noqa
_,idx = np.unique(np.hstack([O,phase.reshape(-1,1)]),return_index=True,axis=0)
idx = np.sort(idx)
else:

View File

@ -257,7 +257,7 @@ class Grid:
@staticmethod
def load_DREAM3D(fname,
feature_IDs=None,cell_data='CellData',
feature_IDs=None,cell_data=None,
phases='Phases',Euler_angles='EulerAngles',
base_group=None):
"""
@ -279,7 +279,8 @@ class Grid:
grain-wise data. Defaults to 'None', in which case cell-wise
data is used.
cell_data : str
Name of the group (folder) containing cell-wise data. Defaults to 'CellData'.
Name of the group (folder) containing cell-wise data. Defaults to
None in wich case it is automatically detected.
phases : str
Name of the dataset containing the phase ID. It is not used for
grain-wise data, i.e. when feature_IDs is not None.
@ -296,19 +297,21 @@ class Grid:
"""
b = util.DREAM3D_base_group(fname) if base_group is None else base_group
c = util.DREAM3D_cell_data_group(fname) if cell_data is None else cell_data
f = h5py.File(fname, 'r')
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')][()]
if feature_IDs is None:
phase = f[os.path.join(b,cell_data,phases)][()].reshape(-1,1)
O = Rotation.from_Euler_angles(f[os.path.join(b,cell_data,Euler_angles)]).as_quaternion().reshape(-1,4) # noqa
phase = f[os.path.join(b,c,phases)][()].reshape(-1,1)
O = Rotation.from_Euler_angles(f[os.path.join(b,c,Euler_angles)]).as_quaternion().reshape(-1,4) # noqa
unique,unique_inverse = np.unique(np.hstack([O,phase]),return_inverse=True,axis=0)
ma = np.arange(cells.prod()) if len(unique) == cells.prod() else \
np.arange(unique.size)[np.argsort(pd.unique(unique_inverse))][unique_inverse]
else:
ma = f[os.path.join(b,cell_data,feature_IDs)][()].flatten()
ma = f[os.path.join(b,c,feature_IDs)][()].flatten()
return Grid(ma.reshape(cells,order='F'),size,origin,util.execution_stamp('Grid','load_DREAM3D'))

View File

@ -29,7 +29,7 @@ __all__=[
'execution_stamp',
'shapeshifter', 'shapeblender',
'extend_docstring', 'extended_docstring',
'DREAM3D_base_group'
'DREAM3D_base_group', 'DREAM3D_cell_data_group'
]
####################################################################################################
@ -379,14 +379,52 @@ def extended_docstring(f,extra_docstring):
def DREAM3D_base_group(fname):
"""
Determine the base group of a DREAM.3D file.
The base group is defined as the group (folder) that contains
a 'SPACING' dataset in a '_SIMPL_GEOMETRY' group.
Parameters
----------
fname : str
Filename of the DREAM.3D (HDF5) file.
"""
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
raise ValueError('Could not determine base group in file {fname}.')
return base_group
def DREAM3D_cell_data_group(fname):
"""
Determine the cell data group of a DREAM.3D file.
The cell data group is defined as the group (folder) that contains
a dataset in the base group whose length matches the total number
of points as specified in '_SIMPL_GEOMETRY/DIMENSIONS'.
Parameters
----------
fname : str
Filename of the DREAM.3D (HDF5) file.
"""
base_group = DREAM3D_base_group(fname)
with h5py.File(fname,'r') as f:
N_points = np.prod(f[os.path.join(base_group,'_SIMPL_GEOMETRY','DIMENSIONS')])
cell_data_group = f[base_group].visititems(lambda path,obj: path.split('/')[0] \
if isinstance(obj,h5py._hl.dataset.Dataset) and np.prod(np.shape(obj)) == N_points \
else None)
if cell_data_group is None:
raise ValueError('Could not determine cell data group in file {fname}.')
return cell_data_group
####################################################################################################
# Classes
####################################################################################################

View File

@ -133,7 +133,7 @@ class TestConfigMaterial:
def test_load_DREAM3D_reference(self,tmp_path,ref_path,update):
config = ConfigMaterial.load_DREAM3D(ref_path/'measured.dream3d',cell_data='EBSD Scan Data')
config = ConfigMaterial.load_DREAM3D(ref_path/'measured.dream3d')
config.save(tmp_path/'material.yaml')
if update:
config.save(ref_path/'measured.material_yaml')

View File

@ -442,7 +442,7 @@ class TestGrid:
def test_load_DREAM3D_reference(self,ref_path,update):
current = Grid.load_DREAM3D(ref_path/'measured.dream3d',cell_data='EBSD Scan Data')
current = Grid.load_DREAM3D(ref_path/'measured.dream3d')
reference = Grid.load(ref_path/'measured')
if update:
current.save(ref_path/'measured.vtr')