amfe.element module

Element Module in which the Finite Elements are described on Element level.

This Module is arbitrarily extensible. The idea is to use the basis class Element which provides the functionality for an efficient solution of a time integration by only once calling the internal tensor computation and then extracting the tangential stiffness matrix and the internal force vector in one run.

Some remarks resulting in the observations of the profiler: Most of the time is spent with python-functions, when they are used. For instance the kron-function in order to build the scattered geometric stiffness matrix or the trace function are very inefficient. They can be done better when using direct functions.

If some element things are really time critical, it is recommended to port the heavy computation to FORTRAN. This can be achieved by using the provided f2py routines and reprogram the stuff for the own use.

class amfe.element.Element(material=None)[source]

Bases: object

Anonymous baseclass for all elements. It contains the methods needed for the computation of the element stuff.

material

instance of amfe.HyperelasticMaterial – Class containing the material behavior.

name

str – Name for the postprocessing tool to identify the characteristics of the element

f_int(X, u, t=0)[source]

Returns the internal element restoring force f_int

Parameters:
  • X (ndarray) – nodal coordinates given in Voigt notation (i.e. a 1-D-Array of type [x_1, y_1, z_1, x_2, y_2, z_2 etc.])
  • u (ndarray) – nodal displacements given in Voigt notation
  • t (float, optional) – time, default value: 0.
Returns:

f_int – The nodal force vector (numpy.ndarray of dimension (ndim,))

Return type:

ndarray

k_and_f_int(X, u, t=0)[source]

Returns the tangential stiffness matrix and the internal nodal force of the Element.

Parameters:
  • X (ndarray) – nodal coordinates given in Voigt notation (i.e. a 1-D-Array of type [x_1, y_1, z_1, x_2, y_2, z_2 etc.])
  • u (ndarray) – nodal displacements given in Voigt notation
  • t (float) – time
Returns:

  • k_int (ndarray) – The tangential stiffness matrix (ndarray of dimension (ndim, ndim))
  • f_int (ndarray) – The nodal force vector (ndarray of dimension (ndim,))

Examples

TODO

k_f_S_E_int(X, u, t=0)[source]

Returns the tangential stiffness matrix, the internal nodal force, the strain and the stress tensor (voigt-notation) of the Element.

Parameters:
  • X (ndarray) – nodal coordinates given in Voigt notation (i.e. a 1-D-Array of type [x_1, y_1, z_1, x_2, y_2, z_2 etc.])
  • u (ndarray) – nodal displacements given in Voigt notation
  • t (float) – time
Returns:

  • K (ndarray) – The tangential stiffness matrix (ndarray of dimension (ndim, ndim))
  • f (ndarray) – The nodal force vector (ndarray of dimension (ndim,))
  • S (ndarray) – The stress tensor (ndarray of dimension (no_of_nodes, 6))
  • E (ndarray) – The strain tensor (ndarray of dimension (no_of_nodes, 6))

Examples

TODO

k_int(X, u, t=0)[source]

Returns the tangential stiffness matrix of the Element.

Parameters:
  • X (ndarray) – nodal coordinates given in Voigt notation (i.e. a 1-D-Array of type [x_1, y_1, z_1, x_2, y_2, z_2 etc.])
  • u (ndarray) – nodal displacements given in Voigt notation
  • t (float) – time
Returns:

k_int – The tangential stiffness matrix (numpy.ndarray of type ndim x ndim)

Return type:

ndarray

m_and_vec_int(X, u, t=0)[source]

Returns mass matrix of the Element and zero vector of size X.

Parameters:
  • X (ndarray) – nodal coordinates given in Voigt notation (i.e. a 1-D-Array of type [x_1, y_1, z_1, x_2, y_2, z_2 etc.])
  • u (ndarray) – nodal displacements given in Voigt notation
  • t (float, optional) – time, default value: 0.
Returns:

  • m_int (ndarray) – The consistent mass matrix of the element (numpy.ndarray of dimension (ndim,ndim))
  • vec (ndarray) – vector (containing zeros) of dimension (ndim,)

m_int(X, u, t=0)[source]

Returns the mass matrix of the element.

Parameters:
  • X (ndarray) – nodal coordinates given in Voigt notation (i.e. a 1-D-Array of type [x_1, y_1, z_1, x_2, y_2, z_2 etc.])
  • u (ndarray) – nodal displacements given in Voigt notation
  • t (float, optional) – time, default value: 0.
Returns:

m_int – The consistent mass matrix of the element (numpy.ndarray of dimension (ndim,ndim))

Return type:

ndarray

name = None
class amfe.element.Tri3(*args, **kwargs)[source]

Bases: amfe.element.Element

Element class for a plane triangle element in Total Lagrangian formulation. The displacements are given in x- and y-coordinates;

Notes

The Element assumes constant strain and stress over the whole element. Thus the approximation quality is very moderate.

References

Basis for this implementation is the Monograph of Ted Belytschko: Nonlinear Finite Elements for Continua and Structures. pp. 201 and 207.

name = 'Tri3'
plane_stress = True
class amfe.element.Tri6(*args, **kwargs)[source]

Bases: amfe.element.Element

6 node second order triangle Triangle Element with 6 dofs; 3 dofs at the corner, 3 dofs in the intermediate point of every face.

name = 'Tri6'
plane_stress = True
class amfe.element.Quad4(*args, **kwargs)[source]

Bases: amfe.element.Element

Quadrilateral 2D element with bilinear shape functions.

        ^ eta
4       |       3
  o_____|_____o
  |     |     |
  |     |     |
--------+---------->
  |     |     |   xi
  |     |     |
  o_____|_____o
