37. Parameter Dependent Problems#
Many problems from mechanics are of the form
with the bilinear-form depending on a small parameter \(\varepsilon\):
where \(a(.,.)\) and \(c(.,.)\) are non-negative bilinear-forms on Hilbert spaces \(V\) and \(Q\), and
is a continuous, linear operator.
We introduce a new variable
Using this new variable in the first row, and writing down the definition of \(p\) in weak form we get the system
We rename the off-diagonal forms as
This structure is called a mixed formulation with penalty term.
In the mixed formulation we can set the parameter \(\varepsilon = 0\).
We are interested in the convergence
37.1. Example: Dirichlet boundary condition by penalty#
Here,
37.2. Example: Nearly incompressible materials#
with \(\varepsilon(u) = \tfrac{1}{2} (\nabla u + (\nabla u)^T)\), and \(\lambda \gg \mu\). Introduce \(p := \lambda \operatorname{div} u\)
from ngsolve import *
from ngsolve.webgui import Draw
mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
fes = VectorH1(mesh,order=2, dirichlet="bottom|top")
u,v = fes.TnT()
mu,lam = 1, 1e5
a = BilinearForm(2*mu*InnerProduct (Sym(Grad(u)), Sym(Grad(v))) * dx \
+ lam*div(u)*div(v)*dx).Assemble()
f = LinearForm(v[0]*ds("left")).Assemble()
gfu = GridFunction(fes)
gfu.vec.data = a.mat.Inverse(fes.FreeDofs())*f.vec
The deformation looks reasonable, however the pressure looks strange:
Draw (gfu, deformation=True);
Draw (lam*div(gfu), mesh);
Now, we replace the term
by introducing an \(L_2\)-projection into elements-wise constants:
fesP = L2(mesh, order=0)
a = BilinearForm(2*mu*InnerProduct (Sym(Grad(u)), Sym(Grad(v))) * dx \
+ lam*Interpolate(div(u),fesP)*div(v)*dx).Assemble()
f = LinearForm(v[0]*ds("left")).Assemble()
gfu = GridFunction(fes)
gfu.vec.data = a.mat.Inverse(fes.FreeDofs())*f.vec
Now, also the pressure looks good (as good as it can be in \(P^0\)):
Draw (gfu, deformation=True);
Draw (lam*Interpolate(div(gfu), fesP), mesh);
How to analyze that ?
Instead of solving the parameter dependent problem for \(u\), we define the pressure \(p\):
Inserting \(p\), and using the defining equation for \(p\) we obtain the saddle-point problem:
Now, the large parameter \(\lambda\) got a small paramter \(\lambda^{-1}\). The limit problem is a Stokes equation.
Discretize it with a stable method for Stokes leads to the linear system
The meaning of the second discrete equation is
If we use discontinuous finite elements for the pressure, the matrix \(C\) is block-diagonal, and thus cheaply invertible. This allows to eliminate \(p\) as
and insert it into the first equation to obtain
We obtain exactly the same system by inserting the projector into the bilinear-form.