v0.6.9
Time dependent nonlinear mix formulation (unsaturated flow)

Table of Contents

In this tutorial, we show how to solve a strongly nonlinear time dependent equation using mix-finite-elements. We focus attention on construction of finite elements, finite elements users data operators (UDO), setting up discrete manager (DM) and time solver (TS). Issues related to the nonlinear solver and line searchers are briefly discussed. Also, we show how to use automatic differentiation (ADOL-C) to calculate derivatives of constitutive model.

Problem solved here is similar to the solution of nonlinear Poisson's equation described in A nonlinear Poisson equation and is a generalization of the problem described in Mix formulation and integration on skeleton (h-adaptivity). You might look to those tutorials if some implementation aspects of presented methodology are unclear. Also, you are encouraged to pose questions on our Q&A forum.

Presented methodology is applied to the problem of unsaturated water transport. We have no intention to give a detailed description of transport in unsaturated materials, here we only focus on problem discretization, equations linearization and implementation of finite element. Also, we do not show here how to do hp-adaptivity, this is part of other tutorials. This problem can be solved with block solver, improving its efficiency and robustness, see example Using fieldsplit solver and DM sub problem.. If you like to extend this work, adding block solver or hp-adaptivity, we offer our guidance and help.

The application for unsaturated flow is implemented in unsaturated_transport.cpp, material models are implemented in MaterialUnsaturatedFlow.hpp, problem definition and implementation is implemented in UnsaturatedFlow.hpp.

The tutorial is organized as follows. First, we explain how to run the application. Next, we show how to material model is implemented, after that we explain the mathematical and numerical model. Finally, we explain implementation details.

Note
Not all information presented here is essential, some parts you can skip and come back to them later. For example, if you like to add your material model you have to focus attention on section Physical equations.

Running code

The easiest installation is with Docker containers; please see youtube video how to run installation link. You can also look here Installation with Docker. Open terminal and run Docker container

docker run --rm=true -it --volumes-from mofem_build -v $HOME/mofem-cephas/mofem:/mofem -v $HOME:$HOME -e HOSTHOME=$HOME mofem_build /bin/bash

You are now in Docker container, go to directory where users modules are located and compile code

cd /mofem_build/um
make -j2

It can take some time when you compile it the first time since all parts of the code are compiled.

Numerical example 1

You can imagine that we do an experiment, such that we have wet silt and dry clay, we make three layers, with silt in the middle, and clay on the top and below. Such column of soil we put in the box, for example, made for plexiglass. Thus all sides are impermeable. We will observe water content and capillary pressure for two days, that is long enough to approach equilibrium.

Mesh

First, we create a mesh; we can make it in Salome, gMesh, Tetgen, or any other code. MoFEM uses MOAB for mesh database with can read various mesh file formats. Here we using Cubit to make mesh, journal file this problem is as follows

reset
# create body volume
brick x 0.01 y 0.002 z 0.1
move Volume all x 0 y 0 z -0.05 include_merged
brick x 0.1 y 0.06 z 0.02
move Volume 2  location volume 1  include_merged
chop volume 1 with volume 2
imprint volume all
merge volume all
# make block for clay
block 1 volume 4 5
block 1 name 'SOIL_CLAY1'
# make block for silt
block 2 volume 3
block 2 name 'SOIL_SILT2'
# make mesh
volume all size auto factor 10
volume all scheme Tetmesh
mesh volume all
# refine mesh at interfaces between clay & silt
refine surface 13 14 size 0.0025 bias 1.5 depth 1 smooth

Note that we do not set any boundary conditions here, in that case, is assumed that fluxes on the boundary are zero, that is boundary is impermeable.

Model parameters

We need to make a configuration file, where all problem specific parameters are set. We discuss details how of that file in mix_us_flow_physical_equation_impl

[mat_block_1]
# VanGenuchten
# SimpleDarcy
material_name=VanGenuchten
# Default material for VanGenuchten is Clay (no need to set parameters if clay is used)
# Assumed units are in meters and days.

# Model parameters controlling convergence
sCale=1e6       # Scale mass conservation equation
ePsilon1=1e-5   # Minimal capacity when removing spurious oscillations in time

Ah=-9 # head suction at z = 0
AhZ=1 # gradient of initial head suction
AhZZ=0 # initial head suction is calculated from equation h = Ah + AhZZ*z + AhZZ*zz

[mat_block_2]
material_name=VanGenuchten
# Material parameters for Silt from Vogel., van Genuchten Cislerova 2001
sCale=1e6
ePsilon1=1e-5
Ah=-0.09
AhZ=1
thetaS=0.46
thetaM=0.46
thetaR=0.034
alpha=1.6
n=1.37
Ks=0.006
hS=0

The material parameters are taken from [44], the initial parameters for clay (mat_block_1) and silt (mat_block_2) are set with the \(A_h\) and \(A_h^z\), such that initial capillary pressure head is given by \( h = A_h+A_h^z z \) where

