Navier-Stokes equation

We consider the incompressible Navier-Stokes equation in a domain Ω=(0,Lx)×(0,Ly)\Omega=(0,L_x)\times(0,L_y),

ϱ(tu+(u)u)(2νD(u))p=f, in Ωu=0 \varrho(\partial_t\mathbf{u} + (\mathbf{u}\cdot\nabla)\mathbf{u}) - \nabla\cdot\big(2\nu\mathbf{D}(\mathbf{u})\big) - \nabla p = \mathbf{f},\qquad\text{ in }\Omega \\ \nabla\cdot\mathbf{u} = 0

with velocity u\mathbf{u}, pressure pp, density ϱ\varrho, dynamic viscosity ν\nu, symmetric rate-of-strain tensor D(u)=12(u+u)\mathbf{D}(\mathbf{u})=\frac{1}{2}(\nabla\mathbf{u}+\nabla\mathbf{u}^\top), and external volume force f\mathbf{f}, w.r.t. to boundary conditions u=g\mathbf{u}=\mathbf{g} on Ω\partial\Omega. For an overview, see also Wikipedia.

Discretization

Weak formulation

In a weak variational formulation, we try to find u(t,)Vg:={vH1(Ω)d:trΩv=g}\mathbf{u}(t,\cdot)\in \mathbf{V}_\mathbf{g}:=\{\mathbf{v}\in H^1(\Omega)^d\,:\, \operatorname{tr}_{\partial\Omega}\mathbf{v} = \mathbf{g}\} and p(t,)Q:=L02(Ω)p(t,\cdot)\in Q:=L_0^2(\Omega), such that

Ωϱ(tu+(u)u)v+νu:v+pvdx=Ωfvdx,vV0,Ωqudx=0,qQ. \int_\Omega \varrho(\partial_t\mathbf{u} + (\mathbf{u}\cdot\nabla)\mathbf{u})\cdot\mathbf{v} + \nu\nabla\mathbf{u}:\nabla\mathbf{v} + p\,\nabla\cdot\mathbf{v}\,\text{d}\mathbf{x} = \int_\Omega \mathbf{f}\cdot\mathbf{v}\,\text{d}\mathbf{x},\qquad\forall \mathbf{v}\in \mathbf{V}_0, \\ \int_\Omega q\,\nabla\cdot\mathbf{u}\,\text{d}\mathbf{x} = 0,\qquad \forall q\in Q.

for t[0,T]t\in[0,T], given u(0,)u0\mathbf{u}(0,\cdot)\equiv\mathbf{u}^0.

The pair (u(t,),p(t,))\big(\mathbf{u}(t,\cdot), p(t,\cdot)\big) is thus in the product space Vg×Q\mathbf{V}_\mathbf{g}\times Q where Vg\mathbf{V}_\mathbf{g} is also a product of H1H^1 spaces, i.e., Vg=Vg0×Vg1×Vg2[Vgi]d\mathbf{V}_\mathbf{g}=V_{g_0}\times V_{g_1}\times V_{g_2}\simeq [V_{g_i}]^d. Later, this space will be approximated with a Taylor-Hood space of piecewise quadratic and linear functions.

Time-discretization

A linearized time-discretization with a semi-implicit backward Euler scheme then reads:

For k=0,1,2,,Nk=0,1,2,\ldots, N find uk+1Vg,pQ\mathbf{u}^{k+1}\in\mathbf{V}_\mathbf{g},\, p\in Q, s.t.

