Triangle class

Coordinate-free formulation for a finite deformation bi-linear triangle.

Linear Theory (small displacements)

The element uses an iso-parametric mapping for both the reference and the spatial configuration. In the reference configuration,

(1)\[\mathbf{X}=\mathbf{X}_I \xi^I \quad I=0,1,2\equiv u,s,t \quad \xi^I \ge 0 ~\&~ \xi^0+\xi^1+\xi^2 = 1\]

where \(\xi^I\) are standard triangle coordinates. Alternatively, the mapping may be expressed using just two coordinates, \(\xi^1\equiv\xi^s\) and \(\xi^2\equiv\xi^t\), as

(2)\[\mathbf{X}=\mathbf{X}_0+\mathbf{G}_i \xi^i \quad i=1,2 \quad \xi^i \ge 0 ~\&~ \xi^1+\xi^2 \leq 1\]

The covariant base vectors follow from (1) and (2) as

(3)\[\mathbf{G}_i = \mathbf{X}_{,i} = \mathbf{X}_i - \mathbf{X}_0 \quad i=1,2\]

We can express differential \(d\xi^i\) in terms of \(d\mathbf{X}\) by means of the dual base as

\[d\mathbf{X}=\mathbf{G}_i d \xi^i \rightarrow d \xi^i=\mathbf{G}^i \cdot d \mathbf{X} \quad w/ ~~ \mathbf{G}^i \mathbf{G}_j=\delta_j^i\]

The displacement field is represented as

\[\mathbf{u} = \mathbf{u}_I \xi^I \quad I=0,1,2\]

Using (3) we obtain the displacement gradient as

\[d\mathbf{u}=\frac{\partial\mathbf{u}}{\xi^i} d \xi^i=\left(\mathbf{u}_I \otimes \mathbf{G}^I\right) \cdot d \mathbf{X} =: \nabla\mathbf{u} \cdot d\mathbf{\bar{X}} \Longrightarrow \nabla\mathbf{u}=\mathbf{u}_I \otimes \mathbf{G}^I\]

where \(d\xi^0+d\xi^1+d\xi^2=0\) and \(\mathbf{G}^0=-\mathbf{G}^1-\mathbf{G}^2\) were used.

The strain tensors follows as

\[\symbf{\varepsilon} = \frac{1}{2} (\nabla\mathbf{u}^t + \nabla\mathbf{u}) = \frac{1}{2} (\mathbf{G}^I \otimes \mathbf{u}_I + \mathbf{u}_I \otimes \mathbf{G}^I)\]

The variation of strain is needed for the principle of virtual displacements. It follows as

\[\delta \symbf{\varepsilon} = \frac{1}{2} (\mathbf{G}^I \otimes \delta \mathbf{u}_I + \delta \mathbf{u}_I \otimes \mathbf{G}^I) \qquad I=0,1,2\]

Application of virtual work

The principle of virtual displacements (PVD) takes the form

\[\int\!\!\!\int\limits_{\mathcal{A}} \delta \symbf{\varepsilon}:\symbf{\sigma}\: d\mathcal{A} = \int\!\!\!\int\limits_{\mathcal{A}} \delta\mathbf{u}\cdot\mathbf{b} \: d\mathcal{A} + \int\limits_{\partial\mathcal{A}_\sigma} \delta\mathbf{u}\cdot\bar{\mathbf{T}} \: d\mathcal{s}\]

where \(\symbf{\sigma}=\mathbb{C}:d\symbf{\varepsilon}\) is the Cauchy stress tensor, and \(\mathbb{C}\) is the material stiffness tensor.

The left hand term in the PVD yields

\[\delta\symbf{\varepsilon}:\symbf{\sigma} = \frac{1}{2} (\mathbf{G}^I \otimes \delta \mathbf{u}_I + \delta \mathbf{u}_I \otimes \mathbf{G}^I):\symbf{\sigma} = \delta \mathbf{u}_I \cdot \symbf{\sigma} \cdot\mathbf{G}^I = \delta \mathbf{u}_I \cdot \mathbf{T}^I\]

where \(\mathbf{T}^I= \symbf{\sigma} \cdot \mathbf{G}^I\) is the traction on the edge opposing node I.

Expressing the stress as \(\symbf{\sigma}=\mathbb{C}:\symbf{\varepsilon}\) yields under consideration of symmetries in \(\mathbb{C}\)

\[\delta\symbf{\varepsilon}:\mathbb{C}:\symbf{\varepsilon} = \delta\mathbf{u}_I \left( \mathbb{C}^{IKLJ} \mathbf{G}_K\otimes\mathbf{G}_L \right) \mathbf{u}_J =: \delta\mathbf{u}_I \: \mathbb{K}^{IJ} \:\mathbf{u}_J\]

