Merge branch 'load-DREAM3d-improvements' into 'development'
Load dream3d improvements See merge request damask/DAMASK!821
This commit is contained in:
commit
57d8a6c438
|
@ -106,15 +106,12 @@ class ConfigMaterial(Config):
|
|||
Load DREAM.3D (HDF5) file.
|
||||
|
||||
Data in DREAM.3D files can be stored per cell ('CellData')
|
||||
and/or per grain ('Grain Data'). Per default, cell-wise data
|
||||
is assumed.
|
||||
|
||||
damask.Grid.load_DREAM3D allows to get the corresponding geometry
|
||||
for the grid solver.
|
||||
and/or per grain ('Grain Data'). Per default, i.e. if
|
||||
'grain_data' is None, cell-wise data is assumed.
|
||||
|
||||
Parameters
|
||||
----------
|
||||
fname : str
|
||||
fname : str or pathlib.Path
|
||||
Filename of the DREAM.3D (HDF5) file.
|
||||
grain_data : str
|
||||
Name of the group (folder) containing grain-wise data. Defaults
|
||||
|
@ -140,36 +137,43 @@ class ConfigMaterial(Config):
|
|||
and grain- or cell-wise data. Defaults to None, in which case
|
||||
it is set as the path that contains _SIMPL_GEOMETRY/SPACING.
|
||||
|
||||
Notes
|
||||
-----
|
||||
Homogenization and phase entries are emtpy and need to be defined separately.
|
||||
|
||||
Returns
|
||||
-------
|
||||
loaded : damask.ConfigMaterial
|
||||
Material configuration from file.
|
||||
|
||||
Notes
|
||||
-----
|
||||
damask.Grid.load_DREAM3D gives the corresponding geometry for
|
||||
the grid solver.
|
||||
|
||||
For cell-wise data, only unique combinations of
|
||||
orientation and phase are considered.
|
||||
|
||||
Homogenization and phase entries are emtpy and need to be
|
||||
defined separately.
|
||||
|
||||
"""
|
||||
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')
|
||||
with h5py.File(fname, 'r') as f:
|
||||
b = util.DREAM3D_base_group(f) if base_group is None else base_group
|
||||
c = util.DREAM3D_cell_data_group(f) if cell_data is None else cell_data
|
||||
|
||||
if grain_data is None:
|
||||
phase = f['/'.join([b,c,phases])][()].flatten()
|
||||
O = Rotation.from_Euler_angles(f['/'.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:
|
||||
phase = f['/'.join([b,grain_data,phases])][()]
|
||||
O = Rotation.from_Euler_angles(f['/'.join([b,grain_data,Euler_angles])]).as_quaternion() # noqa
|
||||
idx = np.arange(phase.size)
|
||||
if grain_data is None:
|
||||
phase = f['/'.join([b,c,phases])][()].flatten()
|
||||
O = Rotation.from_Euler_angles(f['/'.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:
|
||||
phase = f['/'.join([b,grain_data,phases])][()]
|
||||
O = Rotation.from_Euler_angles(f['/'.join([b,grain_data,Euler_angles])]).as_quaternion() # noqa
|
||||
idx = np.arange(phase.size)
|
||||
|
||||
if cell_ensemble_data is not None and phase_names is not None:
|
||||
try:
|
||||
names = np.array([s.decode() for s in f['/'.join([b,cell_ensemble_data,phase_names])]])
|
||||
phase = names[phase]
|
||||
except KeyError:
|
||||
pass
|
||||
if cell_ensemble_data is not None and phase_names is not None:
|
||||
try:
|
||||
names = np.array([s.decode() for s in f['/'.join([b,cell_ensemble_data,phase_names])]])
|
||||
phase = names[phase]
|
||||
except KeyError:
|
||||
pass
|
||||
|
||||
|
||||
base_config = ConfigMaterial({'phase':{k if isinstance(k,int) else str(k): None for k in np.unique(phase)},
|
||||
|
|
|
@ -358,14 +358,14 @@ class Grid:
|
|||
"""
|
||||
Load DREAM.3D (HDF5) file.
|
||||
|
||||
Data in DREAM.3D files can be stored per cell ('CellData') and/or
|
||||
per grain ('Grain Data'). Per default, cell-wise data is assumed.
|
||||
Data in DREAM.3D files can be stored per cell ('CellData')
|
||||
and/or per grain ('Grain Data'). Per default, i.e. if
|
||||
'feature_IDs' is None, cell-wise data is assumed.
|
||||
|
||||
damask.ConfigMaterial.load_DREAM3D gives the corresponding material definition.
|
||||
|
||||
Parameters
|
||||
----------
|
||||
fname : str or or pathlib.Path
|
||||
fname : str or pathlib.Path
|
||||
Filename of the DREAM.3D (HDF5) file.
|
||||
feature_IDs : str, optional
|
||||
Name of the dataset containing the mapping between cells and
|
||||
|
@ -392,23 +392,31 @@ class Grid:
|
|||
loaded : damask.Grid
|
||||
Grid-based geometry from file.
|
||||
|
||||
Notes
|
||||
-----
|
||||
damask.ConfigMaterial.load_DREAM3D gives the corresponding
|
||||
material definition.
|
||||
|
||||
For cell-wise data, only unique combinations of
|
||||
orientation and phase are considered.
|
||||
|
||||
"""
|
||||
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')
|
||||
with h5py.File(fname, 'r') as f:
|
||||
b = util.DREAM3D_base_group(f) if base_group is None else base_group
|
||||
c = util.DREAM3D_cell_data_group(f) if cell_data is None else cell_data
|
||||
|
||||
cells = f['/'.join([b,'_SIMPL_GEOMETRY','DIMENSIONS'])][()]
|
||||
size = f['/'.join([b,'_SIMPL_GEOMETRY','SPACING'])] * cells
|
||||
origin = f['/'.join([b,'_SIMPL_GEOMETRY','ORIGIN'])][()]
|
||||
cells = f['/'.join([b,'_SIMPL_GEOMETRY','DIMENSIONS'])][()]
|
||||
size = f['/'.join([b,'_SIMPL_GEOMETRY','SPACING'])] * cells
|
||||
origin = f['/'.join([b,'_SIMPL_GEOMETRY','ORIGIN'])][()]
|
||||
|
||||
if feature_IDs is None:
|
||||
phase = f['/'.join([b,c,phases])][()].reshape(-1,1)
|
||||
O = Rotation.from_Euler_angles(f['/'.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['/'.join([b,c,feature_IDs])][()].flatten()
|
||||
if feature_IDs is None:
|
||||
phase = f['/'.join([b,c,phases])][()].reshape(-1,1)
|
||||
O = Rotation.from_Euler_angles(f['/'.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['/'.join([b,c,feature_IDs])][()].flatten()
|
||||
|
||||
return Grid(material = ma.reshape(cells,order='F'),
|
||||
size = size,
|
||||
|
|
|
@ -698,7 +698,7 @@ def pass_on(keyword: str,
|
|||
return wrapper
|
||||
return decorator
|
||||
|
||||
def DREAM3D_base_group(fname: _Union[str, _Path]) -> str:
|
||||
def DREAM3D_base_group(fname: _Union[str, _Path, _h5py.File]) -> str:
|
||||
"""
|
||||
Determine the base group of a DREAM.3D file.
|
||||
|
||||
|
@ -707,7 +707,7 @@ def DREAM3D_base_group(fname: _Union[str, _Path]) -> str:
|
|||
|
||||
Parameters
|
||||
----------
|
||||
fname : str or pathlib.Path
|
||||
fname : str, pathlib.Path, or _h5py.File
|
||||
Filename of the DREAM.3D (HDF5) file.
|
||||
|
||||
Returns
|
||||
|
@ -716,15 +716,19 @@ def DREAM3D_base_group(fname: _Union[str, _Path]) -> str:
|
|||
Path to the base group.
|
||||
|
||||
"""
|
||||
with _h5py.File(_Path(fname).expanduser(),'r') as f:
|
||||
def get_base_group(f: _h5py.File) -> str:
|
||||
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(f'could not determine base group in file "{fname}"')
|
||||
return base_group
|
||||
|
||||
if base_group is None:
|
||||
raise ValueError(f'could not determine base group in file "{fname}"')
|
||||
if isinstance(fname,_h5py.File):
|
||||
return get_base_group(fname)
|
||||
|
||||
return base_group
|
||||
with _h5py.File(_Path(fname).expanduser(),'r') as f:
|
||||
return get_base_group(f)
|
||||
|
||||
def DREAM3D_cell_data_group(fname: _Union[str, _Path]) -> str:
|
||||
def DREAM3D_cell_data_group(fname: _Union[str, _Path, _h5py.File]) -> str:
|
||||
"""
|
||||
Determine the cell data group of a DREAM.3D file.
|
||||
|
||||
|
@ -734,7 +738,7 @@ def DREAM3D_cell_data_group(fname: _Union[str, _Path]) -> str:
|
|||
|
||||
Parameters
|
||||
----------
|
||||
fname : str or pathlib.Path
|
||||
fname : str, pathlib.Path, or h5py.File
|
||||
Filename of the DREAM.3D (HDF5) file.
|
||||
|
||||
Returns
|
||||
|
@ -743,17 +747,21 @@ def DREAM3D_cell_data_group(fname: _Union[str, _Path]) -> str:
|
|||
Path to the cell data group.
|
||||
|
||||
"""
|
||||
base_group = DREAM3D_base_group(fname)
|
||||
with _h5py.File(_Path(fname).expanduser(),'r') as f:
|
||||
def get_cell_data_group(f: _h5py.File) -> str:
|
||||
base_group = DREAM3D_base_group(f)
|
||||
cells = tuple(f['/'.join([base_group,'_SIMPL_GEOMETRY','DIMENSIONS'])][()][::-1])
|
||||
cell_data_group = f[base_group].visititems(lambda path,obj: path.split('/')[0] \
|
||||
if isinstance(obj,_h5py._hl.dataset.Dataset) and _np.shape(obj)[:-1] == cells \
|
||||
else None)
|
||||
if cell_data_group is None:
|
||||
raise ValueError(f'could not determine cell-data group in file "{fname}/{base_group}"')
|
||||
return cell_data_group
|
||||
|
||||
if cell_data_group is None:
|
||||
raise ValueError(f'could not determine cell-data group in file "{fname}/{base_group}"')
|
||||
if isinstance(fname,_h5py.File):
|
||||
return get_cell_data_group(fname)
|
||||
|
||||
return cell_data_group
|
||||
with _h5py.File(_Path(fname).expanduser(),'r') as f:
|
||||
return get_cell_data_group(f)
|
||||
|
||||
|
||||
def Bravais_to_Miller(*,
|
||||
|
|
|
@ -146,21 +146,25 @@ class TestUtil:
|
|||
assert 'DAMASK' in style('DAMASK')
|
||||
|
||||
@pytest.mark.parametrize('complete',[True,False])
|
||||
def test_D3D_base_group(self,tmp_path,complete):
|
||||
@pytest.mark.parametrize('fhandle',[True,False])
|
||||
def test_D3D_base_group(self,tmp_path,complete,fhandle):
|
||||
base_group = ''.join(random.choices('DAMASK', k=10))
|
||||
with h5py.File(tmp_path/'base_group.dream3d','w') as f:
|
||||
f.create_group('/'.join((base_group,'_SIMPL_GEOMETRY')))
|
||||
if complete:
|
||||
f['/'.join((base_group,'_SIMPL_GEOMETRY'))].create_dataset('SPACING',data=np.ones(3))
|
||||
|
||||
fname = tmp_path/'base_group.dream3d'
|
||||
if fhandle: fname = h5py.File(fname)
|
||||
if complete:
|
||||
assert base_group == util.DREAM3D_base_group(tmp_path/'base_group.dream3d')
|
||||
assert base_group == util.DREAM3D_base_group(fname)
|
||||
else:
|
||||
with pytest.raises(ValueError):
|
||||
util.DREAM3D_base_group(tmp_path/'base_group.dream3d')
|
||||
util.DREAM3D_base_group(fname)
|
||||
|
||||
@pytest.mark.parametrize('complete',[True,False])
|
||||
def test_D3D_cell_data_group(self,tmp_path,complete):
|
||||
@pytest.mark.parametrize('fhandle',[True,False])
|
||||
def test_D3D_cell_data_group(self,tmp_path,complete,fhandle):
|
||||
base_group = ''.join(random.choices('DAMASK', k=10))
|
||||
cell_data_group = ''.join(random.choices('KULeuven', k=10))
|
||||
cells = np.random.randint(1,50,3)
|
||||
|
@ -172,11 +176,13 @@ class TestUtil:
|
|||
if complete:
|
||||
f['/'.join((base_group,cell_data_group))].create_dataset('data',shape=np.append(cells,1))
|
||||
|
||||
fname = tmp_path/'cell_data_group.dream3d'
|
||||
if fhandle: fname = h5py.File(fname)
|
||||
if complete:
|
||||
assert cell_data_group == util.DREAM3D_cell_data_group(tmp_path/'cell_data_group.dream3d')
|
||||
assert cell_data_group == util.DREAM3D_cell_data_group(fname)
|
||||
else:
|
||||
with pytest.raises(ValueError):
|
||||
util.DREAM3D_cell_data_group(tmp_path/'cell_data_group.dream3d')
|
||||
util.DREAM3D_cell_data_group(fname)
|
||||
|
||||
|
||||
@pytest.mark.parametrize('full,reduced',[({}, {}),
|
||||
|
|
Loading…
Reference in New Issue