Ωϱ(1τ(uk+1uk)+(uk)uk+1)v+νuk+1:v+pvdx=Ωf(tk+1)vdx,vV0,Ωquk+1dx=0,qQ. \int_\Omega \varrho\Big(\frac{1}{\tau}(\mathbf{u}^{k+1} - \mathbf{u}^k) + (\mathbf{u}^k\cdot\nabla)\mathbf{u}^{k+1}\Big)\cdot\mathbf{v} + \nu\nabla\mathbf{u}^{k+1}:\nabla\mathbf{v} + p\,\nabla\cdot\mathbf{v}\,\text{d}\mathbf{x} \\ = \int_\Omega \mathbf{f}(t^{k+1})\cdot\mathbf{v}\,\text{d}\mathbf{x},\qquad\forall \mathbf{v}\in \mathbf{V}_0, \\ \int_\Omega q\,\nabla\cdot\mathbf{u}^{k+1}\,\text{d}\mathbf{x} = 0,\qquad \forall q\in Q.

with 0=t0<t1<<tN=T0=t^0<t^1<\ldots<t^N=T, timestep τ=tk+1tk\tau=t^{k+1} - t^k, and given initial velocity u0\mathbf{u}^0.

Finite-Element formulation

Let Ω\Omega be partitioned into quasi-uniform elements by a conforming triangulation Th\mathcal{T}_h. The functions u\mathbf{u} and pp are approximated by functions from the discrete functionspace Th:=[V2,gi]d×V1\mathbb{T}_h := [\mathbb{V}_{2,\,g_i}]^d\times \mathbb{V}_1 with Vp:={vC0(Ω):vTPp(T),TTh}\mathbb{V}_{p} := \{v\in C^0(\Omega)\,:\, v|_T\in\mathbb{P}_p(T),\, \forall T\in\mathcal{T}_h\} and Vp,g:={vVp:trΩv=g}\mathbb{V}_{p,\,g} := \{v\in \mathbb{V}_p\,:\,\operatorname{tr}_{\partial\Omega} v=g\}, where Pp\mathbb{P}_p is the space of polynomials of order at most pp. The space Th\mathbb{T}_h will be denoted as Taylor-Hood space.

Let Grid grid be the Dune grid type representing the triangulation Th\mathcal{T}_h of the domain Ω\Omega, then the basis of the Taylor-Hood space can be written as a composition of Lagrange bases:

using namespace Dune::Functions::BasisFactory;
auto basis = composite(
      power<Grid::dimensionworld>(
        lagrange<2>(), flatInterleaved()),
      lagrange<1>(), flatLexicographic());
At this point, any Dirichlet boundary conditions are ignored in the definition of the bases. We build the Taylor-Hood basis as a composition of a product of d×V2d\times\mathbb{V}_2 bases for the velocity components and a V1\mathbb{V}_1 basis for the pressure. We have to use composite(...) to combine the different types for velocity and pressure space, while we can use power<d>(...) for the d velocity component spaces of the same type.

For efficiency, we use flat indexing strategies flatInterleaved and flatLexicographic to build a global continuous numbering of the basis functions.

See the tutorial Grids and Discrete Functions for more details.

Implementation

We split the implementation into three parts, 1. the time derivative, 2. the stokes operator and 3. any external forces. For addressing the different components of the Taylor-Hood basis, we introduce the treepaths _v = Dune::Indices::_0 and _p = Dune::Indices::_1 for velocity and pressure, respectively.

Problem framework

We have an in-stationary problem consisting of a sequence of stationary equations for each timestep, thus we combine a ProblemInstat with a ProblemStat .

using namespace AMDiS;
ProblemStat prob{"stokes", grid, basis};
prob.initialize(INIT_ALL);

ProblemInstat probInstat{"stokes", prob};
probInstat.initialize(INIT_UH_OLD);
For convenience, we use class template argument deduction here.

Time derivative

For simplicity of this example, we use a backward Euler time stepping and a linearization of the nonlinear advection term:

// define a constant fluid density
double density = 1.0;
Parameters::get("stokes->density", density);

// a reference to 1/tau
auto invTau = probInstat.invTau();

// <1/tau * u, v>
auto opTime = makeOperator(tag::testvec_trialvec{},
  density * invTau);
prob.addMatrixOperator(opTime, _v, _v);

// <1/tau * u^old, v>
auto opTimeOld = makeOperator(tag::testvec{},
  density * invTau * probInstat.oldSolution(_v));
