Netgen/NGSolve

i-Tutorials

  • Interactive NGSolve Tutorial
  • Whetting the Appetite
  • 1. Getting started
  • 2. Advanced Topics
  • 3. Time-dependent and non-linear problems
  • 4. Geometric modeling and mesh generation
  • 5. MPI - Parallelization and CUDA Support
  • 6. Various Topics
  • 7. Shape- and Topology Optimization
  • 8. Unfitted Finite Elements
  • 9. Extending by C++ programming
  • 10. NGSolve and ...
  • 11. Boundary element and Fast Multipole methods
  • Appendix

C++ Tutorials

  • NGSolve Basics
  • Advanced NGSolve Examples
Netgen/NGSolve
  • Docs »
  • Navier Stokes Equation
  • Homepage

Navier Stokes Equation¶

We solve the time-dependent incompressible Navier Stokes Equation. For that

  • we use the P3/P2 Taylor-Hood mixed finite element pairing

  • and perform operator splitting time-integration with the non-linear term explicit, but time-dependent Stokes fully implicit.

The example is from the Schäfer-Turek benchmark <http://www.mathematik.tu-dortmund.de/lsiii/cms/papers/SchaeferTurek1996.pdf> a two-dimensional cylinder, at Reynolds number 100

Download navierstokes.py

from ngsolve import *

# viscosity
nu = 0.001

# timestepping parameters
tau = 0.001
tend = 10

from netgen.geom2d import SplineGeometry
geo = SplineGeometry()
geo.AddRectangle( (0, 0), (2, 0.41), bcs = ("wall", "outlet", "wall", "inlet"))
geo.AddCircle ( (0.2, 0.2), r=0.05, leftdomain=0, rightdomain=1, bc="cyl", maxh=0.02)
mesh = Mesh( geo.GenerateMesh(maxh=0.07))

mesh.Curve(3)

V = VectorH1(mesh,order=3, dirichlet="wall|cyl|inlet")
Q = H1(mesh,order=2)

X = V*Q

u,p = X.TrialFunction()
v,q = X.TestFunction()

stokes = nu*InnerProduct(grad(u), grad(v))+div(u)*q+div(v)*p - 1e-10*p*q
a = BilinearForm(X, symmetric=True)
a += stokes*dx
a.Assemble()

# nothing here ...
f = LinearForm(X)   
f.Assemble()

# gridfunction for the solution
gfu = GridFunction(X)

# parabolic inflow at inlet:
uin = CoefficientFunction( (1.5*4*y*(0.41-y)/(0.41*0.41), 0) )
gfu.components[0].Set(uin, definedon=mesh.Boundaries("inlet"))

# solve Stokes problem for initial conditions:
inv_stokes = a.mat.Inverse(X.FreeDofs())

res = f.vec.CreateVector()
res.data = f.vec - a.mat*gfu.vec
gfu.vec.data += inv_stokes * res


# matrix for implicit Euler 
mstar = BilinearForm(X, symmetric=True)
mstar += (u*v + tau*stokes)*dx
mstar.Assemble()
inv = mstar.mat.Inverse(X.FreeDofs(), inverse="sparsecholesky")

# the non-linear term 
conv = BilinearForm(X, nonassemble = True)
conv += (grad(u) * u) * v * dx

# for visualization
Draw (Norm(gfu.components[0]), mesh, "velocity", sd=3)

# implicit Euler/explicit Euler splitting method:
t = 0
with TaskManager():
    while t < tend:
        print ("t=", t, end="\r")

        conv.Apply (gfu.vec, res)
        res.data += a.mat*gfu.vec
        gfu.vec.data -= tau * inv * res    

        t = t + tau
        Redraw()

The absolute value of velocity:

../_images/navierstokes_solution.jpg

© Copyright 2017, Netgen/NGSolve team.

Built with Sphinx using a theme provided by Read the Docs.
NGSolve v: 6.2.2604-97-gab8e0804e