leading to the stiffness matrix as

\[\mathbf{K}^{IJ} = \int\!\!\!\int_\mathcal{A} \mathbb{K}^{IJ} \: d\mathcal{A}\]

Implementation notes

Writing the kinematic equation in cartesian component form, we obtain

\[\begin{split}\left\{\begin{array}{c} \varepsilon_{xx} \\ \varepsilon_{yy} \\ \gamma_{xy} \end{array}\right\} = \underbrace{ \left[\begin{array}{cc} G^I_x & 0 \\ 0 & G^I_y \\ G^I_y & G^I_x \end{array}\right] }_{=: [\mathbf{B}^I]} \left\{\begin{array}{c} u_{I,x} \\ u_{I,y} \end{array}\right\}\end{split}\]

leading to

\[[ \mathbb{K}^{IJ} ] = [\mathbf{B}^I]^t [\mathbb{C}] [\mathbf{B}^J]\]
Finite deformation theory

The element uses an iso-parametric mapping for both the reference and the spatial configuration. In the reference configuration,

(4)\[\mathbf{X}=\mathbf{X}_I \xi^I \quad I=0,1,2\equiv u,s,t \quad \xi^I \ge 0 ~\&~ \xi^0+\xi^1+\xi^2 = 1\]

where \(\xi^I\) are standard triangle coordinates. Alternatively, the mapping may be expressed using just two coordinates, \(\xi^1\equiv\xi^s\) and \(\xi^2\equiv\xi^t\), as

(5)\[\mathbf{X}=\mathbf{X}_0+\mathbf{G}_i \xi^i \quad i=1,2 \quad \xi^i \ge 0 ~\&~ \xi^1+\xi^2 \leq 1\]

The covariant base vectors follow from (4) and (5) as

(6)\[\mathbf{G}_i = \mathbf{X}_{,i} = \mathbf{X}_i - \mathbf{X}_0 \quad i=1,2\]

We can express differential \(d\xi^i\) in terms of \(d\mathbf{X}\) by means of the dual base as

\[d\mathbf{X}=\mathbf{G}_i d \xi^i \rightarrow d \xi^i=\mathbf{G}^i \cdot d \mathbf{X} \quad w/ ~~ \mathbf{G}^i \mathbf{G}_j=\delta_j^i\]

Similar to the above, we map the spatial configuration as

\[\mathbf{x} = \mathbf{x}_I \xi^I = \mathbf{x}_0+ \mathbf{g}_i \xi^i \quad I=0,1,2 \quad i=1,2\]

Using (6) we obtain the deformation gradient as

\[d\mathbf{x}=\mathbf{g}_i d \xi^i=\left(\mathbf{g}_i \otimes \mathbf{G}^i\right) \cdot d \mathbf{X}=\mathbf{F} \cdot d\mathbf{\bar{x}} \Longrightarrow \mathbf{F}=\mathbf{g}_i \otimes \mathbf{G}^i\]

Various strain tensors can be computed from it. We use the right Cauchy-Green deformation tensor,

\[\mathbf{C} = \mathbf{F}^t \mathbf{F}=g_{i j} \mathbf{G}^i \otimes \mathbf{G}^j\]

or the Green-Lagrange strain tensor

\[\mathbf{E}=\frac{1}{2}\left(\mathbf{F}^t \mathbf{F}-\mathbf{1}\right)=\frac{1}{2}\left(g_{i j}-G_{i j}\right) \mathbf{G}^i \otimes \mathbf{G}^j\]

Variation/linearization as needed for the principle of virtual displacements follows as

\[\begin{split}\delta \mathbf{F}=\delta\mathbf{g}_i \otimes \mathbf{G}^i =\left(\delta \mathbf{x}_i-\delta \mathbf{x}_0\right) \otimes \mathbf{G}^i =\delta \mathbf{x}_I \otimes \mathbf{G}^I \quad \underline{\text { note }}: \quad \begin{array}{l} i=1,2 \\ I=0,1,2 \end{array}\end{split}\]

where \(\mathbf{G}^0=-\mathbf{G}^1-\mathbf{G}^2\) and

