2018-11-17 12:42:12 +05:30
|
|
|
#!/usr/bin/env python3
|
2011-08-25 22:14:36 +05:30
|
|
|
|
2016-06-29 23:39:42 +05:30
|
|
|
import os
|
2019-12-21 19:42:01 +05:30
|
|
|
import sys
|
2019-12-21 22:49:54 +05:30
|
|
|
from io import StringIO
|
2019-06-14 16:33:30 +05:30
|
|
|
from optparse import OptionParser
|
|
|
|
|
2016-06-29 23:39:42 +05:30
|
|
|
import numpy as np
|
2019-06-14 16:33:30 +05:30
|
|
|
|
2014-07-10 14:57:51 +05:30
|
|
|
import damask
|
2011-08-25 22:14:36 +05:30
|
|
|
|
2019-06-14 16:33:30 +05:30
|
|
|
|
2016-01-27 22:36:00 +05:30
|
|
|
scriptName = os.path.splitext(os.path.basename(__file__))[0]
|
|
|
|
scriptID = ' '.join([scriptName,damask.version])
|
2011-08-25 22:14:36 +05:30
|
|
|
|
2016-03-24 18:49:00 +05:30
|
|
|
def volTetrahedron(coords):
|
2016-03-24 17:05:33 +05:30
|
|
|
"""
|
2016-04-15 04:02:30 +05:30
|
|
|
Return the volume of the tetrahedron with given vertices or sides.
|
2020-04-21 01:14:38 +05:30
|
|
|
|
|
|
|
If vertices are given they must be in a NumPy array with shape (4,3): the
|
2016-03-24 17:05:33 +05:30
|
|
|
position vectors of the 4 vertices in 3 dimensions; if the six sides are
|
|
|
|
given, they must be an array of length 6. If both are given, the sides
|
|
|
|
will be used in the calculation.
|
|
|
|
|
|
|
|
This method implements
|
|
|
|
Tartaglia's formula using the Cayley-Menger determinant:
|
|
|
|
|0 1 1 1 1 |
|
|
|
|
|1 0 s1^2 s2^2 s3^2|
|
|
|
|
288 V^2 = |1 s1^2 0 s4^2 s5^2|
|
|
|
|
|1 s2^2 s4^2 0 s6^2|
|
|
|
|
|1 s3^2 s5^2 s6^2 0 |
|
|
|
|
where s1, s2, ..., s6 are the tetrahedron side lengths.
|
|
|
|
|
|
|
|
from http://codereview.stackexchange.com/questions/77593/calculating-the-volume-of-a-tetrahedron
|
|
|
|
"""
|
|
|
|
# The indexes of rows in the vertices array corresponding to all
|
|
|
|
# possible pairs of vertices
|
|
|
|
vertex_pair_indexes = np.array(((0, 1), (0, 2), (0, 3),
|
|
|
|
(1, 2), (1, 3), (2, 3)))
|
|
|
|
|
|
|
|
# Get all the squares of all side lengths from the differences between
|
|
|
|
# the 6 different pairs of vertex positions
|
2020-03-17 16:52:48 +05:30
|
|
|
vertices = np.concatenate((coords[0],coords[1],coords[2],coords[3])).reshape(4,3)
|
2016-03-24 17:05:33 +05:30
|
|
|
vertex1, vertex2 = vertex_pair_indexes[:,0], vertex_pair_indexes[:,1]
|
|
|
|
sides_squared = np.sum((vertices[vertex1] - vertices[vertex2])**2,axis=-1)
|
|
|
|
|
|
|
|
|
|
|
|
# Set up the Cayley-Menger determinant
|
|
|
|
M = np.zeros((5,5))
|
|
|
|
# Fill in the upper triangle of the matrix
|
|
|
|
M[0,1:] = 1
|
|
|
|
# The squared-side length elements can be indexed using the vertex
|
|
|
|
# pair indices (compare with the determinant illustrated above)
|
|
|
|
M[tuple(zip(*(vertex_pair_indexes + 1)))] = sides_squared
|
|
|
|
|
|
|
|
# The matrix is symmetric, so we can fill in the lower triangle by
|
|
|
|
# adding the transpose
|
|
|
|
M = M + M.T
|
2016-03-24 18:49:00 +05:30
|
|
|
return np.sqrt(np.linalg.det(M) / 288)
|
2016-03-24 17:05:33 +05:30
|
|
|
|
|
|
|
|
2016-03-24 18:49:00 +05:30
|
|
|
def volumeMismatch(size,F,nodes):
|
2016-03-24 17:05:33 +05:30
|
|
|
"""
|
2019-12-21 22:49:54 +05:30
|
|
|
Calculates the volume mismatch.
|
2020-04-21 01:14:38 +05:30
|
|
|
|
|
|
|
volume mismatch is defined as the difference between volume of reconstructed
|
2019-12-21 22:49:54 +05:30
|
|
|
(compatible) cube and determinant of deformation gradient at Fourier point.
|
2016-03-24 17:05:33 +05:30
|
|
|
"""
|
2016-03-24 18:49:00 +05:30
|
|
|
coords = np.empty([8,3])
|
2020-04-21 01:14:38 +05:30
|
|
|
vMismatch = np.empty(F.shape[:3])
|
|
|
|
|
2016-03-24 17:05:33 +05:30
|
|
|
#--------------------------------------------------------------------------------------------------
|
|
|
|
# calculate actual volume and volume resulting from deformation gradient
|
2020-04-21 01:14:38 +05:30
|
|
|
for k in range(grid[0]):
|
2016-10-25 00:46:29 +05:30
|
|
|
for j in range(grid[1]):
|
2020-04-21 01:14:38 +05:30
|
|
|
for i in range(grid[2]):
|
2016-06-29 23:39:42 +05:30
|
|
|
coords[0,0:3] = nodes[k, j, i ,0:3]
|
|
|
|
coords[1,0:3] = nodes[k ,j, i+1,0:3]
|
|
|
|
coords[2,0:3] = nodes[k ,j+1,i+1,0:3]
|
|
|
|
coords[3,0:3] = nodes[k, j+1,i ,0:3]
|
|
|
|
coords[4,0:3] = nodes[k+1,j, i ,0:3]
|
|
|
|
coords[5,0:3] = nodes[k+1,j, i+1,0:3]
|
|
|
|
coords[6,0:3] = nodes[k+1,j+1,i+1,0:3]
|
|
|
|
coords[7,0:3] = nodes[k+1,j+1,i ,0:3]
|
|
|
|
vMismatch[k,j,i] = \
|
|
|
|
( abs(volTetrahedron([coords[6,0:3],coords[0,0:3],coords[7,0:3],coords[3,0:3]])) \
|
2016-03-24 18:49:00 +05:30
|
|
|
+ abs(volTetrahedron([coords[6,0:3],coords[0,0:3],coords[7,0:3],coords[4,0:3]])) \
|
|
|
|
+ abs(volTetrahedron([coords[6,0:3],coords[0,0:3],coords[2,0:3],coords[3,0:3]])) \
|
|
|
|
+ abs(volTetrahedron([coords[6,0:3],coords[0,0:3],coords[2,0:3],coords[1,0:3]])) \
|
|
|
|
+ abs(volTetrahedron([coords[6,0:3],coords[4,0:3],coords[1,0:3],coords[5,0:3]])) \
|
2016-06-29 23:39:42 +05:30
|
|
|
+ abs(volTetrahedron([coords[6,0:3],coords[4,0:3],coords[1,0:3],coords[0,0:3]]))) \
|
|
|
|
/np.linalg.det(F[k,j,i,0:3,0:3])
|
2020-04-21 01:14:38 +05:30
|
|
|
return vMismatch/(size.prod()/grid.prod())
|
2016-03-24 18:49:00 +05:30
|
|
|
|
|
|
|
|
|
|
|
def shapeMismatch(size,F,nodes,centres):
|
2016-03-24 17:05:33 +05:30
|
|
|
"""
|
2019-12-21 22:49:54 +05:30
|
|
|
Routine to calculate the shape mismatch.
|
2020-04-21 01:14:38 +05:30
|
|
|
|
2016-04-15 04:02:30 +05:30
|
|
|
shape mismatch is defined as difference between the vectors from the central point to
|
2016-03-24 17:05:33 +05:30
|
|
|
the corners of reconstructed (combatible) volume element and the vectors calculated by deforming
|
2019-12-21 22:49:54 +05:30
|
|
|
the initial volume element with the current deformation gradient.
|
2016-03-24 17:05:33 +05:30
|
|
|
"""
|
2020-04-21 01:14:38 +05:30
|
|
|
sMismatch = np.empty(F.shape[:3])
|
|
|
|
|
2016-03-24 18:49:00 +05:30
|
|
|
#--------------------------------------------------------------------------------------------------
|
|
|
|
# initial positions
|
2020-04-21 01:44:57 +05:30
|
|
|
delta = size/grid*.5
|
|
|
|
coordsInitial = np.vstack((delta * np.array((-1,-1,-1)),
|
|
|
|
delta * np.array((+1,-1,-1)),
|
|
|
|
delta * np.array((+1,+1,-1)),
|
|
|
|
delta * np.array((-1,+1,-1)),
|
|
|
|
delta * np.array((-1,-1,+1)),
|
|
|
|
delta * np.array((+1,-1,+1)),
|
|
|
|
delta * np.array((+1,+1,+1)),
|
|
|
|
delta * np.array((-1,+1,+1))))
|
2020-04-21 01:14:38 +05:30
|
|
|
|
2016-03-24 18:49:00 +05:30
|
|
|
#--------------------------------------------------------------------------------------------------
|
|
|
|
# compare deformed original and deformed positions to actual positions
|
2020-04-21 01:14:38 +05:30
|
|
|
for k in range(grid[0]):
|
2016-10-25 00:46:29 +05:30
|
|
|
for j in range(grid[1]):
|
2020-04-21 01:14:38 +05:30
|
|
|
for i in range(grid[2]):
|
2016-06-29 23:39:42 +05:30
|
|
|
sMismatch[k,j,i] = \
|
|
|
|
+ np.linalg.norm(nodes[k, j, i ,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[0,0:3]))\
|
2020-04-21 01:14:38 +05:30
|
|
|
+ np.linalg.norm(nodes[k+1,j, i ,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[1,0:3]))\
|
|
|
|
+ np.linalg.norm(nodes[k+1,j+1,i ,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[2,0:3]))\
|
2016-06-29 23:39:42 +05:30
|
|
|
+ np.linalg.norm(nodes[k, j+1,i ,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[3,0:3]))\
|
2020-04-21 01:14:38 +05:30
|
|
|
+ np.linalg.norm(nodes[k, j, i+1,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[4,0:3]))\
|
2016-06-29 23:39:42 +05:30
|
|
|
+ np.linalg.norm(nodes[k+1,j, i+1,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[5,0:3]))\
|
|
|
|
+ np.linalg.norm(nodes[k+1,j+1,i+1,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[6,0:3]))\
|
2020-04-21 01:14:38 +05:30
|
|
|
+ np.linalg.norm(nodes[k ,j+1,i+1,0:3] - centres[k,j,i,0:3] - np.dot(F[k,j,i,:,:], coordsInitial[7,0:3]))
|
2016-03-24 18:49:00 +05:30
|
|
|
return sMismatch
|
|
|
|
|
2016-03-24 17:05:33 +05:30
|
|
|
|
2011-08-25 22:14:36 +05:30
|
|
|
# --------------------------------------------------------------------
|
|
|
|
# MAIN
|
|
|
|
# --------------------------------------------------------------------
|
|
|
|
|
2019-02-16 22:11:56 +05:30
|
|
|
parser = OptionParser(option_class=damask.extendableOption, usage='%prog options [ASCIItable(s)]', description = """
|
2015-08-08 00:33:26 +05:30
|
|
|
Add column(s) containing the shape and volume mismatch resulting from given deformation gradient.
|
|
|
|
Operates on periodic three-dimensional x,y,z-ordered data sets.
|
2011-08-25 22:14:36 +05:30
|
|
|
|
2014-08-06 18:57:09 +05:30
|
|
|
""", version = scriptID)
|
|
|
|
|
2011-08-25 22:14:36 +05:30
|
|
|
|
2015-08-08 00:33:26 +05:30
|
|
|
parser.add_option('-c','--coordinates',
|
2016-04-27 02:19:58 +05:30
|
|
|
dest = 'pos',
|
2015-08-08 00:33:26 +05:30
|
|
|
type = 'string', metavar = 'string',
|
|
|
|
help = 'column heading of coordinates [%default]')
|
|
|
|
parser.add_option('-f','--defgrad',
|
|
|
|
dest = 'defgrad',
|
|
|
|
type = 'string', metavar = 'string ',
|
|
|
|
help = 'column heading of deformation gradient [%default]')
|
|
|
|
parser.add_option('--no-shape','-s',
|
|
|
|
dest = 'shape',
|
|
|
|
action = 'store_false',
|
|
|
|
help = 'omit shape mismatch')
|
|
|
|
parser.add_option('--no-volume','-v',
|
|
|
|
dest = 'volume',
|
|
|
|
action = 'store_false',
|
|
|
|
help = 'omit volume mismatch')
|
2016-04-27 02:19:58 +05:30
|
|
|
parser.set_defaults(pos = 'pos',
|
|
|
|
defgrad = 'f',
|
|
|
|
shape = True,
|
|
|
|
volume = True,
|
2015-08-08 00:33:26 +05:30
|
|
|
)
|
2011-08-25 22:14:36 +05:30
|
|
|
|
|
|
|
(options,filenames) = parser.parse_args()
|
2015-08-21 01:12:05 +05:30
|
|
|
if filenames == []: filenames = [None]
|
2015-08-08 00:33:26 +05:30
|
|
|
|
2019-12-21 19:42:01 +05:30
|
|
|
|
2014-07-22 01:25:05 +05:30
|
|
|
for name in filenames:
|
2015-09-24 14:54:42 +05:30
|
|
|
damask.util.report(scriptName,name)
|
2020-04-21 01:14:38 +05:30
|
|
|
|
2019-12-21 19:42:01 +05:30
|
|
|
table = damask.Table.from_ASCII(StringIO(''.join(sys.stdin.read())) if name is None else name)
|
2020-01-13 07:21:49 +05:30
|
|
|
grid,size,origin = damask.grid_filters.cell_coord0_gridSizeOrigin(table.get(options.pos))
|
2016-06-27 23:08:12 +05:30
|
|
|
|
2020-04-21 01:14:38 +05:30
|
|
|
F = table.get(options.defgrad).reshape(tuple(grid)+(-1,),order='F').reshape(tuple(grid)+(3,3))
|
2019-12-21 22:47:04 +05:30
|
|
|
nodes = damask.grid_filters.node_coord(size,F)
|
2020-04-21 01:14:38 +05:30
|
|
|
|
2016-06-27 23:08:12 +05:30
|
|
|
if options.shape:
|
2020-01-12 01:03:29 +05:30
|
|
|
centers = damask.grid_filters.cell_coord(size,F)
|
2020-04-21 01:14:38 +05:30
|
|
|
shapeMismatch = shapeMismatch(size,F,nodes,centers)
|
2019-12-21 19:42:01 +05:30
|
|
|
table.add('shapeMismatch(({}))'.format(options.defgrad),
|
2020-04-21 01:14:38 +05:30
|
|
|
shapeMismatch.reshape(-1,1,order='F'),
|
2019-12-21 19:42:01 +05:30
|
|
|
scriptID+' '+' '.join(sys.argv[1:]))
|
2020-04-21 01:14:38 +05:30
|
|
|
|
2016-06-27 23:08:12 +05:30
|
|
|
if options.volume:
|
2020-04-21 01:14:38 +05:30
|
|
|
volumeMismatch = volumeMismatch(size,F,nodes)
|
2019-12-21 19:42:01 +05:30
|
|
|
table.add('volMismatch(({}))'.format(options.defgrad),
|
2020-04-21 01:14:38 +05:30
|
|
|
volumeMismatch.reshape(-1,1,order='F'),
|
2019-12-21 19:42:01 +05:30
|
|
|
scriptID+' '+' '.join(sys.argv[1:]))
|
2016-06-27 23:08:12 +05:30
|
|
|
|
2020-09-03 20:49:19 +05:30
|
|
|
table.to_file(sys.stdout if name is None else name)
|