prob.addVectorOperator(opTimeOld, _v);

for (int i = 0; i < Grid::dimensionworld; ++i) {
  // <(u^old * nabla)u_i, v_i>
  auto opNonlin = makeOperator(tag::test_gradtrial{},
    density * prob.solution(_v));
  prob.addMatrixOperator(opNonlin, makeTreePath(_v,i), makeTreePath(_v,i));
}

Stokes operator

The stokes part of the Navier-Stokes equation is a common pattern that emerges in several fluid equations.

Ωνu:v+pv+qudx \int_\Omega \nu\nabla\mathbf{u}:\nabla\mathbf{v} + p\,\nabla\cdot\mathbf{v} + q\,\nabla\cdot\mathbf{u}\,\text{d}\mathbf{x}

for u,vV\mathbf{u},\mathbf{v}\in V and p,qQp,q\in Q.

// define a constant fluid viscosity
double viscosity = 1.0;
Parameters::get("stokes->viscosity", viscosity);

for (int i = 0; i < Grid::dimensionworld; ++i) {
  // <viscosity*grad(u_i), grad(v_i)>
  auto opL = makeOperator(tag::gradtest_gradtrial{}, viscosity);
  prob.addMatrixOperator(opL, makeTreePath(_v,i), makeTreePath(_v,i));
}

// <d_i(v_i), p>
auto opP = makeOperator(tag::divtestvec_trial{}, 1.0);
prob.addMatrixOperator(opP, _v, _p);

// <q, d_i(u_i)>
auto opDiv = makeOperator(tag::test_divtrialvec{}, 1.0);
prob.addMatrixOperator(opDiv, _p, _v);
where we have used the identity

u:v=iuivi \nabla\mathbf{u}:\nabla\mathbf{v} = \sum_i \nabla u_i\cdot\nabla v_i

External volume force

The force may be implemented as a vector-valued function of the global coordinates or any other expression that leads to a local force vector at the quadrature points, e.g.

using WorldVector = typename Grid::template Codim<0>::Geometry::GlobalCoordinate;
auto opForce = makeOperator(tag::testvec{},
  [](WorldVector const& x) { return WorldVector{1.0, 0.0}; });
prob.addVectorOperator(opForce, _v);

Numerical example

We consider a lid-driven cavity problem in a square domain Ω\Omega with Lx=Ly=1L_x=L_y=1 no external force, and with Dirichlet boundary conditions

g(x)=(0v1x1(1x1)(1x0)) \mathbf{g}(\mathbf{x}) = \begin{pmatrix}0 \\ v_1 x_1(1 - x_1)(1 - x_0)\end{pmatrix}

Thus, the boundary value is zero on top, right, and bottom boundary and parabolic on the left boundary. The boundary condition is implemented, by first setting boundary ids for parts of the grid boundary and second, setting the Dirichlet value:

double vel = 1.0;
Parameters::get("stokes->boundary velocity", v1);

// define boundary values
auto g = [vel](WorldVector const& x)
{
  return WorldVector{0.0, v1*x[1]*(1.0 - x[1])*(1.0 - x[0])};
};

// set boundary conditions for velocity
prob.boundaryManager()->setBoxBoundary({1,1,1,1});
prob.addDirichletBC(1, _v, _v, g);

Simulation

The time-stepping process with a stationary problem in each iteration can be started using an AdaptInstationary manager class:

// set initial conditions
prob.solution(_v).interpolate(g);

// start simulation
AdaptInfo adaptInfo("adapt");
AdaptInstationary adapt("adapt", prob, adaptInfo, probInstat, adaptInfo);
adapt.adapt();

Complete Example

The full source code of the Navier-Stokes example can be found in the repository at examples/navier_stokes.cc.

Compile with

cmake --build build --target navier_stokes.2d

and run with

./build-cmake/examples/navier_stokes.2d examples/init/navier_stokes.dat.2d

The flow field with density=1, viscosity=1 and boundary velocity=10 will look like