\[ \left\{ \begin{array}{ll} h(z,t=0) = -9+z&\textrm{for clay}\\ h(z,t=0) = -0.09+z&\textrm{for silt} \end{array} \right. \]

Analysis

The analysis is run using script

#!/bin/bash
# Get file name
if [ -z ${1+x} ]; then
MESHFILE=1dTest.cub
else
MESHFILE=$1
fi
# Get config file name
if [ -z ${2+x} ]; then
CONFIGFILE=unsaturated.cfg
else
CONFIGFILE=$2
fi
# Get numbet of processors
if [ -z ${3+x} ]; then
NBPROCS=4
else
NBPROCS=$3
fi
# Get step size
if [ -z ${4+x} ]; then
DT=0.0001
else
DT=$4
fi
# Get final time
if [ -z ${5+x} ]; then
FT=1.0
else
FT=$5
fi
# Get final time
if [ -z ${6+x} ]; then
ORDER=0
else
ORDER=$6
fi
# Partition mesh
../../tools/mofem_part -my_file $MESHFILE -meshsets_config $CONFIGFILE -my_nparts $NBPROCS
# Run code
rm -f out_*.h5m && \
make -j 4 unsaturated_transport && \
mpirun -np $NBPROCS ./unsaturated_transport \
-my_file out.h5m -configure $CONFIGFILE \
-ts_monitor -ts_type beuler \
-ts_dt $DT -ts_final_time $FT \
-ts_monitor -ts_adapt_monitor \
-ts_adapt_type basic \
-ts_adapt_dt_max 0.005 \
-ts_adapt_dt_min 0.0001 \
-ts_rtol 0.01 -ts_atol 0.01 \
-ts_adapt_basic_safety 0.8 \
-ts_error_if_step_fails false \
-ts_theta_adapt false \
-ts_alpha_adapt false \
-ts_max_reject -1 \
-ts_max_snes_failures -1 \
-ksp_type gmres -pc_type lu -pc_factor_mat_solver_package mumps \
-snes_type newtonls \
-snes_linesearch_type l2 \
-snes_linesearch_minlambda 1e-3 -snes_linesearch_damping 1. -snes_linesearch_max_it 3 \
-snes_atol 1e-7 -snes_rtol 1e-4 -snes_stol 0 -snes_max_it 32 \
-snes_converged_reason \
-my_order $ORDER \
-how_often_output 5 \
-my_max_post_proc_ref_level 0 2>&1 | tee log && \
rm -f out_*.vtk && \
../nonlinear_elasticity/do_vtk.sh out_*h5m

. Editing that file you can set range of parameters, controlling time solve, tolerances, stopping criteria, line-searcher. For details look to PETSc manual link.

Note solved problem is strongly nonlinear, and we use here secant line-search in the L2 norm (see link). For stability of results at each Newton iteration, we execute three sub-iterations with line searcher. Also, you have option controlling time step adaptively and linear solver pre-conditioner. All command line options, both PETSc and MoFEM are printed if you add switch -help when you run MoFEM application.

MoFEM specific options are

Finally, we can run simulation, be executing above script

cd /mofem_build/um/basic_finite_elements/mix_transport
./run_uf.sh soil_impermeable.cub unsaturated.cfg 4 0.001 1 0

where

Results

For problem set-up like above, results are shown in figure 1, where the evolution of various parameters over time until equilibrium is reached is shown. You can note that initially, clay is relatively dry, whereas silt is wet. Over time clay suck the water from silt. On the left-hand side observed the distribution of saturation. Initially, high gradients of pressure head are created on the interface between two materials. Those are captured very well, despite coarse mesh, showing capabilities of presented problem formulation. For creating of this figure we used the 2nd order of polynomial to approximate pressure head and 3rd order polynomial to approximate fluxes.

unsaturated_flow_example_fig3.gif
Figure 1. On the bottom axis is depth, on the left axis is pressure head, right axis is water content. Red line is water head, green line is water content. Colormap on the left hand side is for saturation.

Numerical example 2

In the second example, we show how to run downward infiltration problem with three layers, i.e. clay, silt and clay.

Mesh

In that case on the top, we apply zero pressure head, rest of the boundary is impermeable. See journal script to generate such mesh

reset
# create body volume
brick x 0.01 y 0.002 z 0.1
move Volume all x 0 y 0 z -0.05 include_merged
brick x 0.1 y 0.06 z 0.02
move Volume 2  location volume 1  include_merged
chop volume 1 with volume 2
imprint volume all
merge volume all
# make block for clay
block 1 volume 4 5
block 1 name 'SOIL_CLAY1'
# make block for silt
block 2 volume 3
block 2 name 'SOIL_SILT2'
block 3 surface 1
block 3 name 'HEAD1'
block 3 attribute count 1
block 3 attribute index 1 0
# make mesh
volume all size auto factor 10
volume all scheme Tetmesh
mesh volume all
# refine mesh at interfaces between clay & silt
refine surface 1 13 14 size 0.0025 bias 1.5 depth 1 smooth

Here we add third block name HEAD1, where boundary condition is applied. Similarly, non-zero fluxes can be applied by making a block with some surface elements and called it for example "FLUX1".

Model parameters

The material parameters are similar to one given in the previous problem

[mat_block_1]
# VanGenuchten
# SimpleDarcy
material_name=VanGenuchten
# Default material for VanGenuchten is Clay (no need to set parameters if clay is used)
# Assumed units are in meters and days.

# Model parameters controlling convergence
sCale=1e6       # Scale mass conservation equation
ePsilon1=1e-5   # Minimal capacity when removing spurious oscillations in time

Ah=-9 # head suction at z = 0
AhZ=1 # gradient of initial head suction
AhZZ=0 # initial head suction is calculated from equation h = Ah + AhZZ*z + AhZZ*zz

[mat_block_2]
material_name=VanGenuchten
# Material parameters for Silt from Vogel., van Genuchten Cislerova 2001
sCale=1e6
ePsilon1=1e-5
Ah=-9
AhZ=1
thetaS=0.46
thetaM=0.46
thetaR=0.034
alpha=1.6
n=1.37
Ks=0.006
hS=0

# You can change boundary conditions here
#[head_block_3]
#head = 0

Note for this example we set initial pressure head distribution such that clay and silt are in equilibrium.

\[ \left\{ \begin{array}{ll} h(z,t=0) = -9+z&\textrm{for clay}\\ h(z,t=0) = -9+z&\textrm{for silt} \end{array} \right. \]

Moreover at the end we show how to modify boundary condition. As an excrete one can change value of head on the top of the sample.

Analysis

We run script as before

./run_uf.sh upward_infiltration.cub unsaturated.cfg 4 0.001 1 0

Results

For this example we observe an influx of water into the column of solid, we can plot how it is changing over time. At every time step solver monitor printing on the screen following lines

63 TS dt 0.00242774 time 0.0710581
Flux at time 0.07106 3.859e-07
    0 SNES Function norm 4.123926528377e-02
    1 SNES Function norm 6.552371691866e-02
    2 SNES Function norm 7.740255184123e-02
    3 SNES Function norm 7.999178179888e-02
    4 SNES Function norm 7.881110800442e-02
    5 SNES Function norm 7.582982118207e-02
    6 SNES Function norm 7.173427725809e-02
    7 SNES Function norm 6.666200930695e-02
    8 SNES Function norm 6.040633270649e-02
    9 SNES Function norm 5.232159632926e-02
   10 SNES Function norm 4.057267597689e-02
   11 SNES Function norm 1.625005251310e-02
   12 SNES Function norm 4.316837250846e-05
   13 SNES Function norm 8.980432593466e-07
  Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 13
      TSAdapt 'basic': step  63 accepted t=0.0710581  + 2.428e-03 wlte=0.546 family='beuler' scheme=0:'' dt=2.627e-03
64 TS dt 0.00262737 time 0.0734859
Flux at time 0.07349 3.812e-07
    0 SNES Function norm 4.023589003253e-02
    1 SNES Function norm 6.962486297267e-02
    2 SNES Function norm 8.329691099303e-02
    3 SNES Function norm 8.641640743250e-02
    4 SNES Function norm 8.532070311444e-02
    5 SNES Function norm 8.223165652333e-02
    6 SNES Function norm 7.792739113001e-02
    7 SNES Function norm 7.257891287567e-02
    8 SNES Function norm 6.598809018682e-02
    9 SNES Function norm 5.751011198160e-02
   10 SNES Function norm 4.535116560548e-02
   11 SNES Function norm 2.147388609216e-02
   12 SNES Function norm 1.073710142541e-05
   13 SNES Function norm 4.011739514453e-11
  Nonlinear solve converged due to CONVERGED_FNORM_ABS iterations 13
      TSAdapt 'basic': step  64 accepted t=0.0734859  + 2.627e-03 wlte=0.588 family='beuler' scheme=0:'' dt=2.742e-03

where for step 63 and 64, total flux is printed in \(m^3/day\),

Flux at time 0.07106 3.859e-07
Flux at time 0.07349 3.812e-07

We can extract that information using awk/grep as follows

awk '/Flux at/ {print $4,$5/(0.01*0.002)}' log | tee gp

This shell command line takes all lines which have the phrase "Flux at", takes fourth and the fifth column from those lines, whereas fifth column is divided by area of the surface where the head is applied. Finally, results are stored in file gp. You can plot results in Excel or any other software; we choose to use gnuplot. As results, we get figure below.

unsaturated_flow_example_fig5.png
Figure 3. Influx of water for upward infiltration.

On figure below, we can observe how wetting front is moving downwards. Note how the gradient of head differ between clay and silt, and how it moves through the interface between materials.

unsaturated_flow_example_fig4.gif
Figure 2. On the bottom axis is depth, on the left axis is pressure head, right axis is water content. Red line is water head, green line is water content. Colormap on the left hand side is for saturation.

Mathematical model

W considered porous material partially filled with water. In this tutorial, we going to simulate water transport in such system. The primary behaviour of the system depends on the material structure, capillary forces and wetting angle, taken into count by experimentally estimated material parameters. Darcy's law gives the relation between flux and suction head (capillary pressure). Constitutive equations describe the relationship between hydraulic conductivity and water content. It is assumed that state of the material is a function of water content, and using water retention curve one can map between capillary head and water content, i.e. \(\theta=\theta(h)\). The process of water flow has to obey conservation of mass, and consequently, process of unsaturated transport is described by the parabolic differential equation.

PDE

The Richards equation describes problem of unsaturated flow

\[ \frac{\partial \theta}{\partial t}- \nabla \cdot \left[ \left( K(\theta) \nabla ( h - z) \right) \right] = 0 \]

where \(z\) is depth and \(K=K(\theta)\) is hydraulic conductivity. \(\nabla\) and \(\nabla \cdot\) are differential operators for gradient and divergence, respectively. The above equation alternatively can be expressed by the system of first order partial differential equations

\[ \left\{ \begin{array}{l} \frac{1}{K(\theta)}\boldsymbol\sigma + \nabla (h -z) = 0 \\ c(\theta)\frac{\partial h}{\partial t} - \nabla \cdot \boldsymbol\sigma = 0 \end{array} \right. \]

where \(\boldsymbol\sigma\) is flux. First equation expresses constitutive equation; second is conservation law. The second equation is called as well continuity equation, and it does not allow for creation or annihilation of mass. The capacity term \(c\) in the second question is given by

\[ c = \frac{\partial \theta}{\partial h} \]

Physical equations

It exists semi-empirical relation between \(\theta\) and \(h\), called water retention curve, and for unsaturated state one can uniquely map between one and another when water retention hysteresis is not included. For simplicity here, we do not include hysteresis in our equations, and relations need to be established,

\[ \left\{ \begin{array}{l} \theta=\theta(h)\\ K = K(\theta) \end{array} \right. \]

Generic material model is implemented in MixTransport::GenericMaterial, which is an abstract class, thus it instance can not be created. This class sets methods and data used by finite element instances. The is not modified when new material is implemented. This class is located in UnsaturatedFlow.hpp. To implement new material model user need to implement method

Form MixTransport::GenericMaterial, the MixTransport::CommonMaterialData is derived. This class is unchanged if new material is implemented unless new material parameters need to be added to the model. Material parameters are read from text file given with command line option -configure unsaturated.cfg. Example of such fill looks as follows,

[mat_block_1]
# VanGenuchten
# SimpleDarcy
material_name=VanGenuchten
# Default material for VanGenuchten is Clay (no need to set paramters if clay is used)
# Assumed usinits are in meters and days.

# Model parameters controlling convergence
sCale=1e6       # Scale mass conservation equation
ePsilon1=1e-5   # Minimal capacity when removing spurious oscillations in time

Ah=-9 # head suction at z = 0
AhZ=1 # gradient of initial head suction
AhZZ=0 # initial head suction is calculated from equation h = Ah + AhZZ*z + AhZZ*zz

[mat_block_2]
material_name=VanGenuchten
# Material parameters for Silt from Vogel., van Genuchten Cislerova 2001
sCale=1e6
ePsilon1=1e-5
Ah=-0.09
AhZ=1
thetaS=0.46
thetaM=0.46
thetaR=0.034
alpha=1.6
n=1.37
Ks=0.006
hS=0

# You can change boundary conditions here
#[head_block_3]
#head = 0

Note that we have to blocks, which are sets of elements on the mesh. To each block, we can attach material and set material parameters. A user can add more material parameters, by declaring them in MixTransport::CommonMaterialData and add them to method MixTransport::CommonMaterialData::addOptions, such that added material parameter are read from the configuration file.

For proposes of this tutorial we have implemented two material models, i.e. MixTransport::MaterialDarcy and MixTransport::MaterialVanGenuchten, both classes derived from MixTransport::CommonMaterialData.

Adding new material

The simplest method to implement new material is to derive class form MixTransport::MaterialWithAutomaticDifferentiation, and add it to file MaterialUnsaturatedFlow.cpp. MixTransport::MaterialWithAutomaticDifferentiation is a class to simplify implementation of the material model by utilizing advantages of automatic differentiation. Comprehensive description of library ADOL-C for animatic differentiation can be found here (link).

In principle automatic differentiation record mathematical operations in the function. Recorded mathematical operations are represented algorithmically in the tree, where each tree leaf is expanded into Taylor series, with sufficient order to calculate derivatives exactly. Automatic differential is different from symbolic differentiation, such that it can calculate derivatives of algorithm, so you can have loops and other commands in your differentiated function.

In consequence, user need to build a tree to calculate water content density \(\theta=\theta(h)\) and relative hydraulic conductivity \(K_r=K_r(h)\). Using terminology from ADOL-C, one has to record tape with differentiated function. This is done by overloading two functions, MaterialWithAutomaticDifferentiation::recordTheta and MaterialWithAutomaticDifferentiation::recordKr. For example

struct MyMaterialForUnsaturatedFlow: public MaterialWithAutomaticDifferentiation {
static boost::shared_ptr<CommonMaterialData> createMatPtr(const CommonMaterialData &data) {
return boost::shared_ptr<CommonMaterialData>(new MyMaterialForUnsaturatedFlow(data));
}
MyMaterialForUnsaturatedFlow(const CommonMaterialData &data):
MaterialWithAutomaticDifferentiation(data) {
recordTheta();
recordKr();
}
adouble ah;
adouble aTheta;
adouble aKr;
adouble aSe;
adouble aSeStar;
template <typename TYPE>
inline TYPE funTheta(TYPE &h,const double m) {
return thetaR+(thetaS-thetaR)/pow(1+pow(-alpha*h,n),m);
}
virtual void recordTheta() {
trace_on(2*blockId+0,true);
h = -1-hS;
ah <<= h;
const double m = 1-1/n;
aTheta = funTheta(ah,m);
double r_theta;
aTheta >>= r_theta;
trace_off();
}
virtual void recordKr() {
trace_on(2*blockId+1,true);
h = -1-hS;
ah <<= h;
const double m = 1-1/n;
aTheta = funTheta(ah,m);
aSe = funSe(aTheta);
aKr = aSe*pow(1-pow(1-pow(Se,1/m),m),2);
double r_Kr;
aKr >>= r_Kr;
trace_off();
}
};

where method recordTheta, represents equation

\[ \theta = \theta_r + \frac{\theta_s-\theta_r}{(1+(-\alpha h)^n)^m} \]

and method recordKr represents equation

\[ K_r(S_e) = S_e\left[1-(1-S_e^{1/m})^m\right]^2 \]

Note that \(h\) is independent variable, whereas \(\theta\) and \(K_r\) are dependent variables. \(\alpha\), \(n\) and \(m\) are model parameters, see for details here [44]. Dissecting method recordTheta, we have

trace_on(2*blockId+0,true);
/// Code here
trace_off();

where we start and stop tape recording/tracing. Each tape has unique Id, in our case we have a convention that even Id's are used to record operations to calculate \(\theta\) and odd Id's are used to record \(K_r\). Each material block has own set of tapes since it uses unique material parameters and physical model has calculated from block number set on the mesh. In following code

h = -1-hS;
ah <<= h;
// Calculate aTheta here
double r_theta;
aTheta >>= r_theta;

we start by seeding point for which we build tree, \(h\) can be arbitrary but for the unsaturated state. Next, using "<<=" we set indented variable and using operator ">>=" to set dependent variable.

Finally new material model need to be registered, by adding line to the following code as follows

struct RegisterMaterials {
static map<std::string,CommonMaterialData::RegisterHook> mapOfRegistredMaterials;
MoFEMErrorCode operator()() const {
mapOfRegistredMaterials["SimpleDarcy"] = MaterialDarcy::createMatPtr;
mapOfRegistredMaterials["VanGenuchten"] = MaterialVanGenuchten::createMatPtr;
mapOfRegistredMaterials["MyModel"] = MyMaterialForUnsaturatedFlow::createMatPtr;
}
};

Numerical model

In following we apply mix-formulation, building on tutorial Mix formulation and integration on skeleton (h-adaptivity). We will show below that such approximation enables to solve Richards equation efficiently. This enables to address some issue of the regularity of equation difficult to approximate for classical formulation of finite elements.

Mix-formulations allows converging to the solution faster than classical single filed finite element formulation. Mix-elements have built in error estimators into formulation letting for efficient hp-adaptation and least not least separation of nonlinearities, driving the development of effect solvers for strongly nonlinear problems.

Conservation of mass

Our PDE is system of two equations. The second equation is universal law of conservation of mass, which enforce continuity of fluxes, such that for arbitrary surface \(\Gamma\) dividing body on two parts, that is \(\Omega = \Omega^1 \cup \Omega^2\)

\[ \mathbf{n}(\mathbf{x}) \cdot \left( \boldsymbol\sigma^1(\mathbf{x}) - \boldsymbol\sigma^2(\mathbf{x}) \right) = 0,\;\forall\mathbf{x}\in\Gamma \]

Utilizing this observation we choose approximation base which a priori satisfy such continuity on element faces. However, it not enforce continuity of tangential components on the element face, since those do not contribute to mass exchange on element faces. Such approximation base is constructed for Raviart Thomas elements. Here we are using the recipe from [22] to construct arbitrary polynomial base order on elements, which a priori enforce continuity of fluxes.

The mass conservation in discrete form can be expressed as follows

\[ (v,c(\theta) \frac{\partial h}{\partial t})_\Omega - (v,\boldsymbol\sigma)_\Omega = 0 \quad \forall v \]

We obtain this equation by multiplying mass conservation by test function \(v\) and integrating both sides. Tested and test base functions are \(h,v \in L^2(\Omega)\), whereas \(\boldsymbol\sigma \in H(\nabla\cdot,\Omega)\) are divergence integrable functions which property explained in above paragraph. Here, we use notation, suitable for problems with many terms in the variational formulations, requires some explanation. First, we use the shorthand notation

\[ (v,c(\theta) \frac{\partial h}{\partial t})_\Omega = \int_\Omega v c(\theta) \frac{\partial h}{\partial t} \textrm{d}\Omega \]

and analogically for other terms.

Physical equation

The physical equation is the source of strong nonlinearities and problems with convergence of nonlinear solver and irregularities in the solution. We note that continuity of capillary head (capillary pressure) is not derived from conservation law, but dictated by the regularity of hydraulic conductivity \(K=K(\theta)\).

The philosophy of finite element formulation presented here is to use approximation base having a priory continuity demanded by conservation law, whereas not to enforce continuity of capillary pressure head, which is controlled by the regularity of a constitutive law. That allows developing a numerical model for the bigger class of materials, initial and boundary conditions, suitable for analysis interface effects. Here, continuity of capillary pressure head is enforced posterior, as a consequence of assumed physical model.

Taking physical equation and testing it with test function \(\boldsymbol\tau \in H(\nabla\cdot,\Omega)\) we get

\[ (\boldsymbol\tau,\frac{1}{K(\theta)} \boldsymbol\sigma)_\Omega + (\boldsymbol\tau,\nabla (h-z))_\Omega = 0 \quad \forall \boldsymbol\tau \]

next applying integration by parts we get

\[ (\boldsymbol\tau,\frac{1}{K(\theta)} \boldsymbol\sigma)_\Omega - (\nabla \cdot \boldsymbol\tau,(h-z))_\Omega + (\mathbf{n} \cdot \boldsymbol\tau,(\overline{h}-z))_{\partial\Omega} = 0 \]

where \(\overline{h}\) is prescribed capillary pressure on the boundary. Note that

Note
This equation is irregular, hydraulic conductivity can changes sharply, as a function of the capillary head. The problem can be regularized by applying Kirchhoff transform, see details here [34]. However such approach is limited, and can be applied for special cases, in particular when no hysteresis is considered. Applying Kirchhoff transform, for steady case problem transform nonlinear problem to the linear problem, for not steady case problem has still nonlinear capacity term. Such approach enables to calculate large scale problems for unsaturated soils and is highly recommended to be implemented here as an alternative approach.

Semi-discrete linearized system of equations

In MoFEM we do not implement own Newton method, but use one implemented in PETSc. PETSc Newton solver is implemented in a package of nonlinear solvers SNES. In the nonlinear time dependent problem, SNES is called by TS (time solver). Here we focus attention on the calculation of residual vector and tangent matrix.

In this section, we focus attention on discretization in space, without restricting yourself to the particular time integration scheme. That is exclusively managed by TS solver from PETSc. This explains title of this section, the formulation is semi-discrete, since here apply discretization in space, but we assume that function in time is given by some unknown analytical function.

Above we put physical justification for using particular spaces for fluxes, which is \(\boldsymbol\tau,\boldsymbol\sigma \in H(\nabla\cdot,\Omega)\), and capillary pressure head \(u,h \in L^2(\Omega)\). One can show that selection of such approximation spaces, once the problem is discretized (i.e. finite element approximation spaces are used), numerical solution is stable since the choice of approximation spaces comply with De Rham diagram.

To solve problem we define residuals

\[ \left\{ \begin{array}{l} r_\tau = (\boldsymbol\tau,\frac{1}{K(\theta)} \boldsymbol\sigma)_\Omega - (\nabla \cdot \boldsymbol\tau,(h-z))_\Omega + (\mathbf{n} \cdot \boldsymbol\tau,(\overline{h}-z))_{\partial\Omega}\\ r_v = (v,c(\theta) \frac{\partial h}{\partial t})_\Omega - (v,\boldsymbol\sigma)_\Omega \end{array} \right. \]

Applying standard procedure where one solves nonlinear system of equations, we expand both equations into Taylor series, for terms \(\boldsymbol\sigma\), \(h\) and \(\frac{\partial h}{\partial t}\) since we have time dependent problem. Truncating the Taylor series after linear term, we get

\[ \left\{ \begin{array}{l} r^i_\tau + (\boldsymbol\tau,\frac{1}{K^i} \delta\boldsymbol\sigma^{i+1})_\Omega - (\nabla \cdot \boldsymbol\tau,\delta h^{i+1})_\Omega + (\boldsymbol\tau,\frac{\textrm{d}K}{\textrm{d}\theta}|^i\frac{\textrm{d}\theta}{\textrm{d}h}|^i(K^i)^{-2} \boldsymbol\sigma^i) \delta h^{i+1})_\Omega = 0\\ r_v^i + (v,c^i\frac{\partial \delta h}{\partial t}|^{i+1})_\Omega + (v,\frac{\textrm{d}c}{\textrm{d}\theta}|^i\frac{\textrm{d}\theta}{\textrm{d}h}|^i \frac{\partial h}{\partial t}|^i) \delta h^{i+1})_\Omega - (v,\delta\boldsymbol\sigma^{i+1})_\Omega = 0 \end{array} \right. \]

where \(i\) is Newton iteration number. Note that tangent matrices and residuals are evaluetd at known iteration \(i\) and unknowns are solved for iteration \(i+1\). Iterations are repeated until residuals are small in the absolute or relative norm. Once we achieve equilibrium we move to next time step.

Expressing test and tested function by vectors of basis functions and vectors of degrees of freedoms we have

\[ \boldsymbol\sigma^i_n(\mathbf{x}) = \boldsymbol \phi^\textrm{T}(\mathbf{x}) \mathbf{q}^i_n,\; h^i_n(\mathbf{x},t_n) = \boldsymbol \psi^\textrm{T}(\mathbf{x},t_n) \mathbf{p}^i_n \]

where base vectors \(\boldsymbol \phi\) and \(\boldsymbol \psi\) are for H-div and L2 space respectively. \(n\) is time step. Finally, we get system of linear equations

\[ \left[ \begin{array}{cc} \mathbf{C} & \mathbf{G} \\ \mathbf{H} & \mathbf{A} \end{array} \right] \left\{ \begin{array}{c} \delta \mathbf{q}^{i+1}\\ \delta \mathbf{p}^{i+1} \end{array} \right\} = \left[ \begin{array}{c} \mathbf{r}^i_\tau\\ \mathbf{r}^i_v \end{array} \right] \]

where vectors of unknowns are updated as follows

\[ \mathbf{q}^{i+1}_{n+1} = \mathbf{q}_n + \Delta \mathbf{q}_{n+1}^i + \delta \mathbf{q}^{i+1}\quad \mathbf{p}^{i+1}_{n+1} = \mathbf{p}_n + \Delta \mathbf{p}_{n+1}^i + \delta \mathbf{p}^{i+1} \]

Tangent matrices and vectors are implemented in following operators

Time discretization

In MoFEM we do not implement time discretization algorithms unless we like to do something novel and non-standard, and we would rather recommend implementing such an algorithm using PETSc time solver shell rather directly in MoFEM. Detailed description of how to use time solver (TS) in PETSc can be found in link, see chapter 6.

The constructions of time discretized problem are managed transparently for the user. If you are interested in details, see MoFEM::TsCtx class and auxiliary functions like MoFEM::f_TSSetIFunction, MoFEM::f_TSSetIJacobian and MoFEM::f_TSMonitorSet. Those functions are called by the discrete manager (DM) which is linked to time solver.

In the previous section, we transform PDE, by applying space discretization, into the system of linearized ordinary differential equations (ODE). In the following description, we transform the system of ODEs into a system of linear equations. In our the particular problem we need to calculate the rate of capillary head, which is discretized in time and space by

\[ h(\mathbf{x},t) = \boldsymbol \psi(\mathbf{x},t)^\textrm{T} \mathbf{p} = [\boldsymbol \xi(t).\boldsymbol{\tilde{\psi}}(\mathbf{x})]^\textrm{T} \mathbf{p} \]

where \([\xi(t).\boldsymbol{\tilde{\psi}}(\mathbf{x})]\) is element by element product of two vectors. With that at hand we can calculate rate of hydraulic head

\[ \frac{\partial h}{\partial t}(\mathbf{x},t) = \left[ \frac{\textrm{d}\boldsymbol\xi}{\textrm{d}t} (t).\boldsymbol{\tilde{\psi}}(\mathbf{x}) \right]^\textrm{T} \mathbf{p} \]

With above at hand, matrix \(\mathbf{C}\) where time derivative is present, have following form

\[ \mathbf{C} = \left(\boldsymbol\psi^\textrm{T}, \left( c^i \frac{\partial \boldsymbol\xi}{\partial t}+ \frac{\textrm{d}c}{\textrm{d}\theta}|^i\frac{\textrm{d}\theta}{\textrm{d}h}|^i \frac{\partial h}{\partial t}|^i \right) \boldsymbol{\tilde\psi}\right)_\Omega \]

Note for any ODE integration method the approximation of \(\frac{\partial h}{\partial t}\) is by linear function, hence

\[ \frac{\partial \boldsymbol\xi}{\partial t} = (shift). \]

For example, iterative change of capillary pressure head, when backwards Euler’s method is applied, \((shift) = 1/\Delta t\), what can be derived as follows

\[ \left(\frac{\partial \delta h}{\partial t}\right)^{i+1}_{n+1} = \left(\frac{\partial \delta h}{\partial \delta \mathbf{p}}\frac{\partial \delta \mathbf{p}}{\partial t} \right)^{i+1}_{n+1} = \left( \boldsymbol{\tilde \psi} \delta \mathbf{p}^{i+1}\right) \left( 1/\Delta t \right) = \left( \boldsymbol{\tilde \psi} \delta \mathbf{p}^{i+1} \right) (shift) \]

Such approach allows for code generalisation such that finite element formulation which concerns discretization in space, does not depend on discretization in time. That makes it possible to switch between time integration schemes without changing implementation, see list of time integration schemes in link, table 10.

Above is exploited in MixTransport::UnsaturatedFlowElement::OpVU_L2L2::doWork, where shift, as follows, is applied at integration points

for(int gg = 0;gg!=nb_gauss_pts;gg++) {
// get integration weight and multiply by element volume
double alpha = t_w*vol*scale;
// evaluate material model at integration points
// to calculate capacity and tangent of capacity term
block_data->h = t_h;
block_data->h_t = t_h_t;
block_data->x = t_coords(0);
block_data->y = t_coords(1);
block_data->z = t_coords(2);
ierr = block_data->calC(); CHKERRG(ierr);
ierr = block_data->calDiffC(); CHKERRG(ierr);
const double C = block_data->C;
const double diffC = block_data->diffC;
// assemble local entity tangent matrix
FTensor::Tensor0<double*> t_a(&*nN.data().begin());
// iterate base functions on rows
for(int kk = 0;kk!=nb_row;kk++) {
// get first base function on column at integration point gg
FTensor::Tensor0<double*> t_n_col = col_data.getFTensor0N(gg,0);
// iterate base functions on columns
for(int ll = 0;ll!=nb_col;ll++) {
// assemble elements of local matrix
t_a += (alpha*(C*ts_a+diffC*t_h_t))*t_n_row*t_n_col;
++t_n_col; // move to next base function on column
++t_a; // move to next element in local tangent matrix
}
++t_n_row; // move to next base function on row
}
++t_w; // move to next integration weight
++t_coords; // move to next coordinate at integration point
++t_h; // move to next capillary head at integration point
++t_h_t; // move to next capillary head rate at integration point
}

where \(\textrm{ts_a} = (shift)\). Note that evaluation of capillary head DOFs, i.e. \(\frac{\partial \mathbf{p}}{\partial t}\) is managed exclusively by PETSc and depend on used time integration shame.

Implementation

Declaration of fields, finite elements and problem

In method MixTransport::UnsaturatedFlowElement::addFields we declare, define fields on entities and set order of approximation as follows

//Fields
ierr = mField.add_field(fluxes,HDIV,DEMKOWICZ_JACOBI_BASE,1); CHKERRG(ierr);
ierr = mField.add_field(values,L2,AINSWORTH_LEGENDRE_BASE,1); CHKERRG(ierr);
ierr = mField.add_field(values+"_t",L2,AINSWORTH_LEGENDRE_BASE,1); CHKERRG(ierr);
//add entities to field
if(it->getName().compare(0,4,"SOIL")!=0) continue;
ierr = mField.add_ents_to_field_by_type(
dMatMap[it->getMeshsetId()]->tEts,MBTET,fluxes
); CHKERRG(ierr);
ierr = mField.add_ents_to_field_by_type(
dMatMap[it->getMeshsetId()]->tEts,MBTET,values
); CHKERRG(ierr);
ierr = mField.add_ents_to_field_by_type(
dMatMap[it->getMeshsetId()]->tEts,MBTET,values+"_t"
); CHKERRG(ierr);
}
ierr = mField.set_field_order(0,MBTET,fluxes,order+1); CHKERRG(ierr);
ierr = mField.set_field_order(0,MBTRI,fluxes,order+1); CHKERRG(ierr);
ierr = mField.set_field_order(0,MBTET,values,order); CHKERRG(ierr);
ierr = mField.set_field_order(0,MBTET,values+"_t",order); CHKERRG(ierr);

Note that we set field space independently from the approximation base. We can have different recipes to construct approximation base, for example, [1] and [22]. Note, to get solvable system of linear equations, pressure head approximation has one order lower that approximation for fluxes. Here for simplicity, we set homogeneous approximation order, however, implementation can be easily extended to the case where approximation order is locally increased near interfaces or boundary were pressure head is applied.

In similar way we declare and define finite elements in MixTransport::UnsaturatedFlowElement::addFiniteElements, as follows

// Define element "MIX". Note that this element will work with fluxes_name and
// values_name. This reflect bilinear form for the problem
ierr = mField.add_finite_element("MIX",MF_ZERO); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_row("MIX",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_col("MIX",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_row("MIX",values_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_col("MIX",values_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_data("MIX",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_data("MIX",values_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_data("MIX",values_name+"_t"); CHKERRG(ierr);
if(it->getName().compare(0,4,"SOIL")!=0) continue;
ierr = mField.add_ents_to_finite_element_by_type(
dMatMap[it->getMeshsetId()]->tEts,MBTET,"MIX"
); CHKERRG(ierr);
}
// Define element to integrate natural boundary conditions, i.e. set values.
ierr = mField.add_finite_element("MIX_BCVALUE",MF_ZERO); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_row("MIX_BCVALUE",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_col("MIX_BCVALUE",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_data("MIX_BCVALUE",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_data("MIX_BCVALUE",values_name); CHKERRG(ierr);
if(it->getName().compare(0,4,"HEAD")!=0) continue;
ierr = mField.add_ents_to_finite_element_by_type(
bcValueMap[it->getMeshsetId()]->eNts,MBTRI,"MIX_BCVALUE"
); CHKERRG(ierr);
}
// Define element to apply essential boundary conditions.
ierr = mField.add_finite_element("MIX_BCFLUX",MF_ZERO); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_row("MIX_BCFLUX",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_col("MIX_BCFLUX",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_data("MIX_BCFLUX",fluxes_name); CHKERRG(ierr);
ierr = mField.modify_finite_element_add_field_data("MIX_BCFLUX",values_name); CHKERRG(ierr);
if(it->getName().compare(0,4,"FLUX")!=0) continue;
ierr = mField.add_ents_to_finite_element_by_type(
bcFluxMap[it->getMeshsetId()]->eNts,MBTRI,"MIX_BCFLUX"
); CHKERRG(ierr);
}

Note that we create three elements "MIX", "MIX_BCVALUE" and "MIX_BCFLUX", to evaluate fields in volume, on the boundary where natural conditions are applied and where essential conditions are applied, respectively.

Both fields and finite elements are defined on part of the mesh, i.e. entities set, where mesh block name starts with first four letters "SOIL". In case of "MIX_BCVALUE" and "MIX_BCFLUX" where blocks name having first for letters "HEAD" and "FLUX", respectively. The user can define multiple blocks and attach material model, material parameters, initial conditions independently to each of them independently.

In a similar way, other types of elements can be added to the formulation, for example, we can create evaporation element on some part of the boundary, or one can extend formulation with other mechanical and non-mechanical fields. For example, one can add heat conduction and solve coupled problem with freezing and thawing processes.

In method MixTransport::UnsaturatedFlowElement::buildProblem we create problem and Discrete Manager operating on that problem

// create DM instance
ierr = DMCreate(PETSC_COMM_WORLD,&dM);CHKERRG(ierr);
// setting that DM is type of DMMOFEM, i.e. MOFEM implementation manages DM
ierr = DMSetType(dM,"DMMOFEM");CHKERRG(ierr);
// mesh is portioned, each process keeps only part of problem
ierr = DMMoFEMSetIsPartitioned(dM,PETSC_TRUE);
// creates problem in DM, "MIX" is name of problem
ierr = DMMoFEMCreateMoFEM(dM,&mField,"MIX",ref_level); CHKERRG(ierr);
// discretized problem creates square matrix (that makes some optimizations)
ierr = DMMoFEMSetIsPartitioned(dM,PETSC_TRUE); CHKERRG(ierr);
// set DM options from command line
ierr = DMSetFromOptions(dM); CHKERRG(ierr);
// add finite elements
ierr = DMMoFEMAddElement(dM,"MIX"); CHKERRG(ierr);
ierr = DMMoFEMAddElement(dM,"MIX_BCFLUX"); CHKERRG(ierr);
ierr = DMMoFEMAddElement(dM,"MIX_BCVALUE"); CHKERRG(ierr);
// constructor data structures
ierr = DMSetUp(dM); CHKERRG(ierr);

Note that problem in DM is created by setting finite elements on which it operates.

Creating finite element instances

Definition of the problem enables to construct a matrix since matrix connectivity is derived from problem definition. To fill the matrix with the numbers you need implement finite element. In fact, you can have more than one implementation of the same problem.

Once you have problem defined, the developer is released from the need managing complexities related to management of DOFs, iteration over elements, matrix assembly, etc. It can focus on the integration of local entity matrices, as is shown in Semi-discrete linearized system of equations.

MoFEM provides technology with simplifies implementation of the element, enabling the developer to break a complex problem into smaller tasks, which are easier to implement, debug and reuse in different contexts, that has been discussed in several tutorials already. Assembly of residual vectors and tangent matrices is made with the use so called users data operators (UDO), for example, the sequence of users data operators of finite element instance to assemble tangent matrix is shown in figure 3.

unsaturated_flow_fig1.png
Figure 3. Sequence of operators on instance of MIX element

We have already discuss last four operators in sequence for feVolLhs finite element instance, those are used to assemble tangent matrix blocks, i.e. \(\mathbf{C}\), \(\mathbf{G}\), \(\mathbf{H}\) and \(\mathbf{A}\). First, three operators are used to evaluate fields at integration points, i.e. capillary head, flux and rate of capillary head. Note first three operators has been are reused from other tutorials.

In this problem, we using two generic implementations of finite element, MoFEM::VolumeElementForcesAndSourcesCore and MoFEM::FaceElementForcesAndSourcesCore, for integrating fields on volumes and faces. Finite element instances are iterated over a particular element in the problem, and for each entity on finite element entity, the sequence of UDO is called.

We use smart pointers to create finite element instances, which gives code flexibility and safe memory management. Finite element instances are created in method MixTransport::UnsaturatedFlowElement::setFiniteElements, as follows

feFaceBc = boost::shared_ptr<ForcesAndSourcesCore>(new FaceElementForcesAndSourcesCore(mField));
feFaceRhs = boost::shared_ptr<ForcesAndSourcesCore>(new FaceElementForcesAndSourcesCore(mField));
feVolInitialPc = boost::shared_ptr<ForcesAndSourcesCore>(new VolumeElementForcesAndSourcesCore(mField));
feVolLhs = boost::shared_ptr<ForcesAndSourcesCore>(new VolumeElementForcesAndSourcesCore(mField));
feVolRhs = boost::shared_ptr<ForcesAndSourcesCore>(new VolumeElementForcesAndSourcesCore(mField));

Next, we set integration rule, as has been discussed in tutorial Solving the Poisson equation, as follows

// set integration rule to elements
feFaceBc->getRuleHook = face_rule;
feFaceRhs->getRuleHook = face_rule;
feVolInitialPc->getRuleHook = vol_rule;
feVolLhs->getRuleHook = vol_rule;
feVolRhs->getRuleHook = vol_rule;

In addition, we add hooks to functions which are run at the beginning and the end of the loop over a particular element. Those functions are used to perProcess, and postProcess matrix assemble. Those pre-processing and post-processing hooks are added as follows

feVolRhs->preProcessHook = preProcessVol(*this,feVolRhs);
feVolLhs->preProcessHook = preProcessVol(*this,feVolLhs);
feVolRhs->postProcessHook = postProcessVol(*this,feVolRhs);
feVolLhs->postProcessHook = postProcessVol(*this,feVolLhs);

We use classes MixTransport::UnsaturatedFlowElement::preProcessVol and MixTransport::UnsaturatedFlowElement::postProcessVol to enforce essential boundary conditions. The method of enforcing essential boundary conditions is discussed in Mix formulation and integration on skeleton (h-adaptivity). In additional in function MixTransport::UnsaturatedFlowElement::postProcessVol::operator()(), we scatter vector of rate of capillary pressure on field data. Since rate of capillary head is not primary unknown, but linear function of solution at current and previous steps, is evaluated by PETSc for particular integration scheme.

ierr = fePtr->mField.getInterface<VecManager>()->setOtherLocalGhostVector(
fePtr->problemPtr,"VALUES",string("VALUES")+"_t",
ROW,fePtr->ts_u_t,INSERT_VALUES,SCATTER_REVERSE
); CHKERRG(ierr);

Time solver provides vector of DOFs, \(\left(\frac{\textrm{d}\mathbf{p}}{\textrm{d}t}\right)^{i+1}_{n+1}\) at integration \(i+1\) and time step \(n+1\), and above code scatter vector values on mesh, such that values of capillary head at integration points are early evaluated, i.e. \(\frac{\textrm{d}h}{\textrm{d}t}|^{i+1}_{n+1}\). Capillary head DOFs are sub-vector of fePtr->ts_u_t.

Adding operators to elements

Users data operators are added to finite element calculating tangent matrix as follows

// set tangent matrix finite element operators
feVolLhs->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues(string("VALUES")+"_t",headRateAtGaussPts,MBTET)
);
feVolLhs->getOpPtrVector().push_back(new OpValuesAtGaussPts(*this,"VALUES"));
feVolLhs->getOpPtrVector().push_back(new OpFluxDivergenceAtGaussPts(*this,"FLUXES"));
feVolLhs->getOpPtrVector().push_back(new OpTauDotSigma_HdivHdiv(*this,"FLUXES"));
feVolLhs->getOpPtrVector().push_back(new OpVU_L2L2(*this,"VALUES"));
feVolLhs->getOpPtrVector().push_back(new OpVDivSigma_L2Hdiv(*this,"VALUES","FLUXES"));
feVolLhs->getOpPtrVector().push_back(new OpDivTauU_HdivL2(*this,"FLUXES","VALUES"));

and similarly for other elements, evaluating boundary conditions and vector of residuals. For more details we refer to tutorial Hello world.

In addition finite element instances are created to post-process results and to run time stepping monitor, see code in MixTransport::UnsaturatedFlowElement::setFiniteElements.

Adding finite element instances to DM

Once we have created finite element instances, in MixTransport::UnsaturatedFlowElement::setFiniteElements, we add finite instances to DM manager, such that system of equations could be assembled. This is simply done by

ierr = DMMoFEMTSSetIFunction(dM,"MIX_BCVALUE",feFaceRhs,null,null); CHKERRG(ierr);
ierr = DMMoFEMTSSetIFunction(dM,"MIX",feVolRhs,null,null); CHKERRG(ierr);
ierr = DMMoFEMTSSetIJacobian(dM,"MIX",feVolLhs,null,null); CHKERRG(ierr);

Those functions are specific to MoFEM implementation, there were implemented based on PETSc documentation, link, see Chapter 6.

Solving problem using TS

Finally we can solve problem, by creating TS solver and link it to DM, as follows

// Create time solver
TS ts;
ierr = TSCreate(PETSC_COMM_WORLD,&ts); CHKERRG(ierr);
// Use backward Euler method
ierr = TSSetType(ts,TSBEULER); CHKERRG(ierr);
// Not test functions to evalye matrix and rigth hand side, TS ask DM
// for doing that
ierr = TSSetIFunction(ts,F,PETSC_NULL,PETSC_NULL); CHKERRG(ierr);
ierr = TSSetIJacobian(ts,Aij,Aij,PETSC_NULL,PETSC_NULL); CHKERRG(ierr);
// Set final time
double ftime = 1;
ierr = TSSetDuration(ts,PETSC_DEFAULT,ftime); CHKERRG(ierr);
// Setup solver from commabd line
ierr = TSSetFromOptions(ts); CHKERRG(ierr);
// Set DM to TS
ierr = TSSetDM(ts,dM); CHKERRG(ierr);
#if PETSC_VERSION_GE(3,7,0)
ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER); CHKERRG(ierr);
#endif
// Set-up monitor
TsCtx *ts_ctx;
DMMoFEMGetTsCtx(dM,&ts_ctx);
ierr = TSMonitorSet(ts,f_TSMonitorSet,ts_ctx,PETSC_NULL); CHKERRG(ierr);
ierr = TSSolve(ts,D); CHKERRG(ierr);

Note that TS calls finite elements instance by DM, similarly how SNES (Newton solver) calls finite elements in A nonlinear Poisson equation, or KSP (Krylov solver) calls finite elements in Solving the Poisson equation. Above is implemented in MixTransport::UnsaturatedFlowElement::solveProblem, where prior to solution initial conditions are calculated, and statistics at the end printed to the terminal.

Application code

Application code is unsaturated_transport.cpp, which calls sequence of functions

// Initialise
PetscInitialize(&argc,&argv,(char *)0,help);
// Register materials
ierr = RegisterMaterials()(); CHKERRG(ierr);
// Register MOFEM DM implementation
ierr = DMRegister_MoFEM("DMMOFEM"); CHKERRG(ierr); // register MoFEM DM in PETSc
try {
// Create moab instance
moab::Core mb_instance;
moab::Interface& moab = mb_instance;
// Read mesh
...
// Create mofem interface
MoFEM::Core core(moab);
MoFEM::Interface& m_field = core;
// Create unsaturated flow structures
UnsaturatedFlowElement uf(m_field);
\\ Read data from meshsets
...
ierr = uf.addFields("VALUES","FLUXES",order); CHKERRG(ierr);
ierr = uf.addFiniteElements("FLUXES","VALUES"); CHKERRG(ierr);
zero_flux_ents,MBTRI,"MIX_BCFLUX"
); CHKERRG(ierr);
ierr = uf.buildProblem(); CHKERRG(ierr);
ierr = uf.createMatrices(); CHKERRG(ierr);
ierr = uf.setFiniteElements(); CHKERRG(ierr);
ierr = uf.calculateEssentialBc(); CHKERRG(ierr);
ierr = uf.calculateInitialPc(); CHKERRG(ierr);
ierr = uf.solveProblem(); CHKERRG(ierr);
ierr = uf.destroyMatrices(); CHKERRG(ierr);
} catch (MoFEMException const &e) {
SETERRQ(PETSC_COMM_SELF,e.errorCode,e.errorMessage);
}
ierr = PetscFinalize(); CHKERRG(ierr);

To do

Todo:

Add example with hp-adaptivity

Implement block-solver and multigrid

Apply Kirchhoff transform

Add boundary condition for evaporation

Add air phase and mechanical filed

Add diffusion and other physical phenomena

Make example with flux history or head history. That is implemented but no example is added

Add diffusion and other physical phenomenas