2020-04-10 16:00:39 +05:30
|
|
|
from scipy import spatial as _spatial
|
|
|
|
import numpy as _np
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
def _ks(size,grid,first_order=False):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Get wave numbers operator.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
k_sk = _np.where(_np.arange(grid[0])>grid[0]//2,_np.arange(grid[0])-grid[0],_np.arange(grid[0]))/size[0]
|
2019-11-28 22:52:34 +05:30
|
|
|
if grid[0]%2 == 0 and first_order: k_sk[grid[0]//2] = 0 # Nyquist freq=0 for even grid (Johnson, MIT, 2011)
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
k_sj = _np.where(_np.arange(grid[1])>grid[1]//2,_np.arange(grid[1])-grid[1],_np.arange(grid[1]))/size[1]
|
2019-11-28 22:52:34 +05:30
|
|
|
if grid[1]%2 == 0 and first_order: k_sj[grid[1]//2] = 0 # Nyquist freq=0 for even grid (Johnson, MIT, 2011)
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
k_si = _np.arange(grid[2]//2+1)/size[2]
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
kk, kj, ki = _np.meshgrid(k_sk,k_sj,k_si,indexing = 'ij')
|
|
|
|
return _np.concatenate((ki[:,:,:,None],kj[:,:,:,None],kk[:,:,:,None]),axis = 3)
|
2019-11-28 20:16:22 +05:30
|
|
|
|
|
|
|
|
|
|
|
def curl(size,field):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Calculate curl of a vector or tensor field in Fourier space.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
n = _np.prod(field.shape[3:])
|
2020-02-21 03:20:54 +05:30
|
|
|
k_s = _ks(size,field.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
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
field_fourier = _np.fft.rfftn(field,axes=(0,1,2))
|
2020-04-14 14:43:07 +05:30
|
|
|
curl_ = (_np.einsum('slm,ijkl,ijkm ->ijks', e,k_s,field_fourier)*2.0j*_np.pi if n == 3 else # vector, 3 -> 3
|
|
|
|
_np.einsum('slm,ijkl,ijknm->ijksn',e,k_s,field_fourier)*2.0j*_np.pi) # tensor, 3x3 -> 3x3
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
return _np.fft.irfftn(curl_,axes=(0,1,2),s=field.shape[:3])
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
|
|
|
|
|
|
def divergence(size,field):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Calculate divergence of a vector or tensor field in Fourier space.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
n = _np.prod(field.shape[3:])
|
2020-02-21 03:20:54 +05:30
|
|
|
k_s = _ks(size,field.shape[:3],True)
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
field_fourier = _np.fft.rfftn(field,axes=(0,1,2))
|
2020-04-14 14:43:07 +05:30
|
|
|
div_ = (_np.einsum('ijkl,ijkl ->ijk', k_s,field_fourier)*2.0j*_np.pi if n == 3 else # vector, 3 -> 1
|
|
|
|
_np.einsum('ijkm,ijklm->ijkl',k_s,field_fourier)*2.0j*_np.pi) # tensor, 3x3 -> 3
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
return _np.fft.irfftn(div_,axes=(0,1,2),s=field.shape[:3])
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
|
|
|
|
|
|
def gradient(size,field):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Calculate gradient of a vector or scalar field in Fourier space.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
n = _np.prod(field.shape[3:])
|
2020-02-21 03:20:54 +05:30
|
|
|
k_s = _ks(size,field.shape[:3],True)
|
2019-11-28 10:11:53 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
field_fourier = _np.fft.rfftn(field,axes=(0,1,2))
|
2020-04-14 14:43:07 +05:30
|
|
|
grad_ = (_np.einsum('ijkl,ijkm->ijkm', field_fourier,k_s)*2.0j*_np.pi if n == 1 else # scalar, 1 -> 3
|
|
|
|
_np.einsum('ijkl,ijkm->ijklm',field_fourier,k_s)*2.0j*_np.pi) # vector, 3 -> 3x3
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
return _np.fft.irfftn(grad_,axes=(0,1,2),s=field.shape[:3])
|
2019-11-26 14:55:39 +05:30
|
|
|
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
def cell_coord0(grid,size,origin=_np.zeros(3)):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Cell center positions (undeformed).
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
grid : numpy.ndarray
|
|
|
|
number of grid points.
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-21 22:37:04 +05:30
|
|
|
origin : numpy.ndarray, optional
|
2020-03-21 04:17:54 +05:30
|
|
|
physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
|
|
"""
|
2020-03-21 04:17:54 +05:30
|
|
|
start = origin + size/grid*.5
|
|
|
|
end = origin + size - size/grid*.5
|
2020-04-10 16:00:39 +05:30
|
|
|
return _np.mgrid[start[0]:end[0]:grid[0]*1j,start[1]:end[1]:grid[1]*1j,start[2]:end[2]:grid[2]*1j].T
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-11-26 14:55:39 +05:30
|
|
|
|
2019-12-05 23:02:21 +05:30
|
|
|
def cell_displacement_fluct(size,F):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Cell center displacement field from fluctuation part of the deformation gradient field.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-09 00:50:13 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
integrator = 0.5j*size/_np.pi
|
2019-11-28 22:52:34 +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,
|
|
|
|
integrator,
|
2020-04-10 16:00:39 +05:30
|
|
|
) / 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
|
|
|
|
2019-12-05 23:02:21 +05:30
|
|
|
def cell_displacement_avg(size,F):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Cell center displacement field from average part of the deformation gradient field.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-09 00:50:13 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
F_avg = _np.average(F,axis=(0,1,2))
|
|
|
|
return _np.einsum('ml,ijkl->ijkm',F_avg-_np.eye(3),cell_coord0(F.shape[:3][::-1],size))
|
2019-12-03 23:29:59 +05:30
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-21 22:37:04 +05:30
|
|
|
def cell_displacement(size,F):
|
|
|
|
"""
|
|
|
|
Cell center displacement field from deformation gradient field.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-21 22:37:04 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
|
|
|
|
"""
|
|
|
|
return cell_displacement_avg(size,F) + cell_displacement_fluct(size,F)
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
def cell_coord(size,F,origin=_np.zeros(3)):
|
2019-12-21 22:37:04 +05:30
|
|
|
"""
|
|
|
|
Cell center positions.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-21 22:37:04 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
origin : numpy.ndarray, optional
|
2020-03-21 04:17:54 +05:30
|
|
|
physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
|
|
"""
|
|
|
|
return cell_coord0(F.shape[:3][::-1],size,origin) + cell_displacement(size,F)
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2020-01-13 07:21:49 +05:30
|
|
|
def cell_coord0_gridSizeOrigin(coord0,ordered=True):
|
2019-12-08 23:03:43 +05:30
|
|
|
"""
|
|
|
|
Return grid 'DNA', i.e. grid, size, and origin from array of cell positions.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
coord0 : numpy.ndarray
|
2019-12-09 00:50:13 +05:30
|
|
|
array of undeformed cell coordinates.
|
2019-12-08 23:03:43 +05:30
|
|
|
ordered : bool, optional
|
|
|
|
expect coord0 data to be ordered (x fast, z slow).
|
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
coords = [_np.unique(coord0[:,i]) for i in range(3)]
|
|
|
|
mincorner = _np.array(list(map(min,coords)))
|
|
|
|
maxcorner = _np.array(list(map(max,coords)))
|
|
|
|
grid = _np.array(list(map(len,coords)),'i')
|
|
|
|
size = grid/_np.maximum(grid-1,1) * (maxcorner-mincorner)
|
2019-12-08 22:27:02 +05:30
|
|
|
delta = size/grid
|
|
|
|
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
|
2020-04-10 16:00:39 +05:30
|
|
|
size [_np.where(grid==1)] = origin[_np.where(grid==1)]*2.
|
|
|
|
origin[_np.where(grid==1)] = 0.0
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-08 22:27:02 +05:30
|
|
|
if grid.prod() != len(coord0):
|
|
|
|
raise ValueError('Data count {} does not match grid {}.'.format(len(coord0),grid))
|
|
|
|
|
|
|
|
start = origin + delta*.5
|
2019-12-08 22:43:20 +05:30
|
|
|
end = origin - delta*.5 + size
|
2019-12-08 22:27:02 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
if not _np.allclose(coords[0],_np.linspace(start[0],end[0],grid[0])) and \
|
|
|
|
_np.allclose(coords[1],_np.linspace(start[1],end[1],grid[1])) and \
|
|
|
|
_np.allclose(coords[2],_np.linspace(start[2],end[2],grid[2])):
|
2019-12-08 22:27:02 +05:30
|
|
|
raise ValueError('Regular grid spacing violated.')
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
if ordered and not _np.allclose(coord0.reshape(tuple(grid[::-1])+(3,)),cell_coord0(grid,size,origin)):
|
2020-04-14 14:43:07 +05:30
|
|
|
raise ValueError('Input data is not a regular grid.')
|
2019-12-08 23:03:43 +05:30
|
|
|
|
2019-12-08 22:27:02 +05:30
|
|
|
return (grid,size,origin)
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-21 23:33:36 +05:30
|
|
|
def coord0_check(coord0):
|
|
|
|
"""
|
2019-12-21 23:34:29 +05:30
|
|
|
Check whether coordinates lie on a regular grid.
|
2019-12-21 23:33:36 +05:30
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
coord0 : numpy.ndarray
|
|
|
|
array of undeformed cell coordinates.
|
|
|
|
|
|
|
|
"""
|
2020-01-13 07:21:49 +05:30
|
|
|
cell_coord0_gridSizeOrigin(coord0,ordered=True)
|
2019-12-21 23:33:36 +05:30
|
|
|
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
def node_coord0(grid,size,origin=_np.zeros(3)):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Nodal positions (undeformed).
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
grid : numpy.ndarray
|
|
|
|
number of grid points.
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-21 22:37:04 +05:30
|
|
|
origin : numpy.ndarray, optional
|
2020-03-21 04:17:54 +05:30
|
|
|
physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-09 00:50:13 +05:30
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
return _np.mgrid[origin[0]:size[0]+origin[0]:(grid[0]+1)*1j,
|
2020-03-21 23:34:40 +05:30
|
|
|
origin[1]:size[1]+origin[1]:(grid[1]+1)*1j,
|
|
|
|
origin[2]:size[2]+origin[2]:(grid[2]+1)*1j].T
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-03 23:29:59 +05:30
|
|
|
|
2019-12-05 23:02:21 +05:30
|
|
|
def node_displacement_fluct(size,F):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Nodal displacement field from fluctuation part of the deformation gradient field.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-09 00:50:13 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
|
|
|
|
"""
|
2019-12-05 23:02:21 +05:30
|
|
|
return cell_2_node(cell_displacement_fluct(size,F))
|
2019-12-03 23:29:59 +05:30
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-05 23:02:21 +05:30
|
|
|
def node_displacement_avg(size,F):
|
2019-12-09 00:50:13 +05:30
|
|
|
"""
|
|
|
|
Nodal displacement field from average part of the deformation gradient field.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-09 00:50:13 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
F_avg = _np.average(F,axis=(0,1,2))
|
|
|
|
return _np.einsum('ml,ijkl->ijkm',F_avg-_np.eye(3),node_coord0(F.shape[:3][::-1],size))
|
2019-12-03 23:29:59 +05:30
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-21 22:37:04 +05:30
|
|
|
def node_displacement(size,F):
|
|
|
|
"""
|
|
|
|
Nodal displacement field from deformation gradient field.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-21 22:37:04 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
|
|
|
|
"""
|
|
|
|
return node_displacement_avg(size,F) + node_displacement_fluct(size,F)
|
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
def node_coord(size,F,origin=_np.zeros(3)):
|
2019-12-21 22:37:04 +05:30
|
|
|
"""
|
|
|
|
Nodal positions.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
2020-02-21 03:20:54 +05:30
|
|
|
physical size of the periodic field.
|
2019-12-21 22:37:04 +05:30
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field.
|
|
|
|
origin : numpy.ndarray, optional
|
2020-03-21 04:17:54 +05:30
|
|
|
physical origin of the periodic field. Defaults to [0.0,0.0,0.0].
|
2019-12-21 22:37:04 +05:30
|
|
|
|
|
|
|
"""
|
|
|
|
return node_coord0(F.shape[:3][::-1],size,origin) + node_displacement(size,F)
|
2019-11-28 22:52:34 +05:30
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-03 23:29:59 +05:30
|
|
|
def cell_2_node(cell_data):
|
2019-12-21 22:37:04 +05:30
|
|
|
"""Interpolate periodic cell data to nodal data."""
|
2020-04-10 16:00:39 +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
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
return _np.pad(n,((0,1),(0,1),(0,1))+((0,0),)*len(cell_data.shape[3:]),mode='wrap')
|
2019-12-04 00:00:55 +05:30
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-04 00:00:55 +05:30
|
|
|
def node_2_cell(node_data):
|
2019-12-21 22:37:04 +05:30
|
|
|
"""Interpolate periodic nodal data to cell data."""
|
2020-04-10 16:00:39 +05:30
|
|
|
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
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-04 00:00:55 +05:30
|
|
|
return c[:-1,:-1,:-1]
|
2019-12-04 14:47:55 +05:30
|
|
|
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2020-01-13 07:21:49 +05:30
|
|
|
def node_coord0_gridSizeOrigin(coord0,ordered=False):
|
2019-12-08 23:03:43 +05:30
|
|
|
"""
|
|
|
|
Return grid 'DNA', i.e. grid, size, and origin from array of nodal positions.
|
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
coord0 : numpy.ndarray
|
2020-01-23 21:45:02 +05:30
|
|
|
array of undeformed nodal coordinates.
|
2019-12-08 23:03:43 +05:30
|
|
|
ordered : bool, optional
|
|
|
|
expect coord0 data to be ordered (x fast, z slow).
|
|
|
|
|
|
|
|
"""
|
2020-04-10 16:00:39 +05:30
|
|
|
coords = [_np.unique(coord0[:,i]) for i in range(3)]
|
|
|
|
mincorner = _np.array(list(map(min,coords)))
|
|
|
|
maxcorner = _np.array(list(map(max,coords)))
|
|
|
|
grid = _np.array(list(map(len,coords)),'i') - 1
|
2019-12-08 22:27:02 +05:30
|
|
|
size = maxcorner-mincorner
|
|
|
|
origin = mincorner
|
2020-02-21 03:20:54 +05:30
|
|
|
|
2019-12-08 22:27:02 +05:30
|
|
|
if (grid+1).prod() != len(coord0):
|
|
|
|
raise ValueError('Data count {} does not match grid {}.'.format(len(coord0),grid))
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
if not _np.allclose(coords[0],_np.linspace(mincorner[0],maxcorner[0],grid[0]+1)) and \
|
|
|
|
_np.allclose(coords[1],_np.linspace(mincorner[1],maxcorner[1],grid[1]+1)) and \
|
|
|
|
_np.allclose(coords[2],_np.linspace(mincorner[2],maxcorner[2],grid[2]+1)):
|
2019-12-08 22:27:02 +05:30
|
|
|
raise ValueError('Regular grid spacing violated.')
|
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
if ordered and not _np.allclose(coord0.reshape(tuple((grid+1)[::-1])+(3,)),node_coord0(grid,size,origin)):
|
2020-04-14 14:43:07 +05:30
|
|
|
raise ValueError('Input data is not a regular grid.')
|
2019-12-08 23:03:43 +05:30
|
|
|
|
2019-12-08 22:27:02 +05:30
|
|
|
return (grid,size,origin)
|
|
|
|
|
|
|
|
|
2019-12-04 14:47:55 +05:30
|
|
|
def regrid(size,F,new_grid):
|
2020-01-23 21:45:02 +05:30
|
|
|
"""
|
2020-01-29 21:01:05 +05:30
|
|
|
Return mapping from coordinates in deformed configuration to a regular grid.
|
2020-01-23 21:45:02 +05:30
|
|
|
|
|
|
|
Parameters
|
|
|
|
----------
|
|
|
|
size : numpy.ndarray
|
|
|
|
physical size
|
|
|
|
F : numpy.ndarray
|
|
|
|
deformation gradient field
|
|
|
|
new_grid : numpy.ndarray
|
|
|
|
new grid for undeformed coordinates
|
|
|
|
|
|
|
|
"""
|
2019-12-05 23:02:21 +05:30
|
|
|
c = cell_coord0(F.shape[:3][::-1],size) \
|
|
|
|
+ cell_displacement_avg(size,F) \
|
|
|
|
+ cell_displacement_fluct(size,F)
|
2019-12-04 14:47:55 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
outer = _np.dot(_np.average(F,axis=(0,1,2)),size)
|
2019-12-04 14:47:55 +05:30
|
|
|
for d in range(3):
|
2020-04-10 16:00:39 +05:30
|
|
|
c[_np.where(c[:,:,:,d]<0)] += outer[d]
|
|
|
|
c[_np.where(c[:,:,:,d]>outer[d])] -= outer[d]
|
2019-12-04 14:47:55 +05:30
|
|
|
|
2020-04-10 16:00:39 +05:30
|
|
|
tree = _spatial.cKDTree(c.reshape(-1,3),boxsize=outer)
|
2020-01-23 21:45:02 +05:30
|
|
|
return tree.query(cell_coord0(new_grid,outer))[1].flatten()
|