Continuum Mechanics package

Vector calculus module

vector.antisym_grad(A, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Antisymmetric part of the gradient of a vector function A.

A : Matrix (3, 1), list
Vector function to compute the gradient from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
antisym_grad: Matrix (3, 3)
Matrix with the components of the antisymmetric part of the gradient. The position (i, j) has as components diff(A[i], coords[j].
vector.biharmonic(u, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Bilaplacian of the scalar function u.

u : SymPy expression
Scalar function to compute the bilaplacian from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameters, and it takes a cartesian (x, y, z), as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
bilaplacian: Sympy expression
Bilaplacian of u.
vector.curl(A, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Curl of a function vector A.

A : Matrix, List
Vector function to compute the curl from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
curl : Matrix (3, 1)
Column vector with the curl of A.
vector.div(A, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Divergence of the vector function A.

A : Matrix, list
Vector function to compute the divergence from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
divergence: SymPy expression
Divergence of A.
vector.div_tensor(tensor, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Divergence of a (second order) tensor

tensor : Matrix (3, 3)
Tensor function function to compute the divergence from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
divergence: Matrix
Divergence of tensor.
[RICHARDS]Rowland Richards. Principles of Solids Mechanics. CRC Press, 2011.
vector.dual_tensor(vec)[source]

Compute the dual tensor for an axial vector

In index notation, the dual is defined by

\[C_{ij} = \epsilon_{ijk} C_k\]

where \(\epsilon_{ijk}\) is the Levi-Civita symbol.

vec : Matrix (3)
Axial vector.
dual: Matrix (3, 3)
Second order matrix that is dual of vec.
[ARFKEN]George Arfken, Hans J. Weber and Frank Harris. Mathematical methods for physicists, Elsevier, 2013.
vector.dual_vector(tensor)[source]

Compute the dual (axial) vector for an anti-symmetric tensor

In index notation, the dual is defined by

\[C_{i} = \frac{1}{2}\epsilon_{ijk} C_{jk}\]

where \(\epsilon_{ijk}\) is the Levi-Civita symbol.

tensor : Matrix (3, 3)
Second order tensor.
dual: Matrix (3)
Axial vector that is the dual of tensor.
[ARFKEN]George Arfken, Hans J. Weber and Frank Harris. Mathematical methods for physicists, Elsevier, 2013.
vector.grad(u, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Compute the gradient of a scalara function phi.

u : SymPy expression
Scalar function to compute the gradient from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameters, and it takes a cartesian (x, y, z), as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
gradient: Matrix (3, 1)
Column vector with the components of the gradient.
vector.grad_vec(A, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Gradient of a vector function A.

A : Matrix (3, 1), list
Vector function to compute the gradient from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
gradient: Matrix (3, 3)
Matrix with the components of the gradient. The position (i, j) has as components diff(A[i], coords[j].
vector.lap(u, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Laplacian of the scalar function u.

u : SymPy expression
Scalar function to compute the laplacian from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameters, and it takes a cartesian (x, y, z), as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
laplacian: Sympy expression
Laplacian of u.
vector.lap_vec(A, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Laplacian of a vector function A.

A : Matrix, List
Vector function to compute the laplacian from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
laplacian : Matrix (3, 1)
Column vector with the components of the Laplacian.
vector.levi_civita(i, j, k)[source]

Levi-Civita symbol

vector.scale_coeff(r_vec, coords)[source]

Compute scale coefficients for the vector tranform given by r_vec.

r_vec : Matrix (3, 1)
Transform vector (x, y, z) as a function of coordinates u1, u2, u3.
coords : Tuple (3)
Coordinates for the new reference system.
h_vec : Tuple (3)
Scale coefficients.
vector.scale_coeff_coords(coord_sys, coords, a=1, b=1, c=1)[source]

Return scale factors for predefined coordinate system.

coord_sys : string
Coordinate system.
coords : Tuple (3)
Coordinates for the new reference system.
a : SymPy expression, optional
Additional parameter for some coordinate systems.
b : SymPy expression, optional
Additional parameter for some coordinate systems.
c : SymPy expression, optional
Additional parameter for some coordinate systems.
h_vec : Tuple (3)
Scale coefficients.
[ORTHO]Wikipedia contributors, ‘Orthogonal coordinates’, Wikipedia, The Free Encyclopedia, 2019
vector.sym_grad(A, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Symmetric part of the gradient of a vector function A.

A : Matrix (3, 1), list
Vector function to compute the gradient from.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
sym_grad: Matrix (3, 3)
Matrix with the components of the symmetric part of the gradient. The position (i, j) has as components diff(A[i], coords[j].
vector.transform_coords(coord_sys, coords, a=1, b=1, c=1)[source]

Return transformation for predefined coordinate systems.

coord_sys : string
Coordinate system.
coords : Tuple (3)
Coordinates for the new reference system.
a : SymPy expression, optional
Additional parameter for some coordinate systems.
b : SymPy expression, optional
Additional parameter for some coordinate systems.
c : SymPy expression, optional
Additional parameter for some coordinate systems.
h_vec : Tuple (3)
Scale coefficients.
[ORTHO]Wikipedia contributors, ‘Orthogonal coordinates’, Wikipedia, The Free Encyclopedia, 2019
vector.unit_vec_deriv(vec_i, coord_j, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Compute the derivatives of unit vectors with respect to coordinates

The derivative is defined as

\[\begin{split}\frac{\partial{\hat{\mathbf{e}}_i}}{\partial u_j} = \begin{cases} \hat{\mathbf{e}_j} \frac{1}{h_i} \frac{\partial{h_j}}{\partial u_i} &\text{if } i\neq j\\ -\sum_{\substack{k=1\\ k\neq i}}^3 \hat{\mathbf{e}_k} \frac{1}{h_k} \frac{\partial{h_i}}{\partial u_k} &\text{if } i = j\\ \end{cases}\, ,\end{split}\]

as presented in [ARFKEN].

vec_i : int
Number of the unit vector (0, 1, 2).
coord_j : int
Number of the coordinate (0, 1, 2).
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
deriv: Matrix (3, 1)
Derivative of the i-th unit vector with respect to the j-th coordinate.
[ARFKEN]George Arfken, Hans J. Weber and Frank Harris. Mathematical methods for physicists, Elsevier, 2013.

Solid mechanics module

solids.c_cst(u, parameters, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Corrected-Couple-Stress-Theory (C-CST) elasticity operator of a vector function u, as presented in [CST].

u : Matrix (3, 1), list
Vector function to apply the Navier-Cauchy operator.
parameters : tuple

Material parameters in the following order:

lamda : float
Lamé’s first parameter.
mu : float, > 0
Lamé’s second parameter.
eta : float, > 0
Couple stress modulus in C-CST.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
c_cst : Matrix (3, 1)
Components of the C-CST operator applied to the displacement vector
[CST](1, 2, 3) Ali R. Hadhesfandiari, Gary F. Dargush. Couple stress theory for solids. International Journal for Solids and Structures, 2011, 48, 2496-2510.
solids.disp_def_cst(u, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Compute strain measures for C-CST elasticity, as defined in [CST].

u : Matrix (3, 1), list
Displacement vector function to apply the micropolar operator.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
strain : Matrix (3, 3)
Strain tensor.
curvature : Matrix (3, 3)
Curvature tensor.
solids.disp_def_micropolar(u, phi, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Compute strain measures for micropolar elasticity, as defined in [NOW].

u : Matrix (3, 1), list
Displacement vector function to apply the micropolar operator.
phi : Matrix (3, 1), list
Microrrotation (pseudo)-vector function to apply the micropolar operator.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
strain : Matrix (3, 3)
Strain tensor.
curvature : Matrix (3, 3)
Curvature tensor.
solids.micropolar(u, phi, parameters, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Micropolar operator of a vector function u, as defined in [NOW].

u : Matrix (3, 1), list
Displacement vector function to apply the micropolar operator.
phi : Matrix (3, 1), list
Microrrotation (pseudo)-vector function to apply the micropolar operator.
parameters : tuple

Material parameters in the following order:

lamda : float
Lamé’s first parameter.
mu : float, > 0
Lamé’s second parameter.
alpha : float, > 0
Micropolar parameter.
beta : float
Micropolar parameter.
gamma : float, > 0
Micropolar parameter.
epsilon : float, > 0
Micropolar parameter.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
disp_op : Matrix (3, 1)
Displacement components.
rot_op : Matrix (3, 1)
Microrrotational components.
[NOW](1, 2) Witold Nowacki. Theory of micropolar elasticity. International centre for mechanical sciences, Courses and lectures, No. 25. Berlin: Springer, 1972.
solids.navier_cauchy(u, parameters, coords=(x, y, z), h_vec=(1, 1, 1))[source]

Navier-Cauchy operator of a vector function u.

u : Matrix (3, 1), list
Vector function to apply the Navier-Cauchy operator.
parameters : tuple

Material parameters in the following order:

lamda : float
Lamé’s first parameter.
mu : float, > 0
Lamé’s second parameter.
coords : Tuple (3), optional
Coordinates for the new reference system. This is an optional parameter it takes (x, y, z) as default.
h_vec : Tuple (3), optional
Scale coefficients for the new coordinate system. It takes (1, 1, 1), as default.
navier_op : Matrix (3, 1)
Components of the Navier-Cauchy operator applied to the displacement vector.
solids.strain_stress(strain, parameters)[source]

Return the stress tensor for a given strain tensor and material properties lambda and mu.

u : Matrix (3, 3)
Strain tensor.
parameters : tuple

Material parameters in the following order:

lamda : float
Lamé’s first parameter.
mu : float, > 0
Lamé’s second parameter.
disp_op : Matrix (3, 1)
Displacement components.
solids.strain_stress_cst(strain, curvature, parameters)[source]

Return the force-stress and couple-stress tensor for given strain and curvature tensors and material parameters for the Corrected-Couple-Stress-Theory (C-CST), as presented in [CST].

strain : Matrix (3, 3)
Strain tensor.
curvature : Matrix (3, 3)
Curvature tensor.
parameters : tuple

Material parameters in the following order:

lamda : float
Lamé’s first parameter.
mu : float, > 0
Lamé’s second parameter.
eta : float, > 0
Couple stress modulus in C-CST.
solids.strain_stress_micropolar(strain, curvature, parameters)[source]

Return the force-stress and couple-stress tensor for given strain and curvature tensors and material parameters.

strain : Matrix (3, 3)
Strain tensor.
curvature : Matrix (3, 3)
Curvature tensor.
parameters : tuple

Material parameters in the following order:

lamda : float
Lamé’s first parameter.
mu : float, > 0
Lamé’s second parameter.
alpha : float, > 0
Micropolar parameter.
beta : float
Micropolar parameter.
gamma : float, > 0
Micropolar parameter.
epsilon : float, > 0
Micropolar parameter.

Visualization module

Functions to aid the visualization of mathematical entities such as second rank tensors.

visualization.christofel_eig(C, nphi=30, nth=30)[source]

Compute surfaces of eigenvalues for the Christoffel stiffness

For every direction

\[\mathbf{n} =(\sin(\theta)\cos(\phi),\sin(\theta)\sin(\phi),\cos(\theta))\]

computes the eigenvalues of the Christoffel stiffness tensor. These eigencalues are \(\rho v_p^2\), where \(\rho\) is the mass density and \(v_p\) is the phase speed.

C : (6,6) array
Stiffness tensor in Voigt notation.
nphi : (optional) int
Number of partitions in the azimut angle (phi).
nth : (optional) int
Number of partitions in the cenit angle (theta).
V1, V2, V3 : (nphi, nth) arrays
Eigenvalues for the desired discretization.
phi_vec : (nphi) array
Array with azimut angles.
theta_vec : (nth) array
Array with cenit angles.
visualization.mohr2d(stress, ax=None)[source]

Plot Mohr circle for a 2D tensor

stress : ndarray
Stress tensor.
ax : Matplotlib axes, optional
Axes where the plot is going to be added.
[BRAN]Brannon, R. (2003). Mohr’s Circle and more circles. Poslední revize, 29(10).
visualization.mohr3d(stress, ax=None)[source]

Plot 3D Mohr circles

stress : ndarray
Stress tensor.
ax : Matplotlib axes, optional
Axes where the plot is going to be added.
visualization.plot_surf(R, phi, theta, title='Wave speed', **kwargs)[source]

Plot the surface given by R(phi, theta).

R : (m,n) ndarray
Radius function.
phi_vec : (m,n) ndarray
Meshgrid for the azimut angle (phi).
theta_vec : (m,n) ndarray
Meshgrid for the cenit angle (theta).
**kwargs : keyword arguments (optional)
Keyword arguments for mlab.mesh.
surf : mayavi mesh
Mayavi mesh for the surface R(phi, theta).
visualization.traction_circle(stress, npts=48, ax=None)[source]

Visualize a second order tensor as a collection of tractions vectors over a circle.

stress : ndarray
Stress tensor.
npts : int, optional
Number of vector to plot over the circle.
ax : Matplotlib axes, optional
Axes where the plot is going to be added.