-
Notifications
You must be signed in to change notification settings - Fork 48
NSCoupled
NSCoupled is a steady state incompressible Navier-Stokes solver using a Newton iterative method.
The coupled solver uses a mixed function space for velocity and pressure. There are 3 choices of elements
- TaylorHood - Continuous Lagrange elements for both velocity and pressure. Order of velocity space one higher than pressure.
- CR - Crouzeix-Raviart for velocity and 0 order discontinuous Lagrange for pressure.
- MINI - Linear Lagrange for both velocity and pressure. Velocity space enriched with bubble function.
The solver is run using, e.g.,
oasis NSCoupled problem=DrivenCavity element=CR
The DrivenCavity problem requires some special attention, because this problem may be studied both under transient and steady state conditions. For this reason the mesh
and SubDomain
s are defined in a "Base" module at the highest level of the problems
folder. The file problems/DrivenCavity.py
thus contains
from dolfin import UnitSquareMesh
# Create a mesh
def mesh(Nx=50, Ny=50, **params):
m = UnitSquareMesh(Nx, Ny)
return m
noslip = "std::abs(x[0]*x[1]*(1-x[0]))<1e-8"
top = "std::abs(x[1]-1) < 1e-8"
bottom = "std::abs(x[1]) < 1e-8"
This information may be used both by the transient NSfracStep
solver and the current NSCoupled
solver. To complete the formulation of the problem for the coupled case, the remaining (more specific part) is placed in module problems/NSCoupled/DrivenCavity.py
from ..NSCoupled import *
from ..DrivenCavity import *
# Override some problem specific parameters
NS_parameters.update(
nu = 0.001,
max_iter = 100)
# Specify boundary conditions
def create_bcs(VQ, **NS_namespace):
bc0 = DirichletBC(VQ.sub(0), (0, 0), noslip)
bc1 = DirichletBC(VQ.sub(0), (1, 0), top)
return dict(up = [bc0, bc1])
def theend_hook(u_, p_, mesh, **NS_namespace):
plot(u_, title='Velocity')
plot(p_, title='Pressure')
Note that we first import mesh and subdomains from the base module, before we specify parameters and boundary conditions. The implementation of boundary conditions for the fractional step solver is different since the function spaces are different. Otherwise, theend_hook
may look exactly the same, because u_
and p_
are velocity vector and pressure in both solvers.