2020-04-21 01:56:29 +05:30
|
|
|
|
"""
|
|
|
|
|
Filters for operations on regular grids.
|
|
|
|
|
|
|
|
|
|
The grids are defined as (x,y,z,...) where x is fastest and z is slowest.
|
2021-07-16 13:51:06 +05:30
|
|
|
|
This convention is consistent with the layout in grid vti files.
|
2021-04-23 22:50:07 +05:30
|
|
|
|
|
2020-04-21 01:56:29 +05:30
|
|
|
|
When converting to/from a plain list (e.g. storage in ASCII table),
|
2020-04-21 02:17:19 +05:30
|
|
|
|
the following operations are required for tensorial data:
|
2020-04-21 01:56:29 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
- D3 = D1.reshape(cells+(-1,),order='F').reshape(cells+(3,3))
|
|
|
|
|
- D1 = D3.reshape(cells+(-1,)).reshape(-1,9,order='F')
|
2020-04-21 01:56:29 +05:30
|
|
|
|
|
|
|
|
|
"""
|
2021-11-03 12:23:38 +05:30
|
|
|
|
|
2022-01-12 18:48:38 +05:30
|
|
|
|
from typing import Tuple as _Tuple
|
2021-11-02 23:58:54 +05:30
|
|
|
|
|
2021-11-03 12:23:38 +05:30
|
|
|
|
from scipy import spatial as _spatial
|
2020-04-10 16:00:39 +05:30
|
|
|
|
import numpy as _np
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2022-01-12 18:48:38 +05:30
|
|
|
|
from ._typehints import FloatSequence as _FloatSequence, IntSequence as _IntSequence
|
2021-11-03 12:23:38 +05:30
|
|
|
|
|
2022-01-12 18:48:38 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def _ks(size: _FloatSequence,
|
|
|
|
|
cells: _IntSequence,
|
|
|
|
|
first_order: bool = False) -> _np.ndarray:
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
|
|
|
|
Get wave numbers operator.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
cells : sequence of int, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Number of cells.
|
2020-04-21 02:17:19 +05:30
|
|
|
|
first_order : bool, optional
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Correction for first order derivatives, defaults to False.
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2021-04-23 22:50:07 +05:30
|
|
|
|
k_sk = _np.where(_np.arange(cells[0])>cells[0]//2,
|
|
|
|
|
_np.arange(cells[0])-cells[0],_np.arange(cells[0]))/size[0]
|
2020-12-04 02:28:24 +05:30
|
|
|
|
if cells[0]%2 == 0 and first_order: k_sk[cells[0]//2] = 0 # Nyquist freq=0 for even cells (Johnson, MIT, 2011)
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
k_sj = _np.where(_np.arange(cells[1])>cells[1]//2,
|
|
|
|
|
_np.arange(cells[1])-cells[1],_np.arange(cells[1]))/size[1]
|
2020-12-04 02:28:24 +05:30
|
|
|
|
if cells[1]%2 == 0 and first_order: k_sj[cells[1]//2] = 0 # Nyquist freq=0 for even cells (Johnson, MIT, 2011)
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2020-12-04 02:28:24 +05:30
|
|
|
|
k_si = _np.arange(cells[2]//2+1)/size[2]
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2020-04-20 16:11:03 +05:30
|
|
|
|
return _np.stack(_np.meshgrid(k_sk,k_sj,k_si,indexing = 'ij'), axis=-1)
|
2019-11-28 20:16:22 +05:30
|
|
|
|
|
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def curl(size: _FloatSequence,
|
|
|
|
|
f: _np.ndarray) -> _np.ndarray:
|
2021-04-23 22:50:07 +05:30
|
|
|
|
u"""
|
2019-12-09 00:50:13 +05:30
|
|
|
|
Calculate curl of a vector or tensor field in Fourier space.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
f : numpy.ndarray, shape (:,:,:,3) or (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Periodic field of which the curl is calculated.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
∇ × f : numpy.ndarray, shape (:,:,:,3) or (:,:,:,3,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Curl of f.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2021-04-23 22:50:07 +05:30
|
|
|
|
n = _np.prod(f.shape[3:])
|
|
|
|
|
k_s = _ks(size,f.shape[:3],True)
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
|
e = _np.zeros((3, 3, 3))
|
2020-02-21 03:20:54 +05:30
|
|
|
|
e[0, 1, 2] = e[1, 2, 0] = e[2, 0, 1] = +1.0 # Levi-Civita symbol
|
2019-11-26 14:55:39 +05:30
|
|
|
|
e[0, 2, 1] = e[2, 1, 0] = e[1, 0, 2] = -1.0
|
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
f_fourier = _np.fft.rfftn(f,axes=(0,1,2))
|
2022-11-30 02:31:14 +05:30
|
|
|
|
curl_ = (_np.einsum('slm,ijkl,ijkm ->ijks' if n == 3 else
|
|
|
|
|
'slm,ijkl,ijknm->ijksn',e,k_s,f_fourier)*2.0j*_np.pi) # vector 3->3, tensor 3x3->3x3
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
return _np.fft.irfftn(curl_,axes=(0,1,2),s=f.shape[:3])
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def divergence(size: _FloatSequence,
|
|
|
|
|
f: _np.ndarray) -> _np.ndarray:
|
2021-04-23 22:50:07 +05:30
|
|
|
|
u"""
|
2019-12-09 00:50:13 +05:30
|
|
|
|
Calculate divergence of a vector or tensor field in Fourier space.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
f : numpy.ndarray, shape (:,:,:,3) or (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Periodic field of which the divergence is calculated.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
∇ · f : numpy.ndarray, shape (:,:,:,1) or (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Divergence of f.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2021-04-23 22:50:07 +05:30
|
|
|
|
n = _np.prod(f.shape[3:])
|
|
|
|
|
k_s = _ks(size,f.shape[:3],True)
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
f_fourier = _np.fft.rfftn(f,axes=(0,1,2))
|
2022-12-05 21:30:23 +05:30
|
|
|
|
divergence_ = (_np.einsum('ijkl,ijkl ->ijk' if n == 3 else
|
2022-11-30 02:31:14 +05:30
|
|
|
|
'ijkm,ijklm->ijkl', k_s,f_fourier)*2.0j*_np.pi) # vector 3->1, tensor 3x3->3
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2022-11-26 16:28:36 +05:30
|
|
|
|
return _np.fft.irfftn(divergence_,axes=(0,1,2),s=f.shape[:3])
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def gradient(size: _FloatSequence,
|
|
|
|
|
f: _np.ndarray) -> _np.ndarray:
|
2021-04-23 22:50:07 +05:30
|
|
|
|
u"""
|
2022-01-12 18:48:38 +05:30
|
|
|
|
Calculate gradient of a scalar or vector field in Fourier space.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
f : numpy.ndarray, shape (:,:,:,1) or (:,:,:,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Periodic field of which the gradient is calculated.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
∇ f : numpy.ndarray, shape (:,:,:,3) or (:,:,:,3,3)
|
2022-11-26 16:28:36 +05:30
|
|
|
|
Gradient of f.
|
2021-04-23 22:50:07 +05:30
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2021-04-23 22:50:07 +05:30
|
|
|
|
n = _np.prod(f.shape[3:])
|
|
|
|
|
k_s = _ks(size,f.shape[:3],True)
|
2019-11-28 10:11:53 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
f_fourier = _np.fft.rfftn(f,axes=(0,1,2))
|
2022-12-05 21:30:23 +05:30
|
|
|
|
gradient_ = (_np.einsum('ijkl,ijkm->ijkm' if n == 1 else
|
2022-11-30 02:31:14 +05:30
|
|
|
|
'ijkl,ijkm->ijklm',f_fourier,k_s)*2.0j*_np.pi) # scalar 1->3, vector 3->3x3
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2022-11-26 16:28:36 +05:30
|
|
|
|
return _np.fft.irfftn(gradient_,axes=(0,1,2),s=f.shape[:3])
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
|
|
|
|
|
2022-01-12 18:48:38 +05:30
|
|
|
|
def coordinates0_point(cells: _IntSequence,
|
|
|
|
|
size: _FloatSequence,
|
|
|
|
|
origin: _FloatSequence = _np.zeros(3)) -> _np.ndarray:
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
|
|
|
|
Cell center positions (undeformed).
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
cells : sequence of int, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Number of cells.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
origin : sequence of float, len(3), optional
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
x_p_0 : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Undeformed cell center coordinates.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size_ = _np.array(size,float)
|
2022-03-27 02:36:03 +05:30
|
|
|
|
start = origin + size_/_np.array(cells,_np.int64)*.5
|
|
|
|
|
end = origin + size_ - size_/_np.array(cells,_np.int64)*.5
|
2020-04-20 16:11:03 +05:30
|
|
|
|
|
2020-12-04 02:28:24 +05:30
|
|
|
|
return _np.stack(_np.meshgrid(_np.linspace(start[0],end[0],cells[0]),
|
|
|
|
|
_np.linspace(start[1],end[1],cells[1]),
|
|
|
|
|
_np.linspace(start[2],end[2],cells[2]),indexing = 'ij'),
|
2020-04-20 16:11:03 +05:30
|
|
|
|
axis = -1)
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def displacement_fluct_point(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray) -> _np.ndarray:
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
|
|
|
|
Cell center displacement field from fluctuation part of the deformation gradient field.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
u_p_fluct : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Fluctuating part of the cell center displacements.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2020-02-21 03:20:54 +05:30
|
|
|
|
k_s = _ks(size,F.shape[:3],False)
|
2020-04-10 16:00:39 +05:30
|
|
|
|
k_s_squared = _np.einsum('...l,...l',k_s,k_s)
|
2019-12-03 23:29:59 +05:30
|
|
|
|
k_s_squared[0,0,0] = 1.0
|
2019-11-28 22:52:34 +05:30
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
|
displacement = -_np.einsum('ijkml,ijkl,l->ijkm',
|
|
|
|
|
_np.fft.rfftn(F,axes=(0,1,2)),
|
2019-11-28 22:52:34 +05:30
|
|
|
|
k_s,
|
2023-06-05 03:16:35 +05:30
|
|
|
|
_np.array([0.5j/_np.pi]*3),
|
|
|
|
|
) / k_s_squared[...,_np.newaxis]
|
2019-12-03 23:29:59 +05:30
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
|
return _np.fft.irfftn(displacement,axes=(0,1,2),s=F.shape[:3])
|
2019-12-03 23:29:59 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def displacement_avg_point(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray) -> _np.ndarray:
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
|
|
|
|
Cell center displacement field from average part of the deformation gradient field.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
u_p_avg : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Average part of the cell center displacements.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
|
F_avg = _np.average(F,axis=(0,1,2))
|
2020-12-04 03:30:49 +05:30
|
|
|
|
return _np.einsum('ml,ijkl->ijkm',F_avg - _np.eye(3),coordinates0_point(F.shape[:3],size))
|
2019-12-03 23:29:59 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def displacement_point(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray) -> _np.ndarray:
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
|
|
|
|
Cell center displacement field from deformation gradient field.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
u_p : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Cell center displacements.
|
|
|
|
|
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
2020-12-04 03:30:49 +05:30
|
|
|
|
return displacement_avg_point(size,F) + displacement_fluct_point(size,F)
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def coordinates_point(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray,
|
|
|
|
|
origin: _FloatSequence = _np.zeros(3)) -> _np.ndarray:
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
|
|
|
|
Cell center positions.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
origin : sequence of float, len(3), optional
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
x_p : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Cell center coordinates.
|
|
|
|
|
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
2020-12-04 03:30:49 +05:30
|
|
|
|
return coordinates0_point(F.shape[:3],size,origin) + displacement_point(size,F)
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2021-11-02 23:58:54 +05:30
|
|
|
|
def cellsSizeOrigin_coordinates0_point(coordinates0: _np.ndarray,
|
2022-01-12 18:48:38 +05:30
|
|
|
|
ordered: bool = True) -> _Tuple[_np.ndarray,_np.ndarray,_np.ndarray]:
|
2019-12-08 23:03:43 +05:30
|
|
|
|
"""
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Return grid 'DNA', i.e. cells, size, and origin from 1D array of point positions.
|
2019-12-08 23:03:43 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
coordinates0 : numpy.ndarray, shape (:,3)
|
|
|
|
|
Undeformed cell center coordinates.
|
2019-12-08 23:03:43 +05:30
|
|
|
|
ordered : bool, optional
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Expect coordinates0 data to be ordered (x fast, z slow).
|
2020-12-07 22:19:37 +05:30
|
|
|
|
Defaults to True.
|
2019-12-08 23:03:43 +05:30
|
|
|
|
|
2021-04-24 18:17:52 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2021-04-24 21:30:57 +05:30
|
|
|
|
cells, size, origin : Three numpy.ndarray, each of shape (3)
|
|
|
|
|
Information to reconstruct grid.
|
2021-04-24 18:17:52 +05:30
|
|
|
|
|
2019-12-08 23:03:43 +05:30
|
|
|
|
"""
|
2020-12-04 03:30:49 +05:30
|
|
|
|
coords = [_np.unique(coordinates0[:,i]) for i in range(3)]
|
2020-04-10 16:00:39 +05:30
|
|
|
|
mincorner = _np.array(list(map(min,coords)))
|
|
|
|
|
maxcorner = _np.array(list(map(max,coords)))
|
2022-03-27 02:36:03 +05:30
|
|
|
|
cells = _np.array(list(map(len,coords)),_np.int64)
|
2020-12-04 02:28:24 +05:30
|
|
|
|
size = cells/_np.maximum(cells-1,1) * (maxcorner-mincorner)
|
|
|
|
|
delta = size/cells
|
2019-12-08 22:27:02 +05:30
|
|
|
|
origin = mincorner - delta*.5
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2019-12-24 20:29:09 +05:30
|
|
|
|
# 1D/2D: size/origin combination undefined, set origin to 0.0
|
2022-12-05 21:30:23 +05:30
|
|
|
|
size [_np.where(cells == 1)] = origin[_np.where(cells == 1)]*2.
|
|
|
|
|
origin[_np.where(cells == 1)] = 0.0
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2020-12-04 03:30:49 +05:30
|
|
|
|
if cells.prod() != len(coordinates0):
|
2022-02-22 21:12:05 +05:30
|
|
|
|
raise ValueError(f'data count {len(coordinates0)} does not match cells {cells}')
|
2019-12-08 22:27:02 +05:30
|
|
|
|
|
|
|
|
|
start = origin + delta*.5
|
2019-12-08 22:43:20 +05:30
|
|
|
|
end = origin - delta*.5 + size
|
2020-04-24 23:31:57 +05:30
|
|
|
|
|
|
|
|
|
atol = _np.max(size)*5e-2
|
2020-12-04 02:28:24 +05:30
|
|
|
|
if not (_np.allclose(coords[0],_np.linspace(start[0],end[0],cells[0]),atol=atol) and \
|
|
|
|
|
_np.allclose(coords[1],_np.linspace(start[1],end[1],cells[1]),atol=atol) and \
|
|
|
|
|
_np.allclose(coords[2],_np.linspace(start[2],end[2],cells[2]),atol=atol)):
|
2022-02-22 21:12:05 +05:30
|
|
|
|
raise ValueError('non-uniform cell spacing')
|
2019-12-08 22:27:02 +05:30
|
|
|
|
|
2020-12-04 03:30:49 +05:30
|
|
|
|
if ordered and not _np.allclose(coordinates0.reshape(tuple(cells)+(3,),order='F'),
|
2021-11-02 23:58:54 +05:30
|
|
|
|
coordinates0_point(list(cells),size,origin),atol=atol):
|
2022-02-22 21:12:05 +05:30
|
|
|
|
raise ValueError('input data is not ordered (x fast, z slow)')
|
2019-12-08 23:03:43 +05:30
|
|
|
|
|
2020-12-04 02:28:24 +05:30
|
|
|
|
return (cells,size,origin)
|
2019-12-08 22:27:02 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2022-01-12 18:48:38 +05:30
|
|
|
|
def coordinates0_node(cells: _IntSequence,
|
|
|
|
|
size: _FloatSequence,
|
|
|
|
|
origin: _FloatSequence = _np.zeros(3)) -> _np.ndarray:
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
|
|
|
|
Nodal positions (undeformed).
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
cells : sequence of int, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Number of cells.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
origin : sequence of float, len(3), optional
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
x_n_0 : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Undeformed nodal coordinates.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2020-12-04 02:28:24 +05:30
|
|
|
|
return _np.stack(_np.meshgrid(_np.linspace(origin[0],size[0]+origin[0],cells[0]+1),
|
|
|
|
|
_np.linspace(origin[1],size[1]+origin[1],cells[1]+1),
|
|
|
|
|
_np.linspace(origin[2],size[2]+origin[2],cells[2]+1),indexing = 'ij'),
|
2020-04-20 16:11:03 +05:30
|
|
|
|
axis = -1)
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2019-12-03 23:29:59 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def displacement_fluct_node(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray) -> _np.ndarray:
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
|
|
|
|
Nodal displacement field from fluctuation part of the deformation gradient field.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
u_n_fluct : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Fluctuating part of the nodal displacements.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2020-12-07 22:19:37 +05:30
|
|
|
|
return point_to_node(displacement_fluct_point(size,F))
|
2019-12-03 23:29:59 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def displacement_avg_node(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray) -> _np.ndarray:
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
|
|
|
|
Nodal displacement field from average part of the deformation gradient field.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
u_n_avg : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Average part of the nodal displacements.
|
|
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
|
F_avg = _np.average(F,axis=(0,1,2))
|
2020-12-04 03:30:49 +05:30
|
|
|
|
return _np.einsum('ml,ijkl->ijkm',F_avg - _np.eye(3),coordinates0_node(F.shape[:3],size))
|
2019-12-03 23:29:59 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def displacement_node(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray) -> _np.ndarray:
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
|
|
|
|
Nodal displacement field from deformation gradient field.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2023-02-21 20:57:06 +05:30
|
|
|
|
u_n : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Nodal displacements.
|
|
|
|
|
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
2020-12-04 03:30:49 +05:30
|
|
|
|
return displacement_avg_node(size,F) + displacement_fluct_node(size,F)
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def coordinates_node(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray,
|
|
|
|
|
origin: _FloatSequence = _np.zeros(3)) -> _np.ndarray:
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
|
|
|
|
Nodal positions.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical size of the periodic field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Deformation gradient field.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
origin : sequence of float, len(3), optional
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
x_n : numpy.ndarray, shape (:,:,:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Nodal coordinates.
|
|
|
|
|
|
2019-12-21 22:37:04 +05:30
|
|
|
|
"""
|
2020-12-04 03:30:49 +05:30
|
|
|
|
return coordinates0_node(F.shape[:3],size,origin) + displacement_node(size,F)
|
2019-11-28 22:52:34 +05:30
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2021-11-02 23:58:54 +05:30
|
|
|
|
def cellsSizeOrigin_coordinates0_node(coordinates0: _np.ndarray,
|
2022-01-12 18:48:38 +05:30
|
|
|
|
ordered: bool = True) -> _Tuple[_np.ndarray,_np.ndarray,_np.ndarray]:
|
2019-12-08 23:03:43 +05:30
|
|
|
|
"""
|
2020-12-04 02:28:24 +05:30
|
|
|
|
Return grid 'DNA', i.e. cells, size, and origin from 1D array of nodal positions.
|
2019-12-08 23:03:43 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
coordinates0 : numpy.ndarray, shape (:,3)
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Undeformed nodal coordinates.
|
2019-12-08 23:03:43 +05:30
|
|
|
|
ordered : bool, optional
|
2020-12-04 03:30:49 +05:30
|
|
|
|
Expect coordinates0 data to be ordered (x fast, z slow).
|
2020-12-07 22:19:37 +05:30
|
|
|
|
Defaults to True.
|
2019-12-08 23:03:43 +05:30
|
|
|
|
|
2021-04-24 18:17:52 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
2021-04-24 21:30:57 +05:30
|
|
|
|
cells, size, origin : Three numpy.ndarray, each of shape (3)
|
|
|
|
|
Information to reconstruct grid.
|
2021-04-24 18:17:52 +05:30
|
|
|
|
|
2019-12-08 23:03:43 +05:30
|
|
|
|
"""
|
2020-12-04 03:30:49 +05:30
|
|
|
|
coords = [_np.unique(coordinates0[:,i]) for i in range(3)]
|
2020-04-10 16:00:39 +05:30
|
|
|
|
mincorner = _np.array(list(map(min,coords)))
|
|
|
|
|
maxcorner = _np.array(list(map(max,coords)))
|
2022-03-27 02:36:03 +05:30
|
|
|
|
cells = _np.array(list(map(len,coords)),_np.int64) - 1
|
2019-12-08 22:27:02 +05:30
|
|
|
|
size = maxcorner-mincorner
|
|
|
|
|
origin = mincorner
|
2020-02-21 03:20:54 +05:30
|
|
|
|
|
2020-12-04 03:30:49 +05:30
|
|
|
|
if (cells+1).prod() != len(coordinates0):
|
2022-02-22 21:12:05 +05:30
|
|
|
|
raise ValueError(f'data count {len(coordinates0)} does not match cells {cells}')
|
2019-12-08 22:27:02 +05:30
|
|
|
|
|
2020-04-24 23:31:57 +05:30
|
|
|
|
atol = _np.max(size)*5e-2
|
2020-12-04 02:28:24 +05:30
|
|
|
|
if not (_np.allclose(coords[0],_np.linspace(mincorner[0],maxcorner[0],cells[0]+1),atol=atol) and \
|
|
|
|
|
_np.allclose(coords[1],_np.linspace(mincorner[1],maxcorner[1],cells[1]+1),atol=atol) and \
|
|
|
|
|
_np.allclose(coords[2],_np.linspace(mincorner[2],maxcorner[2],cells[2]+1),atol=atol)):
|
2022-02-22 21:12:05 +05:30
|
|
|
|
raise ValueError('non-uniform cell spacing')
|
2019-12-08 22:27:02 +05:30
|
|
|
|
|
2020-12-04 03:30:49 +05:30
|
|
|
|
if ordered and not _np.allclose(coordinates0.reshape(tuple(cells+1)+(3,),order='F'),
|
2021-11-02 23:58:54 +05:30
|
|
|
|
coordinates0_node(list(cells),size,origin),atol=atol):
|
2022-02-22 21:12:05 +05:30
|
|
|
|
raise ValueError('input data is not ordered (x fast, z slow)')
|
2019-12-08 23:03:43 +05:30
|
|
|
|
|
2020-12-04 02:28:24 +05:30
|
|
|
|
return (cells,size,origin)
|
2019-12-08 22:27:02 +05:30
|
|
|
|
|
|
|
|
|
|
2021-11-02 23:58:54 +05:30
|
|
|
|
def point_to_node(cell_data: _np.ndarray) -> _np.ndarray:
|
2021-04-23 22:50:07 +05:30
|
|
|
|
"""
|
|
|
|
|
Interpolate periodic point data to nodal data.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
cell_data : numpy.ndarray, shape (:,:,:,...)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Data defined on the cell centers of a periodic grid.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
node_data : numpy.ndarray, shape (:,:,:,...)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Data defined on the nodes of a periodic grid.
|
2021-04-24 18:17:52 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
"""
|
|
|
|
|
n = ( cell_data + _np.roll(cell_data,1,(0,1,2))
|
|
|
|
|
+ _np.roll(cell_data,1,(0,)) + _np.roll(cell_data,1,(1,)) + _np.roll(cell_data,1,(2,))
|
|
|
|
|
+ _np.roll(cell_data,1,(0,1)) + _np.roll(cell_data,1,(1,2)) + _np.roll(cell_data,1,(2,0)))*0.125
|
|
|
|
|
|
|
|
|
|
return _np.pad(n,((0,1),(0,1),(0,1))+((0,0),)*len(cell_data.shape[3:]),mode='wrap')
|
|
|
|
|
|
|
|
|
|
|
2021-11-02 23:58:54 +05:30
|
|
|
|
def node_to_point(node_data: _np.ndarray) -> _np.ndarray:
|
2021-04-23 22:50:07 +05:30
|
|
|
|
"""
|
|
|
|
|
Interpolate periodic nodal data to point data.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
node_data : numpy.ndarray, shape (:,:,:,...)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Data defined on the nodes of a periodic grid.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
cell_data : numpy.ndarray, shape (:,:,:,...)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Data defined on the cell centers of a periodic grid.
|
|
|
|
|
|
|
|
|
|
"""
|
|
|
|
|
c = ( node_data + _np.roll(node_data,1,(0,1,2))
|
|
|
|
|
+ _np.roll(node_data,1,(0,)) + _np.roll(node_data,1,(1,)) + _np.roll(node_data,1,(2,))
|
|
|
|
|
+ _np.roll(node_data,1,(0,1)) + _np.roll(node_data,1,(1,2)) + _np.roll(node_data,1,(2,0)))*0.125
|
|
|
|
|
|
|
|
|
|
return c[1:,1:,1:]
|
|
|
|
|
|
|
|
|
|
|
2021-11-02 23:58:54 +05:30
|
|
|
|
def coordinates0_valid(coordinates0: _np.ndarray) -> bool:
|
2021-04-23 22:50:07 +05:30
|
|
|
|
"""
|
2021-04-24 22:42:44 +05:30
|
|
|
|
Check whether coordinates form a regular grid.
|
2021-04-23 22:50:07 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
coordinates0 : numpy.ndarray, shape (:,3)
|
2021-04-23 22:50:07 +05:30
|
|
|
|
Array of undeformed cell coordinates.
|
|
|
|
|
|
2021-04-24 18:17:52 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
valid : bool
|
2021-04-24 22:42:44 +05:30
|
|
|
|
Whether the coordinates form a regular grid.
|
2021-04-24 18:17:52 +05:30
|
|
|
|
|
2021-04-23 22:50:07 +05:30
|
|
|
|
"""
|
|
|
|
|
try:
|
|
|
|
|
cellsSizeOrigin_coordinates0_point(coordinates0,ordered=True)
|
|
|
|
|
return True
|
|
|
|
|
except ValueError:
|
|
|
|
|
return False
|
|
|
|
|
|
|
|
|
|
|
2022-11-09 00:22:08 +05:30
|
|
|
|
def unravel_index(idx: _np.ndarray) -> _np.ndarray:
|
|
|
|
|
"""
|
|
|
|
|
Convert flat indices to coordinate indices.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
idx : numpy.ndarray, shape (:,:,:)
|
|
|
|
|
Grid of flat indices.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
unravelled : numpy.ndarray, shape (:,:,:,3)
|
|
|
|
|
Grid of coordinate indices.
|
|
|
|
|
|
|
|
|
|
Examples
|
|
|
|
|
--------
|
|
|
|
|
Unravel a linearly increasing sequence of material indices on a 3 × 2 × 1 grid.
|
|
|
|
|
|
|
|
|
|
>>> import numpy as np
|
|
|
|
|
>>> import damask
|
|
|
|
|
>>> seq = np.arange(6).reshape((3,2,1),order='F')
|
|
|
|
|
>>> (coord_idx := damask.grid_filters.unravel_index(seq))
|
|
|
|
|
array([[[[0, 0, 0]],
|
|
|
|
|
[[0, 1, 0]]],
|
|
|
|
|
[[[1, 0, 0]],
|
|
|
|
|
[[1, 1, 0]]],
|
|
|
|
|
[[[2, 0, 0]],
|
|
|
|
|
[[2, 1, 0]]]])
|
|
|
|
|
>>> coord_idx[1,1,0]
|
|
|
|
|
array([1, 1, 0])
|
|
|
|
|
|
|
|
|
|
"""
|
|
|
|
|
cells = idx.shape
|
|
|
|
|
idx_ = _np.expand_dims(idx,3)
|
|
|
|
|
return _np.block([ idx_ %cells[0],
|
|
|
|
|
(idx_//cells[0]) %cells[1],
|
|
|
|
|
((idx_//cells[0])//cells[1])%cells[2]])
|
|
|
|
|
|
|
|
|
|
def ravel_index(idx: _np.ndarray) -> _np.ndarray:
|
|
|
|
|
"""
|
|
|
|
|
Convert coordinate indices to flat indices.
|
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
|
|
|
|
idx : numpy.ndarray, shape (:,:,:,3)
|
|
|
|
|
Grid of coordinate indices.
|
|
|
|
|
|
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
ravelled : numpy.ndarray, shape (:,:,:)
|
|
|
|
|
Grid of flat indices.
|
|
|
|
|
|
|
|
|
|
Examples
|
|
|
|
|
--------
|
|
|
|
|
Ravel a reversed sequence of coordinate indices on a 2 × 2 × 1 grid.
|
|
|
|
|
|
|
|
|
|
>>> import numpy as np
|
|
|
|
|
>>> import damask
|
|
|
|
|
>>> (rev := np.array([[1,1,0],[0,1,0],[1,0,0],[0,0,0]]).reshape((2,2,1,3)))
|
|
|
|
|
array([[[[1, 1, 0]],
|
|
|
|
|
[[0, 1, 0]]],
|
|
|
|
|
[[[1, 0, 0]],
|
|
|
|
|
[[0, 0, 0]]]])
|
|
|
|
|
>>> (flat_idx := damask.grid_filters.ravel_index(rev))
|
|
|
|
|
array([[[3],
|
|
|
|
|
[2]],
|
|
|
|
|
[[1],
|
|
|
|
|
[0]]])
|
|
|
|
|
|
|
|
|
|
"""
|
|
|
|
|
cells = idx.shape[:3]
|
|
|
|
|
return idx[:,:,:,0] \
|
|
|
|
|
+ idx[:,:,:,1]*cells[0] \
|
|
|
|
|
+ idx[:,:,:,2]*cells[0]*cells[1]
|
|
|
|
|
|
|
|
|
|
|
2022-01-31 03:06:30 +05:30
|
|
|
|
def regrid(size: _FloatSequence,
|
|
|
|
|
F: _np.ndarray,
|
|
|
|
|
cells: _IntSequence) -> _np.ndarray:
|
2020-01-23 21:45:02 +05:30
|
|
|
|
"""
|
2022-11-09 00:22:08 +05:30
|
|
|
|
Map a deformed grid A back to a rectilinear grid B.
|
|
|
|
|
|
|
|
|
|
The size of grid B is chosen as the average deformed size of grid A.
|
2020-01-23 21:45:02 +05:30
|
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
|
----------
|
2022-01-12 18:48:38 +05:30
|
|
|
|
size : sequence of float, len (3)
|
2022-11-09 00:22:08 +05:30
|
|
|
|
Physical size of grid A.
|
|
|
|
|
F : numpy.ndarray, shape (:,:,:,3,3)
|
|
|
|
|
Deformation gradient field on grid A.
|
2022-01-12 18:48:38 +05:30
|
|
|
|
cells : sequence of int, len (3)
|
2022-11-09 00:22:08 +05:30
|
|
|
|
Cell count along x,y,z of grid B.
|
2020-01-23 21:45:02 +05:30
|
|
|
|
|
2022-11-09 00:22:08 +05:30
|
|
|
|
Returns
|
|
|
|
|
-------
|
|
|
|
|
idx : numpy.ndarray of int, shape (cells)
|
|
|
|
|
Flat index of closest point on deformed grid A for each point on grid B.
|
2019-12-04 14:47:55 +05:30
|
|
|
|
|
2022-11-09 00:22:08 +05:30
|
|
|
|
"""
|
|
|
|
|
box = _np.dot(_np.average(F,axis=(0,1,2)),size)
|
|
|
|
|
c = coordinates_point(size,F)%box
|
|
|
|
|
tree = _spatial.cKDTree(c.reshape((-1,3),order='F'),boxsize=box)
|
|
|
|
|
return tree.query(coordinates0_point(cells,box))[1]
|