\[\begin{split}\begin{align} \delta \mathbf{E} &= \frac{1}{2} \left( \delta \mathbf{F}^t \mathbf{F} + \mathbf{F}^t \delta \mathbf{F} \right) \\ &= \frac{1}{2} \left( \mathbf{G}^I \otimes \delta \mathbf{x}_I \mathbf{g}_i \otimes \mathbf{G}^i + \mathbf{G}^i \otimes \mathbf{g}_i \delta \mathbf{x}_I \otimes \mathbf{G}^I \right) \\ &= \left( \delta \mathbf{x}_I \mathbf{g}_i \right) \frac{1}{2} \left( \mathbf{G}^I \otimes \mathbf{G}^i + \mathbf{G}^i \otimes \mathbf{G}^I \right) \qquad\qquad I=0,1,2 \quad i=1,2 \\ &= \left( \delta \mathbf{x}_I \mathbf{x}_i \right) \frac{1}{2} \left( \mathbf{G}^I \otimes \mathbf{G}^i + \mathbf{G}^i \otimes \mathbf{G}^I \right) - \left( \delta \mathbf{x}_I \mathbf{x}_0 \right) \frac{1}{2} \left( \mathbf{G}^I \otimes (\mathbf{G}^1 + \mathbf{G}^2) + (\mathbf{G}^1 + \mathbf{G}^2) \otimes \mathbf{G}^I \right) \\ &= \left( \delta \mathbf{x}_I \mathbf{x}_J \right) \frac{1}{2} \left( \mathbf{G}^I \otimes \mathbf{G}^J + \mathbf{G}^J \otimes \mathbf{G}^I \right) \qquad\qquad I,J=0,1,2 \end{align}\end{split}\]

Application of virtual work

The incremental form of the principle of virtual displacements (PVD) takes the form

\[\int\!\!\!\int\limits_{\mathcal{A}}\left( \delta\mathbf{E}:\mathbb{C}:d\mathbf{E} + d\delta \mathbf{E}:\mathbf{S} \right) d\mathcal{A} = \int\!\!\!\int\limits_{\mathcal{A}} \delta\mathbf{u}\cdot\mathbf{b} \: d\mathcal{A} + \int\limits_{\partial\mathcal{A}_\sigma} \delta\mathbf{u}\cdot\bar{\mathbf{T}} \: d\mathcal{s}\]

where \(\mathbf{S}\) is the \(2^{nd}\) Piola-Kirchhoff stress tensor, \(\mathbb{C}=\partial\mathbf{S}/\partial\mathbf{E}\) is the material stiffness tensor, and

\[\begin{align} d\delta \mathbf{E} &= \left( \delta \mathbf{x}_I \:d\mathbf{x}_J \right) \frac{1}{2} \left( \mathbf{G}^I \otimes \mathbf{G}^J + \mathbf{G}^J \otimes \mathbf{G}^I \right) \qquad I,J=0,1,2 \end{align}\]

The left hand term in the PVD yields

\[\delta\mathbf{E}:\mathbb{C}:d\mathbf{E} = \delta\mathbf{x}_I \left( \mathbb{C}^{IKLJ} \mathbf{x}_K\otimes\mathbf{x}_L \right) d\mathbf{x}_J =: \delta\mathbf{x}_I \: \mathbb{K}^{IJ}_{m} \:d\mathbf{x}_J\]

and

\[d\delta\mathbf{E}:\mathbf{S} = \delta\mathbf{x}_I \: S^{IJ} \mathbf{1} \: d\mathbf{x}_J\]

leading to the tangent stiffness matrix as

\[\mathbf{K}^{IJ} = \int\!\!\!\int_\mathcal{A} \left( \mathbb{K}^{IJ}_{m} + S^{IJ} \mathbf{1} \right) d\mathcal{A}\]

Implementation notes

\[\delta\mathbf{F} =\left(\mathbf{G}^J \otimes \delta\mathbf{x}_J\right)\left(\mathbf{g}_i \otimes \mathbf{G}^i\right) =\delta\mathbf{x}_J \cdot \mathbf{g}_i \otimes \mathbf{G}^J \otimes \mathbf{G}^i =:\delta\mathbf{x}_J \cdot{\mathbf{B^{J}}^t}\]

Parent class

See also

Class doc

class femedu.elements.linear.Triangle.Triangle(node0, node1, node2, material, label=None)

class: representing a 3-noded plane triangle

computeSurfaceLoads()

compute surface loads using faces

This method should be called during updateState() by every element supporting surface loads

mapGaussPoints(var)

Initiate mapping of Gauss-point values to nodes. This method is an internal method and should not be called by the user. Calling that method explicitly will cause faulty nodal values.

Parameters:

var – variable code for a variable to be mapped from Gauss-points to nodes

resetLoads()

default implementation for resetting element loads.

setSurfaceLoad(face, pn, ps=0)

face ID

nodes defining that face

0

node 0 to node 1

1

node 1 to node 2

2

node 2 to node 0

Parameters:
  • face – face ID for the laoded face

  • pn – magnitude of distributed normal load per area. Tension on a surface is positive.

  • ps – magnitude of distributed shear load per area. Positive direction is defined as shown in the above table.

updateState()