![]() |
v0.13.2 |
The solution of the linear elasticity problem is presented in this tutorial. Lets consider an isotropic elastic cantilever beam with prescribed gravity load as presented in Figure 1.
In order to compute displacement vector field \(\mathbf{u}\), every point in the beam has to satisfy balance of linear momentum and boundary conditions as follows:
\[ \begin{align} \label{eq:momentum_balance} \begin{cases} \nabla \cdot \boldsymbol \sigma \left(\mathbf u(\mathbf{x}) \right) + \mathbf b =0 & \text { in } \Omega \\ \mathbf u(\mathbf{x}) = \bar{\mathbf{u}} & \text { on } \partial\Omega^\mathbf{u} \\ \boldsymbol \sigma \cdot \mathbf n = \bar {\mathbf t} & \text { on } \partial\Omega^\mathbf{t} \end{cases} \end{align} \]
where \(\mathbf{b}\) are body forces, \( \bar{\mathbf{u}} \) are fixed displacements, \( \bar{\mathbf{t}} \) is traction vector and \(\mathbf{n}\) is the unit normal vector on the boundary and \( \boldsymbol \sigma \) is the Cauchy stress tensor. In case of linear elasticity, small strains are assumed whereby the stress measure is defined as:
\[ \begin{align} \boldsymbol \sigma =& \mathbb D : \frac{1}{2} (\nabla \mathbf u^\text{T} + \nabla \mathbf u) \end{align} \]
where \( \mathbb D \) is the 4th order elasticity tensor.
Following a standard Galerkin approach as in SCL-1: Poisson's equation (homogeneous BC), the governing equations are discretised into weak form as follows:
\[ \begin{align} \int_{\Omega} \nabla \delta \mathbf u : \boldsymbol \sigma(\mathbf u) \, \mathrm d \Omega= - \int_{\Omega} \delta \mathbf u \cdot \mathbf b \, \mathrm d \Omega + \int_{\partial{\Omega}^\mathbf t} \delta \mathbf u \cdot \bar{\mathbf t} \, \mathrm d S \quad \forall\, \delta \mathbf u \in \mathbf{H}^1_0(\Omega) \end{align} \]
To approximate displacements \( \mathbf u\) and variance of displacements \( \delta \mathbf u\), once again we will use MoFEM's hierarchical shape functions:
\[ \begin{align} \mathbf u \approx \mathbf u^{h}=\sum_{j=0}^{N-1} \phi_{j} \bar{\mathbf{u}}_{j}, \quad \quad \delta \mathbf u \approx \delta \mathbf u^{h}=\sum_{j=0}^{N-1} \phi_{j} \delta \bar{\mathbf{u}}_{j} \end{align} \]
By substituting the above definitions into the weak form, the global system \( \mathbf{K U}=\mathbf{F}\) can be assembled, where elements contributions from every finite elements are computed as follows:
\[ \begin{align} \mathbf K_{i j}^{e}&=\int_{\Omega^{e}} \nabla \phi_{i} \mathbb D \nabla \phi_{j} \, \mathrm d \Omega \\ \mathbf F_{i}^{e}&=\int_{\partial \Omega^{e}} \phi_{i} \mathbf t \, \mathrm d S + \int_{\Omega^{e}} \phi_{i} \mathbf b \, \mathrm d \Omega \end{align} \]
The Gauss quadrature is utilised to calculate the integrals.
In this section we will focus on the procedure of setting up and solving vector-valued problems like elasticity in MoFEM using MoFEM::Simple and Form Integrators, which are template UDOs used for common operations.
The code developed for this example can be compiled for both 2D and 3D cases by changing the following constant variable:
The constant variable SPACE_DIM will be used throughout the code to determine dimension of the executable. It is also possible to setup CMakeFiles.txt to compile both, 2D and 3D executables (see plastic.cpp example).
Note that based on SPACE_DIM an appropriate type of elements is defined in the following snippet:
The main program consists of a series of functions when invoking the function Example::runProblem():
where Example::readMesh() executes the mesh reading procedure.
The function Example::setupProblem() involves adding the displacement field U on the entities of the mesh and determine the base function space, base, field rank and base order. We are adding U as both domain and a boundary field.
The space chosen is H1
, the base is AINSWORTH_LEGENDRE_BASE
and dimension is SPACE_DIM that corresponds to a vector field with SPACE_DIM number of coefficients. Additionally, we define a command line parameter -order that specifies the order of our base, 2nd order is chosen as a default.
The next function in our main program: Example::createCommonData() defines and initialises the data passed between UDOs. In case of elasticity, we compute elasticity tensor \( \mathbb D\) using the formula below:
\[ \begin{align} \mathbb D_{ijkl} = G \left[ \delta_{ik} \delta_{jl} + \delta_{il} \delta_{jk} \right] + A (K - \frac{2}{3} G) \left[ \delta_{ij} \delta_{kl} \right] \end{align} \]
where \( K\) and \( G\) are bulk and shear modulus, respectively. The coefficient \( A\) depends on the dimension of the problem, for 3D cases and plane strain formulation it is simply \( A = 1\), whereas for plane stress it takes the following form:
\[ \begin{align} A=\frac{2 G}{K+\frac{4}{3} G} \end{align} \]
Using the FTensor::Index type these formulas can be directly implemented using Einstein's summation convention as shown below:
Note that depending on SPACE_DIM our code will appropriately sum over all the indices and set coefficient \( A \) for plane stress or 3D case. In the next lines of the program, we define a vector that will be passed to the UDOs and specify the gravity load. For this example, we will hardcode the load to be always 1 in the Y direction \( \mathbf{ b } = [ 0,-1,0 ]\)
It is also possible to implement passing that data from command line parameter or read it from the mesh blocksets.
Finally, Example::createCommonData() we initialise containers that will be storing data between different UDOs. For elasticity problem we will be storing gradient of deformation, small strain tensor, stress tensor, elasticity tensor and body force vector.
In the last two lines we call previously defined lambda functions for calculating elasticity tensor and body force vector.
In order to apply homogeneous Dirichlet boundary conditions for our problem, we will utilise MoFEM's utility to remove degrees of freedom from a problem.
The first lambda function in the snippet above fix_disp searches for entities in a meshset with a given name. Subsequently, these entities (Range) are passed to another lambda function - remove_ents that removes them from the problem. This operation is equivalent of setting all the degrees of freedom to zero on associated entities. For the examples considered in this tutorial, the provided meshes contain only one meshset FIX_ALL that will essentially fix displacements in all directions to zero. In Figure 2 we show selected triangles on the face of triangle that are used to define the FIX_ALL meshset. Note that in case of 2D mesh, we select edge elements.
For natural boundary conditions we will only specify simple gravity load on the entire domain. We have to push into our Domain RHS pipeline an operator that calculates the following integral:
\[ \begin{align} \mathbf F_{i}^{e}&=\int_{\Omega^{e}} \phi_{i} \mathbf b \, \mathrm d V \end{align} \]
This operation can be done by pushing OpBodyForce
User Data Operator into RHS domain pipeline. Note that arguments of this operator are: name of the field U, container bodyForceMatPtr for vector of body force and lambda function which is defines a coefficient for our vector. The function takes three arguments, which are the coordinates of the integration points. However, in this example, we assume homogeneous body, therefore the coefficient can by simply set to 1.
The OpBodyForce
is an alias of the specialisation of the form integrator template OpBaseTimesVector
declared in the beginning of the file
The first parameter in the above template specifies the dimension of the base (we use scalar base, so set to 1), second parameter requires specification of the number of the coefficients of we field we are operating on (our vector field of displacements has SPACE_DIM coefficients). The last parameter specifies whether the passed vector is constant for all integration points. In this case bodyForceMatPtr
body force vector is assumed constant, (has the same values for all gauss points). However, if the considered body is heterogeneous, then the last parameter would be set to 1.
Another function of the main program Example::assembleSystem() is responsible for defining pipelines used to evaluate linear and bilinear forms of the system and assembly to the global matrices and vectors.
In case of linear elasticity, after pushing boundary and source term operators (gravity load), we only have to define UDOs for LHS pipeline of the domain. In case of 2D problem we push operators (MoFEM::OpCalculateHOJacForFace) that calculate and set the inverse of Jacobian for 2D elements.
The last operator OpK
integrates the bilinear form:
\[ \begin{align} \mathbf K^e_{ij} = \int_{\Omega^e} \nabla \phi_i \mathbb D \nabla \phi_j \, \mathrm d \Omega \end{align} \]
and subsequently assembles it into global system matrix \( \mathbf K \). Note that once again, a new operator does not have to be defined as OpK
is simply an alias to existing operator from MoFEM's repertoire. To that operator we are passing previously defined container that stores elasticity tensor \( \mathbb D \).
The utilised OpGradSymTensorGrad
form integrator template of Bilinear Form has four parameters: rank of the base, rank of row field, rank of the column field and finally, the increment of elasticity tensor, which in case of homogeneous material can be set to 0.
More info about available Form Integrators can be found in User data operators table.
With appropriately defined pipelines for assembling the global system of equations \( \mathbf{K} \mathbf{U} = \mathbf{F} \), the next functions Example::solveSystem() sets up PETSc solver and calculates it directly or iteratively based on specified settings.
Note that in case of simple linear problems like elasticity considered therein the above snippet will look very similar. A detailed description of each function can be found e.g. in SCL-1: Poisson's equation (homogeneous BC).
Finally, once the solver calculations are completed and solution vector is obtained, it is necessary save the data on the visualisation mesh and compute quantities of interest like strains and stresses. In function Example::outputResults() we will define pipeline to postprocess the results and save the calculated values to the output mesh.
We can reuse LHS domain pipeline, for the purpose of postprocessing. The following line resets (cleans) the operators pushed previously
In the next two lines we create a postprocessig element (object) and use its function to generate a reference mesh.
Next, a series of operators is pushed:
OpCalculateHOJacForFace
computes and sets the Jacobian in case of 2D elementsOpCalculateVectorFieldGradient
computes the gradient of the obtained displacement field \( \nabla \mathbf u^{h}=\sum_{m=0}^{N-1} \nabla \phi_m \bar{\mathbf u}_m \)OpSymmetrizeTensor
symmetrises the previously computed gradient to compute small strain tensor \( \boldsymbol \varepsilon = \frac{1}{2} (\nabla \mathbf u^\text{T} + \nabla \mathbf u) \)OpTensorTimesSymmetricTensor
is another form integrator that computes Cauchy stress: \( \boldsymbol \sigma = \mathbb D : \boldsymbol \varepsilon \)In order to run the program that we have been discussing in this tutorial, you will do the following steps
elastic
is located. Depending on how you install MoFEM shown in this page Installation, going to the directory would be something similar to thisThis will produce out_elastic.h5m file which we can convert using
The file out_elastic_3D.vtk can be opened in Paraview. Using the filter WarpByScalar
allows to visualise the deformation as shown in Figure 3.
We can also compute 2D case of the beam. In order to do that, change the below variable:
to the value equal to 2 and compile the code (in your build directory) using make
command. You should be able then to compute 2D variant as follows
and generate vtk file using:
Note that the beams, both for 2D and 3D cases experience spurious deformation at the free end. This is a known feature of small strain kinematic description. To mitigate this issue, a large strain geometrically nonlinear formulation can be used. The extension of the above tutorial to such nonlinear case can be found in VEC-2: Nonlinear elastic. We propose there a simple Hencky strain measure that improves the behaviour of the material experiencing large rotations.
The full source code of the main program can be accessed below: elastic.cpp