1       |       2
name = 'Quad4'
class amfe.element.Quad8(*args, **kwargs)[source]

Bases: amfe.element.Element

Plane Quadrangle with quadratic shape functions and 8 nodes. 4 nodes are at every corner, 4 nodes on every face.

name = 'Quad8'
class amfe.element.Tet4(*args, **kwargs)[source]

Bases: amfe.element.Element

Tetraeder-Element with 4 nodes

name = 'Tet4'
class amfe.element.Tet10(*args, **kwargs)[source]

Bases: amfe.element.Element

Tet10 solid element; Node numbering is done like in ParaView and in [R1]

The node numbering is as follows:

References

[R1]Felippa, Carlos: Advanced Finite Element Methods (ASEN 6367), Spring 2013. Online Source
name = 'Tet10'
class amfe.element.Hexa8(*args, **kwargs)[source]

Bases: amfe.element.Element

Eight point Hexahedron element.

         eta
3----------2
|\     ^   |        | \    |   |         |  \   |   |          |   7------+---6
|   |  +-- |-- | -> xi
0---+---\--1   |
 \  |    \  \  |
  \ |     \  \ |
   \|   zeta  \|
    4----------5
name = 'Hexa8'
class amfe.element.Hexa20(*args, **kwargs)[source]

Bases: amfe.element.Element

20-node brick element.

# # v #3———-2 3—-10—-2 3—-13—-2 #| ^ |\ |\ |\ |\ |#| | | | 19 | 18 |15 24 | 14 #| | | 11 9 9 20 11 #| 7——+—6 | 7—-14+—6 | 7—-19+—6 #| | +– |– | -> u | | | | |22 | 26 | 23| #0—+—–1 | 0—+-8—-1 | 0—+-8—-1 | # | | 15 13 17 25 18 # | | 16 | 17| 10 | 21 12| # | w | | | | | # 4———-5 4—-12—-5 4—-16—-5

name = 'Hexa20'
class amfe.element.Bar2Dlumped(*args, **kwargs)[source]

Bases: amfe.element.Element

Bar-Element with 2 nodes and lumped stiffness matrix

foo()[source]
class amfe.element.BoundaryElement(val, ndof, direct='normal', time_func=None, shadow_area=False)[source]

Bases: amfe.element.Element

Class for the application of Neumann Boundary Conditions.

time_func

func – function returning a value between {-1, 1} which is time dependent storing the time dependency of the Neumann Boundary condition. Example for constant function:

>>> def func(t):
>>>    return 1
f_proj

func – function producing the nodal force vector from the given nodal force vector in normal direction.

val

float

direct

{‘normal’, numpy.array} – ‘normal’ means that force acts normal to the current surface numpy.array is a vector in which the force should act (with global coordinate system as reference)

f

numpy.array – local external force vector of the element

K

numpy.array – tangent stiffness matrix of the element (for boundary elements typically zero)

M

numpy.array – mass matrix of the element (for boundary elements typically zero)

class amfe.element.Tri3Boundary(val, direct, time_func=None, shadow_area=False)[source]

Bases: amfe.element.BoundaryElement

Class for application of Neumann Boundary Conditions.

class amfe.element.Tri6Boundary(val, direct, time_func=None, shadow_area=False)[source]

Bases: amfe.element.BoundaryElement

Boundary element with variatonally consistent boundary forces.

Notes

This function has been updated to give a variationally consistent integrated skin element.

gauss_points = ((0.16666666666666666, 0.16666666666666666, 0.6666666666666666, 0.3333333333333333), (0.16666666666666666, 0.6666666666666666, 0.16666666666666666, 0.3333333333333333), (0.6666666666666666, 0.16666666666666666, 0.16666666666666666, 0.3333333333333333))
class amfe.element.Quad4Boundary(val, direct, time_func=None, shadow_area=False)[source]

Bases: amfe.element.BoundaryElement

Quad4 boundary element for 3D-Problems.

g1 = 0.57735026918962584
gauss_points = ((-0.57735026918962584, -0.57735026918962584, 1.0), (0.57735026918962584, -0.57735026918962584, 1.0), (0.57735026918962584, 0.57735026918962584, 1.0), (-0.57735026918962584, 0.57735026918962584, 1.0))
class amfe.element.Quad8Boundary(val, direct, time_func=None, shadow_area=False)[source]

Bases: amfe.element.BoundaryElement

Quad8 boundary element for 3D-Problems.

g = 0.7745966692414834
gauss_points = ((-0.7745966692414834, -0.7745966692414834, 0.308641975308642), (0.7745966692414834, -0.7745966692414834, 0.308641975308642), (0.7745966692414834, 0.7745966692414834, 0.308641975308642), (-0.7745966692414834, 0.7745966692414834, 0.308641975308642), (0, -0.7745966692414834, 0.49382716049382713), (0.7745966692414834, 0, 0.49382716049382713), (0, 0.7745966692414834, 0.49382716049382713), (-0.7745966692414834, 0, 0.49382716049382713), (0, 0, 0.7901234567901234))
w = 0.5555555555555556
w0 = 0.8888888888888888
class amfe.element.LineLinearBoundary(val, direct, time_func=None, shadow_area=False)[source]

Bases: amfe.element.BoundaryElement

Line Boundary element for 2D-Problems.

N = array([ 0.5, 0.5])
rot_mat = array([[ 0, 1], [-1, 0]])
class amfe.element.LineQuadraticBoundary(val, direct, time_func=None, shadow_area=False)[source]

Bases: amfe.element.BoundaryElement

Quadratic line boundary element for 2D problems.

N = array([ 0.16666667, 0.16666667, 0.66666667])
rot_mat = array([[ 0, 1], [-1, 0]])