Chapter 9
Mathematical Models
Summary
This chapter goes more in depth into a number of problems that FreeFem++ can solve. It is a complement to
chapter 3 which was only an introduction. Users are invited to contribute to make this data base of problem
solutions grow.
9.1 Static Problems
9.1.1 Soap Film
Our starting point here will be the mathematical model to find the shape of soap film which is glued to the ring on
the xy-plane
We assume the shape of the film is described as the graph (x,y,u(x,y)) of the vertical displacement
u(x,y)(x2 + y2 < 1) under a vertical pressure p in terms of force per unit area and an initial tension μ in terms of
force per unit length. Consider “small plane” ABCD, A:(x,y,u(x,y)), B:(x,y,u(x + δx,y)), C:(x,y,u(x + δx,y + δy))
and D:(x,y,u(x,y + δy)). Let us denote by
(x,y) = (nx(x,y),ny(x,y),nz(x,y)) the normal vector of the surface
z = u(x,y). We see that the vertical force due to the tension μ acting along the edge AD is -μnx(x,y)δy and the the
vertical force acting along the edge AD is
Similarly, for the edges AB and DC we have
The force in the vertical direction on the surface ABCD due to the tension μ is given by
Assuming small displacements, we have
Letting δx → dx,δy → dy, we have the equilibrium of the vertical displacement of soap film on ABCD by
p
Using the Laplace operator Δ = ∂2∕∂x2 + ∂2∕∂y2, we can find the virtual displacement write the following
 | (9.1) |
where f = p∕μ, Ω = {(x,y); x2 + y2 < 1}. Poisson’s equation (2.1) appear also in electrostatics taking the form of
f = ρ∕ϵ where ρ is the charge density, ϵ the dielectric constant and u is named as electrostatic potential. The soap
film is glued to the ring ∂Ω = C, then we have the boundary condition
 | (9.2) |
If the force is gravity, for simplify, we assume that f = -1.
Example 9.1 (a_tutorial.edp)
1 : border a(t=0,2⋆pi){ x = cos(t); y = sin(t);label=1;};
2 :
3 : mesh disk = buildmesh(a(50));
4 : plot(disk);
5 : fespace femp1(disk,P1);
6 : femp1 u,v;
7 : func f = -1;
8 : problem laplace(u,v) =
9 : int2d(disk)( dx(u)⋆dx(v) + dy(u)⋆dy(v) ) // bilinear
form
10 : - int2d(disk)( f⋆v ) // linear
form
11 : + on(1,u=0) ; // boundary
condition
12 : func ue = (x^2+y^2-1)/4; // ue: exact
solution
13 : laplace;
14 : femp1 err = u - ue;
15 :
16 : plot (u,ps="aTutorial.eps",value=true,wait=true);
17 : plot(err,value=true,wait=true);
18 :
19 : cout << "error L2=" << sqrt(int2d(disk)( err^2) )<< endl;
20 : cout << "error H10=" << sqrt( int2d(disk)((dx(u)-x/2)^2)
21 : + int2d(disk)((dy(u)-y/2)^2))<< endl;
22 :
23 : disk = adaptmesh(disk,u,err=0.01);
24 : plot(disk,wait=1);
25 :
26 : laplace;
27 :
28 : plot (u,value=true,wait=true);
29 : err = u - ue; // become FE-function on adapted
mesh
30 : plot(err,value=true,wait=true);
31 : cout << "error L2=" << sqrt(int2d(disk)( err^2) )<< endl;
32 : cout << "error H10=" << sqrt(int2d(disk)((dx(u)-x/2)^2)
33 : + int2d(disk)((dy(u)-y/2)^2))<< endl;
In 19th line, the L2-error estimation between the exact solution ue,
and
from 20th line to 21th line, the H1-error seminorm estimation
are
done on the initial mesh. The results are ∥uh -ue∥0,Ω = 0.000384045,|uh -ue|1,Ω = 0.0375506.
After the adaptation, we hava ∥uh -ue∥0,Ω = 0.000109043,|uh -ue|1,Ω = 0.0188411. So the numerical solution is
improved by adaptation of mesh.
9.1.2 Electrostatics
We assume that there is no current and a time independent charge distribution. Then the electric field
satisfy
where ρ is the charge density and ϵ is called the permittivity of free space. From the second equation in (9.3), we can
introduce the electrostatic potential such that
= -∇ϕ. Then we have Poisson equation -Δϕ = f, f = -ρ∕ϵ. We
now obtain the equipotential line which is the level curve of ϕ, when there are no charges except conductors
{Ci}1,
,K. Let us assume K conductors C1,
,CK within an enclosure C0. Each one is held at an electrostatic
potential φi. We assume that the enclosure C0 is held at potential 0. In order to know φ(x) at any point x of the
domain Ω, we must solve
 | (9.4) |
where Ω is the interior of C0 minus the conductors Ci, and Γ is the boundary of Ω, that is ∑
i=0NCi.
Here g is any function of x equal to φi on Ci and to 0 on C0. The boundary equation is a reduced form
for:
 | (9.5) |
Example 9.2 First we give the geometrical informations; C0 = {(x,y); x2 + y2 = 52}, C1 = {(x,y) :
(x -2)2 +
y2 = 1}, C2 = {(x,y) :
(x + 2)2 +
y2 = 1}. Let Ω be the disk enclosed by C0 with the
elliptical holes enclosed by C1 and C2. Note that C0 is described counterclockwise, whereas the elliptical
holes are described clockwise, because the boundary must be oriented so that the computational domain is
to its left.
// a circle with center at (0 ,0) and radius
5
border C0(t=0,2⋆pi) { x = 5 ⋆ cos(t); y = 5 ⋆ sin(t); }
border C1(t=0,2⋆pi) { x = 2+0.3 ⋆ cos(t); y = 3⋆sin(t); }
border C2(t=0,2⋆pi) { x = -2+0.3 ⋆ cos(t); y = 3⋆sin(t); }
mesh Th = buildmesh(C0(60)+C1(-50)+C2(-50));
plot(Th,ps="electroMesh"); // figure
9.3
fespace Vh(Th,P1); // P1
FE-space
Vh uh,vh; // unknown and test
function.
problem Electro(uh,vh) = // definition of the
problem
int2d(Th)( dx(uh)⋆dx(vh) + dy(uh)⋆dy(vh) ) //
bilinear
+ on(C0,uh=0) // boundary condition on
C0
+ on(C1,uh=1) // +1 volt on
C1
+ on(C2,uh=-1) ; // -1 volt on
C2
Electro; // solve the problem, see figure 9.4 for the
solution
plot(uh,ps="electro.eps",wait=true); // figure
9.4
9.1.3 Aerodynamics
Let us consider a wing profile S in a uniform flow. Infinity will be represented by a large circle Γ∞. As previously,
we must solve
 | (9.6) |
where Ω is the area occupied by the fluid, u∞ is the air speed at infinity, c is a constant to be determined so that ∂nφ
is continuous at the trailing edge P of S (so-called Kutta-Joukowski condition). Lift is proportional to c. To find c
we use a superposition method. As all equations in (9.6) are linear, the solution φc is a linear function of
c
 | (9.7) |
where φ0 is a solution of (9.6) with c = 0 and φ1 is a solution with c = 1 and zero speed at infinity. With these two
fields computed, we shall determine c by requiring the continuity of ∂φ∕∂n at the trailing edge. An equation for the
upper surface of a NACA0012 (this is a classical wing profile in aerodynamics; the rear of the wing is called the
trailing edge) is:
 | (9.8) |
Taking an incidence angle α such that tan α = 0.1, we must solve
 | (9.9) |
where Γ2 is the wing profile and Γ1 is an approximation of infinity. One finds c by solving:
The solution φ = φ0 + cφ1 allows us to find c by writing that ∂nφ has no jump at the trailing edge P = (1,0). We have
∂nφ-(φ(P+) -φ(P))∕δ where P+ is the point just above P in the direction normal to the profile at a distance δ. Thus
the jump of ∂nφ is (φ0|P+ + c(φ1|P+ -1)) + (φ0|P- + c(φ1|P- -1)) divided by δ because the normal changes sign
between the lower and upper surfaces. Thus
 | (9.12) |
which can be programmed as:
 | (9.13) |
Example 9.3
// Computation of the potential flow around a NACA0012
airfoil.
// The method of decomposition is used to apply the Joukowski
condition
// The solution is seeked in the form psi0 + beta psi1 and beta
is
// adjusted so that the pressure is continuous at the trailing
edge
border a(t=0,2⋆pi) { x=5⋆cos(t); y=5⋆sin(t); }; // approximates
infinity
border upper(t=0,1) { x = t;
y = 0.17735⋆sqrt(t)-0.075597⋆t
- 0.212836⋆(t^2)+0.17363⋆(t^3)-0.06254⋆(t^4); }
border lower(t=1,0) { x = t;
y= -(0.17735⋆sqrt(t)-0.075597⋆t
-0.212836⋆(t^2)+0.17363⋆(t^3)-0.06254⋆(t^4)); }
border c(t=0,2⋆pi) { x=0.8⋆cos(t)+0.5; y=0.8⋆sin(t); }
wait = true;
mesh Zoom = buildmesh(c(30)+upper(35)+lower(35));
mesh Th = buildmesh(a(30)+upper(35)+lower(35));
fespace Vh(Th,P2); // P1 FE
space
Vh psi0,psi1,vh; // unknown and test
function.
fespace ZVh(Zoom,P2);
solve Joukowski0(psi0,vh) = // definition of the
problem
int2d(Th)( dx(psi0)⋆dx(vh) + dy(psi0)⋆dy(vh) ) // bilinear
form
+ on(a,psi0=y-0.1⋆x) // boundary condition
form
+ on(upper,lower,psi0=0);
plot(psi0);
solve Joukowski1(psi1,vh) = // definition of the
problem
int2d(Th)( dx(psi1)⋆dx(vh) + dy(psi1)⋆dy(vh) ) // bilinear
form
+ on(a,psi1=0) // boundary condition
form
+ on(upper,lower,psi1=1);
plot(psi1);
// continuity of pressure at trailing
edge
real beta = psi0(0.99,0.01)+psi0(0.99,-0.01);
beta = -beta / (psi1(0.99,0.01)+ psi1(0.99,-0.01)-2);
Vh psi = beta⋆psi1+psi0;
plot(psi);
ZVh Zpsi=psi;
plot(Zpsi,bw=true);
ZVh cp = -dx(psi)^2 - dy(psi)^2;
plot(cp);
ZVh Zcp=cp;
plot(Zcp,nbiso=40);
9.1.4 Error estimation
There are famous estimation between the numerical result uh and the exact solution u of the problem 2.1 and 2.2: If
triangulations {
h}h↓0 is regular (see Section 5.4), then we have the estimates
with constants C1,C2 independent of h, if u is in H2(Ω). It is known that u
H2(Ω) if Ω is convex.
In this section we check (9.14) and (9.15). We will pick up numericall error if we use the numerical derivative, so we
will use the following for (9.14).
The constants C1,C2 are depend on
h and f, so we will find them by FreeFem++ . In general, we cannot get the
solution u as a elementary functions (see Section 4.7) even if spetical functions are added. Instead of the exact
solution, here we use the approximate solution u0 in Vh(
h,P2),h ~ 0.
Example 9.4
1 : mesh Th0 = square(100,100);
2 : fespace V0h(Th0,P2);
3 : V0h u0,v0;
4 : func f = x⋆y; //
sin(pi⋆x)⋆cos(pi⋆y);
5 :
6 : solve Poisson0(u0,v0) =
7 : int2d(Th0)( dx(u0)⋆dx(v0) + dy(u0)⋆dy(v0) ) // bilinear
form
8 : - int2d(Th0)( f⋆v0 ) // linear
form
9 : + on(1,2,3,4,u0=0) ; // boundary
condition
10 :
11 : plot(u0);
12 :
13 : real[int] errL2(10), errH1(10);
14 :
15 : for (int i=1; i<=10; i++) {
16 : mesh Th = square(5+i⋆3,5+i⋆3);
17 : fespace Vh(Th,P1);
18 : fespace Ph(Th,P0);
19 : Ph h = hTriangle; // get the size of all
triangles
20 : Vh u,v;
21 : solve Poisson(u,v) =
22 : int2d(Th)( dx(u)⋆dx(v) + dy(u)⋆dy(v) ) // bilinear
form
23 : - int2d(Th)( f⋆v ) // linear
form
24 : + on(1,2,3,4,u=0) ; // boundary
condition
25 : V0h uu = u;
26 : errL2[i-1] = sqrt( int2d(Th0)((uu - u0)^2) )/h[].max^2;
27 : errH1[i-1] = sqrt( int2d(Th0)( f⋆(u0-2⋆uu+uu) ) )/h[].max;
28 : }
29 : cout << "C1 = " << errL2.max <<"("<<errL2.min<<")"<< endl;
30 : cout << "C2 = " << errH1.max <<"("<<errH1.min<<")"<< endl;
We can guess that C1 = 0.0179253(0.0173266) and C2 = 0.0729566(0.0707543), where the numbers inside the
parentheses are minimum in calculation.
9.1.5 Periodic
We now solve the Poisson equation
on the
square ]0,2π[2 under bi-periodic boundary condition u(0,y) = u(2π,y) for all y and u(x,0) = u(x,2π)
for all x. These boundary conditions are achieved from the definition of the periodic finite element
space.
Example 9.5 (periodic.edp)
mesh Th=square(10,10,[2⋆x⋆pi,2⋆y⋆pi]);
// defined the fespacewith periodic
condition
// label : 2 and 4 are left and right side with y the curve
abscissa
// 1 and 2 are bottom and upper side with x the curve
abscissa
fespace Vh(Th,P2,periodic=[[2,y],[4,y],[1,x],[3,x]]);
Vh uh,vh; // unknown and test
function.
func f=sin(x+pi/4.)⋆cos(y+pi/4.); // right hand side
function
problem laplace(uh,vh) = // definion of the
problem
int2d(Th)( dx(uh)⋆dx(vh) + dy(uh)⋆dy(vh) ) // bilinear
form
+ int2d(Th)( -f⋆vh ) // linear
form
;
laplace; // solve the problem plot(uh); // to see the
result
plot(uh,ps="period.eps",value=true);
The periodic condition does not necessarily require parallel to the axis. Example 9.6 give such example.
Example 9.6 (periodic4.edp)
real r=0.25;
// a diamond with a
hole
border a(t=0,1){x=-t+1; y=t;label=1;};
border b(t=0,1){ x=-t; y=1-t;label=2;};
border c(t=0,1){ x=t-1; y=-t;label=3;};
border d(t=0,1){ x=t; y=-1+t;label=4;};
border e(t=0,2⋆pi){ x=r⋆cos(t); y=-r⋆sin(t);label=0;};
int n = 10;
mesh Th= buildmesh(a(n)+b(n)+c(n)+d(n)+e(n));
plot(Th,wait=1);
real r2=1.732;
func abs=sqrt(x^2+y^2);
// warning for periodic condition:
// side a and c
// on side a (label 1) x
[0, 1] or x - y
[-1, 1]
// on side c (label 3) x
[-1, 0] or x - y
[-1, 1]
// so the common abscissa can be respectively x and
x + 1
// or you can can try curviline abscissa x - y and
x - y
// 1 first way
// fespace Vh(Th,P2,periodic=[[2,1+x],[4,x],[1,x],[3,1+x]]);
// 2 second way
fespace Vh(Th,P2,periodic=[[2,x+y],[4,x+y],[1,x-y],[3,x-y]]);
Vh uh,vh;
func f=(y+x+1)⋆(y+x-1)⋆(y-x+1)⋆(y-x-1);
real intf = int2d(Th)(f);
real mTh = int2d(Th)(1);
real k = intf/mTh;
problem laplace(uh,vh) =
int2d(Th)( dx(uh)⋆dx(vh) + dy(uh)⋆dy(vh) ) + int2d(Th)( (k-f)⋆vh ) ;
laplace;
plot(uh,wait=1,ps="perio4.eps");
A other example with no equal border, just to see if the code works.
Example 9.7 (periodic4bis.edp)
// irregular boundary
condition.
// to build border
AB
macro LINEBORDER(A,B,lab) border A#B(t=0,1){real t1=1.-t;
x=A#x⋆t1+B#x⋆t;y=A#y⋆t1+B#y⋆t;label=lab;} //
EOM
// compute ||AB|| a=(ax,ay) et B
=(bx,by)
macro dist(ax,ay,bx,by) sqrt(square((ax)-(bx))+ square((ay)-(by))) //
EOM
macro Grad(u) [dx(u),dy(u)] //
EOM
real Ax=0.9,Ay=1; real Bx=2,By=1;
real Cx=2.5,Cy=2.5; real Dx=1,Dy=2;
real gx = (Ax+Bx+Cx+Dx)/4.; real gy = (Ay+By+Cy+Dy)/4.;
LINEBORDER(A,B,1)
LINEBORDER(B,C,2)
LINEBORDER(C,D,3)
LINEBORDER(D,A,4)
int n=10;
real l1=dist(Ax,Ay,Bx,By);
real l2=dist(Bx,By,Cx,Cy);
real l3=dist(Cx,Cy,Dx,Dy);
real l4=dist(Dx,Dy,Ax,Ay);
func s1=dist(Ax,Ay,x,y)/l1; // absisse on AB =
||AX||/||AB||
func s2=dist(Bx,By,x,y)/l2; // absisse on BC =
||BX||/||BC||
func s3=dist(Cx,Cy,x,y)/l3; // absisse on CD =
||CX||/||CD||
func s4=dist(Dx,Dy,x,y)/l4; // absisse on DA =
||DX||/||DA||
mesh Th=buildmesh(AB(n)+BC(n)+CD(n)+DA(n),fixeborder=1); //
verbosity=6; // to see the abscisse value pour the periodic
condition.
fespace Vh(Th,P1,periodic=[[1,s1],[3,s3],[2,s2],[4,s4]]);
verbosity=1;
Vh u,v;
real cc=0;
cc= int2d(Th)((x-gx)⋆(y-gy)-cc)/Th.area;
cout << " compatibility =" << int2d(Th)((x-gx)⋆(y-gy)-cc) <<endl;
solve Poission(u,v)=int2d(Th)(Grad(u)'⋆Grad(v)+ 1e-10⋆u⋆v)
-int2d(Th)(10⋆v⋆((x-gx)⋆(y-gy)-cc));
plot(u,wait=1,value=1);
Example 9.8 (Period-Poisson-cube-ballon.edp)
verbosity=1;
load "msh3"
load "tetgen"
load "medit"
bool buildTh=0;
mesh3 Th;
try { // a way to build one time the mesh an read if the file
exist.
Th=readmesh3("Th-hex-sph.mesh");
}
catch(...) { buildTh=1;}
if( buildTh ){
...
put the code ??cube-ballon} example page //
225
without the first line
}
fespace Ph(Th,P0);
verbosity=50;
fespace Vh(Th,P1,periodic=[[3,x,z],[4,x,z],[1,y,z],[2,y,z],[5,x,y],[6,x,y]]);
// back and front
verbosity=1;
Ph reg=region;
cout << " centre = " << reg(0,0,0) << endl;
cout << " exterieur = " << reg(0,0,0.7) << endl;
macro Grad(u) [dx(u),dy(u),dz(u)] //
EOM
Vh uh,vh;
real x0=0.3,y0=0.4,z0=06;
func f= sin(x⋆2⋆pi+x0)⋆sin(y⋆2⋆pi+y0)⋆sin(z⋆2⋆pi+z0);
real gn = 1.;
real cf= 1;
problem P(uh,vh)=
int3d(Th,1)( Grad(uh)'⋆Grad(vh)⋆100)
+ int3d(Th,2)( Grad(uh)'⋆Grad(vh)⋆2)
+ int3d(Th) (vh⋆f)
;
P;
plot(uh,wait=1, nbiso=6);
medit(" uh ",Th, uh);
9.1.6 Poisson with mixed boundary condition
Here we consider the Poisson equation with mixed boundary value problems: For given functions f and g, find u
such that
where ΓD is a part of the boundary Γ and ΓN = Γ \ΓD. The solution u has the singularity at the points
{γ1,γ2} = ΓD ∩ΓN. When Ω = {(x,y); -1 < x < 1,0 < y < 1}, ΓN = {(x,y); -1 ≤ x < 0,y = 0}, ΓD = ∂Ω \ΓN, the
singularity will appear at γ1 = (0,0),γ2(-1,0), and u has the expression
with a
constants Ki. Here uS = rj1∕2 sin(θj∕2) by the local polar coordinate (rj,θj at γj such that (r1,θ1) = (r,θ). Instead of
poler coordinate system (r,θ), we use that r = sqrt( x
+y
)and θ = atan2(y,x)in FreeFem++
.
Example 9.9 Assume that f = -2×30(x2 +y2) and g = ue = 10(x2 +y2)1∕4 sin
+30(x2y2),
where ueS is the exact solution.
1 : border N(t=0,1) { x=-1+t; y=0; label=1; };
2 : border D1(t=0,1){ x=t; y=0; label=2;};
3 : border D2(t=0,1){ x=1; y=t; label=2; };
4 : border D3(t=0,2){ x=1-t; y=1; label=2;};
5 : border D4(t=0,1) { x=-1; y=1-t; label=2; };
6 :
7 : mesh T0h = buildmesh(N(10)+D1(10)+D2(10)+D3(20)+D4(10));
8 : plot(T0h,wait=true);
9 : fespace V0h(T0h,P1);
10 : V0h u0, v0;
11 :
12 : func f=-2⋆30⋆(x^2+y^2); // given
function
13 : // the singular term of the solution is K⋆us (K:
constant)
14 : func us = sin(atan2(y,x)/2)⋆sqrt( sqrt(x^2+y^2) );
15 : real K=10.;
16 : func ue = K⋆us + 30⋆(x^2⋆y^2);
17 :
18 : solve Poisson0(u0,v0) =
19 : int2d(T0h)( dx(u0)⋆dx(v0) + dy(u0)⋆dy(v0) ) // bilinear
form
20 : - int2d(T0h)( f⋆v0 ) // linear
form
21 : + on(2,u0=ue) ; // boundary
condition
22 :
23 : // adaptation by the singular
term
24 : mesh Th = adaptmesh(T0h,us);
25 : for (int i=0;i< 5;i++)
26 : {
27 : mesh Th=adaptmesh(Th,us);
28 : } ;
29 :
30 : fespace Vh(Th, P1);
31 : Vh u, v;
32 : solve Poisson(u,v) =
33 : int2d(Th)( dx(u)⋆dx(v) + dy(u)⋆dy(v) ) // bilinear
form
34 : - int2d(Th)( f⋆v ) // linear
form
35 : + on(2,u=ue) ; // boundary
condition
36 :
37 : /⋆ plot the solution ⋆/
38 : plot(Th,ps="adaptDNmix.ps");
39 : plot(u,wait=true);
40 :
41 : Vh uue = ue;
42 : real H1e = sqrt( int2d(Th)( dx(uue)^2 + dy(uue)^2 + uue^2 ) );
43 :
44 : /⋆ calculate the H1 Sobolev norm ⋆/
45 : Vh err0 = u0 - ue;
46 : Vh err = u - ue;
47 : Vh H1err0 = int2d(Th)( dx(err0)^2+dy(err0)^2+err0^2 );
48 : Vh H1err = int2d(Th)( dx(err)^2+dy(err)^2+err^2 );
49 : cout <<"Relative error in first mesh "<< int2d(Th)(H1err0)/H1e<<endl;
50 : cout <<"Relative error in adaptive mesh "<< int2d(Th)(H1err)/H1e<<endl;
From 24th line to 28th, adaptation of meshes are done using the base of singular term. In 42th line, H1e=∥ue∥1,Ω is
calculated. In last 2 lines, the relative errors are calculated, that is,
where uh0 is the numerical solution in T0hand uha is uin this program.
9.1.7 Poisson with mixte finite element
Here we consider the Poisson equation with mixed boundary value problems: For given functions f , gd, gn, find p
such that
where ΓD is a part of the boundary Γ and ΓN = Γ \ΓD.
The mixte formulation is: find p and u such that
where gn is a vector such that gn.n = gn.
The variationnal formulation is,
where the functionnal space are:
and
To write, the FreeFem++ example, we have just to choose the finites elements spaces. here V space
is discretize with Raviart-Thomas finite element RT0and ℙ is discretize by constant finite element
P0.
Example 9.10 (LaplaceRT.edp)
mesh Th=square(10,10);
fespace Vh(Th,RT0);
fespace Ph(Th,P0);
func gd = 1.;
func g1n = 1.;
func g2n = 1.;
Vh [u1,u2],[v1,v2];
Ph p,q;
problem laplaceMixte([u1,u2,p],[v1,v2,q],
solver=GMRES,eps=1.0e-10,
tgv=1e30,dimKrylov=150)
=
int2d(Th)( p⋆q⋆1e-15 // this term is here to be
sur
// that all sub matrix are inversible (LU
requirement)
+ u1⋆v1 + u2⋆v2 + p⋆(dx(v1)+dy(v2)) + (dx(u1)+dy(u2))⋆q )
+ int2d(Th) ( q)
- int1d(Th,1,2,3)( gd⋆(v1⋆N.x +v2⋆N.y)) // on
ΓD
+ on(4,u1=g1n,u2=g2n); // on
ΓN
laplaceMixte;
plot([u1,u2],coef=0.1,wait=1,ps="lapRTuv.eps",value=true);
plot(p,fill=1,wait=1,ps="laRTp.eps",value=true);
9.1.8 Metric Adaptation and residual error indicator
We do metric mesh adaption and compute the classical residual error indicator ηT on the element T for the Poisson
problem.
Example 9.11 (adaptindicatorP2.edp) First, we solve the same problem as in a previous example.
1 : border ba(t=0,1.0){x=t; y=0; label=1;}; // see
Fig,5.13
2 : border bb(t=0,0.5){x=1; y=t; label=2;};
3 : border bc(t=0,0.5){x=1-t; y=0.5;label=3;};
4 : border bd(t=0.5,1){x=0.5; y=t; label=4;};
5 : border be(t=0.5,1){x=1-t; y=1; label=5;};
6 : border bf(t=0.0,1){x=0; y=1-t;label=6;};
7 : mesh Th = buildmesh (ba(6) + bb(4) + bc(4) +bd(4) + be(4) + bf(6));
8 : savemesh(Th,"th.msh");
9 : fespace Vh(Th,P2);
10 : fespace Nh(Th,P0);
11 : Vh u,v;
12 : Nh rho;
13 : real[int] viso(21);
14 : for (int i=0;i<viso.n;i++)
15 : viso[i]=10.^(+(i-16.)/2.);
16 : real error=0.01;
17 : func f=(x-y);
18 : problem Probem1(u,v,solver=CG,eps=1.0e-6) =
19 : int2d(Th,qforder=5)( u⋆v⋆1.0e-10+ dx(u)⋆dx(v) + dy(u)⋆dy(v))
20 : + int2d(Th,qforder=5)( -f⋆v);
21 : /⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆
Now, the local error indicator ηT is:
where hT is the longest’s edge of T,
T is the set of T edge not on Γ = ∂Ω, nT is the outside unit normal to
K, he is the length of edge e, [g] is the jump of the function g across edge (left value minus right value).
Of coarse, we can use a variational form to compute ηT2, with test function constant function in each
triangle.
29 : ⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆⋆/
30 :
31 : varf indicator2(uu,chiK) =
32 : intalledges(Th)(chiK⋆lenEdge⋆square(jump(N.x⋆dx(u)+N.y⋆dy(u))))
33 : +int2d(Th)(chiK⋆square(hTriangle⋆(f+dxx(u)+dyy(u))) );
34 : for (int i=0;i< 4;i++)
35 : {
36 : Probem1;
37 : cout << u[].min << " " << u[].max << endl;
38 : plot(u,wait=1);
39 : cout << " indicator2 " << endl;
40 :
41 : rho[] = indicator2(0,Nh);
42 : rho=sqrt(rho);
43 : cout << "rho = min " << rho[].min << " max=" << rho[].max << endl;
44 : plot(rho,fill=1,wait=1,cmm="indicator density ",ps="rhoP2.eps",
value=1,viso=viso,nbiso=viso.n);
45 : plot(Th,wait=1,cmm="Mesh ",ps="ThrhoP2.eps");
46 : Th=adaptmesh(Th,[dx(u),dy(u)],err=error,anisomax=1);
47 : plot(Th,wait=1);
48 : u=u;
49 : rho=rho;
50 : error = error/2;
51 : } ;
If the method is correct, we expect to look the graphics by an almost constant function η on your computer as in Fig.
9.11.
9.1.9 Adaptation using residual error indicator
In the previous example we compute the error indicator, now we use it, to adapt the mesh.
The new mesh size is given by the following formulae:
where
ηn(x) is the level of error at point x given by the local error indicator, hn is the previous “mesh size” field, and fn is a
user function define by fn = min(3,max(1∕3,ηn∕ηn*)) where ηn* = mean(ηn)c, and c is an user coefficient generally
close to one.
Example 9.12 (AdaptResidualErrorIndicator.edp)
First a macro MeshSizecomputationto get a P1 mesh size as the average of edge lenght.
// macro the get the current mesh size
// parameter
// in: Th the mesh
// Vh P1 fespace on Th
// out :
// h: the Vh finite element finite set to the current mesh size
macro MeshSizecomputation(Th,Vh,h)
{ /⋆ Th mesh Vh P1 finite element space
h the P1 mesh size value ⋆/
real[int] count(Th.nv);
/⋆ mesh size (lenEdge = integral(e) 1 ds) ⋆/
varf vmeshsizen(u,v)=intalledges(Th,qfnbpE=1)(v);
/⋆ number of edge / par vertex ⋆/
varf vedgecount(u,v)=intalledges(Th,qfnbpE=1)(v/lenEdge);
/⋆
computation of the mesh size
----------------------------- ⋆/
count=vedgecount(0,Vh);
h[]=0.;
h[]=vmeshsizen(0,Vh);
cout << " count min = "<< count.min << " " << count.max << endl;
h[]=h[]./count;
cout << " -- bound meshsize = " <<h[].min << " " << h[].max << endl;
} // end of macro MeshSizecomputation
A second macro to remesh according to the new mesh size.
// macro to remesh according the de residual indicator
// in:
// Th the mesh
// Ph P0 fespace on Th
// Vh P1 fespace on Th
// vindicator the varf of to evaluate the indicator to 2
// coef on etameam ..
// ------
macro ReMeshIndicator(Th,Ph,Vh,vindicator,coef)
{
Vh h=0;
/⋆evalutate the mesh size ⋆/
MeshSizecomputation(Th,Vh,h);
Ph etak;
etak[]=vindicator(0,Ph);
etak[]=sqrt(etak[]);
real etastar= coef⋆(etak[].sum/etak[].n);
cout << " etastar = " << etastar << " sum=" << etak[].sum << " " << endl;
/⋆ here etaK is discontinous
we use the P1 L2 projection with mass lumping . ⋆/
Vh fn,sigma;
varf veta(unused,v)=int2d(Th)(etak⋆v);
varf vun(unused,v)=int2d(Th)(1⋆v);
fn[] = veta(0,Vh);
sigma[]= vun(0,Vh);
fn[]= fn[]./ sigma[];
fn = max(min(fn/etastar,3.),0.3333) ;
/⋆ new mesh size ⋆/
h = h / fn ;
/⋆ plot(h,wait=1); ⋆/
/⋆ build the new mesh ⋆/
Th=adaptmesh(Th,IsMetric=1,h,splitpbedge=1,nbvx=10000);
}
We skip the mesh construction, see the previous example,
// FE space definition
---
fespace Vh(Th,P1); // for the mesh size and
solution
fespace Ph(Th,P0); // for the error
indicator
real hinit=0.2; // initial mesh
size
Vh h=hinit; // the FE function for the mesh
size
// to build a mesh with a given mesh size :
meshsize
Th=adaptmesh(Th,h,IsMetric=1,splitpbedge=1,nbvx=10000);
plot(Th,wait=1,ps="RRI-Th-init.eps");
Vh u,v;
func f=(x-y);
problem Poisson(u,v) =
int2d(Th,qforder=5)( u⋆v⋆1.0e-10+ dx(u)⋆dx(v) + dy(u)⋆dy(v))
- int2d(Th,qforder=5)( f⋆v);
varf indicator2(unused,chiK) =
intalledges(Th)(chiK⋆lenEdge⋆square(jump(N.x⋆dx(u)+N.y⋆dy(u))))
+int2d(Th)(chiK⋆square(hTriangle⋆(f+dxx(u)+dyy(u))) );
for (int i=0;i< 10;i++)
{
u=u;
Poisson;
plot(Th,u,wait=1);
real cc=0.8;
if(i>5) cc=1;
ReMeshIndicator(Th,Ph,Vh,indicator2,cc);
plot(Th,wait=1);
}
9.2 Elasticity
Consider an elastic plate with undeformed shape Ω×] -h,h[ in
3, Ω ⊂
2. By the deformation of the plate, we
assume that a point P(x1,x2,x3) moves to
(ξ1,ξ2,ξ3). The vector
= (u1,u2,u3) = (ξ1 -x1,ξ2 -x2,ξ3 -x3) is
called displacement vector. By the deformation, the line segment x,x + τΔx moves approximately to
x + u(x),x + τΔx + u(x + τΔx) for small τ, where x = (x1,x2,x3),Δx = (Δx1,Δx2,Δx3). We now calculate the ratio
between two segments
then we have (see e.g. [16, p.32])
where νi = Δxi|Δx|-1. If the deformation is small, then we may consider that
and the following is called small strain tensor
The tensor eij is called finite strain tensor.
Consider the small plane ΔΠ(x) centered at x with the unit normal direction
= (n1,n2,n3), then the surface on
ΔΠ(x) at x is
where σij(x) is called stress tensor at x. Hooke’s law is the assumption of a linear relation between σij and εij such
as
with the symmetry cijkl = cjikl,cijkl = cijlk,cijkl = cklij.
If Hooke’s tensor cijkl(x) do not depend on the choice of coordinate system, the material is called isotropic at x. If
cijkl is constant, the material is called homogeneous. In homogeneous isotropic case, there is Lamé constants λ,μ
(see e.g. [16, p.43]) satisfying
where δij is Kronecker’s delta. We assume that the elastic plate is fixed on ΓD×] -h,h[,ΓD ⊂ ∂Ω. If the body force
f = (f1,f2,f3) is given in Ω×] -h,h[ and surface force g is given in ΓN×] -h,h[,ΓN = ∂Ω \ΓD, then the equation
of equilibrium is given as follows:
We now explain the plain elasticity.
-
Plain strain:
- On the end of plate, the contact condition u3 = 0,g3 = is satisfied. In this case, we can
suppose that f3 = g3 = u3 = 0 and
(x1,x2,x3) = u(x1,x2) for all -h < x3 < h.
-
Plain stress:
- The cylinder is assumed to be very thin and subjected to no load on the ends x3 = ±h, that is,
The assumption leads that σ3i = 0 in Ω×] -h,h[ and
(x1,x2,x3) = u(x1,x2) for all -h < x3 < h.
-
Generalized plain stress:
- The cylinder is subjected to no load on the ends x3 = ±h. Introducing the mean
values with respect to thickness,
and we derive u3 ≡ 0. Similarly we define the mean values f,g of the body force and surface force
as well as the mean values εij and σij of the components of stress and strain, respectively.
In what follows we omit the overlines of u,f,g,εij and εij. Then we obtain similar equation of equilibrium
given in (9.21) replacing Ω×] -h,h[ with Ω and changing i = 1,2. In the case of plane stress,
σij = λ*δijdivu + 2μεij,λ* = (2λμ)∕(λ + μ).
The equations of elasticity are naturally written in variational form for the displacement vector u(x)
V as
where
V is the linear closed subspace of H1(Ω)2.
Example 9.13 (Beam.edp) Consider elastic plate with the undeformed rectangle shape [0,10] ×[0,2].
The body force is the gravity force
and the boundary force
is zero on lower and upper side. On the two
vertical sides of the beam are fixed.
// a weighting beam sitting on
a
int bottombeam = 2;
border a(t=2,0) { x=0; y=t ;label=1;}; // left
beam
border b(t=0,10) { x=t; y=0 ;label=bottombeam;}; // bottom of
beam
border c(t=0,2) { x=10; y=t ;label=1;}; // rigth
beam
border d(t=0,10) { x=10-t; y=2; label=3;}; // top
beam
real E = 21.5;
real sigma = 0.29;
real mu = E/(2⋆(1+sigma));
real lambda = E⋆sigma/((1+sigma)⋆(1-2⋆sigma));
real gravity = -0.05;
mesh th = buildmesh( b(20)+c(5)+d(20)+a(5));
fespace Vh(th,[P1,P1]);
Vh [uu,vv], [w,s];
cout << "lambda,mu,gravity ="<<lambda<< " " << mu << " " << gravity << endl;
// deformation of a beam under its own
weight
real sqrt2=sqrt(2.); // see lame.edp example
3.9
macro epsilon(u1,u2) [dx(u1),dy(u2),(dy(u1)+dx(u2))/sqrt2] //
EOM
macro div(u,v) ( dx(u)+dy(v) ) //
EOM
solve bb([uu,vv],[w,s])=
int2d(th)(
lambda⋆div(w,s)⋆div(uu,vv)
+2.⋆mu⋆( epsilon(w,s)'⋆epsilon(uu,vv) )
)
+ int2d(th) (-gravity⋆s)
+ on(1,uu=0,vv=0)
;
plot([uu,vv],wait=1);
plot([uu,vv],wait=1,bb=[[-0.5,2.5],[2.5,-0.5]]);
mesh th1 = movemesh(th, [x+uu, y+vv]);
plot(th1,wait=1);
9.2.1 Fracture Mechanics
Consider the plate with the crack whose undeformed shape is a curve Σ with the two edges γ1,γ2. We assume
the stress tensor σij is the state of plate stress regarding (x,y)
ΩΣ = Ω \Σ. Here Ω stands for the
undeformed shape of elastic plate without crack. If the part ΓN of the boundary ∂Ω is fixed and a load
= (
,
)
L2(Ω)2 ×L2(ΓN)2 is given, then the displacement
is the minimizer of the potential energy
functional
over the functional space V(ΩΣ),
where w(x,
) = σij(
)εij(
)∕2,
If the elasticity is homogeneous isotropic, then the displacement
(x) is decomposed in an open neighborhood Uk
of γk as in (see e.g. [17])
 | (9.23) |
with
k,R
H2(ΩΣ ∩Uk)2, where Uk,k = 1,2 are open neighborhoods of γk such that ∂L1 ∩U1 = γ1,∂Lm ∩U2 = γ2,
and
where μ is the shear modulus of elasticity, κ = 3 -4ν (ν is the Poisson’s ratio) for plane strain and κ =
for plane
stress.
The coefficients K1(γi) and K2(γi), which are important parameters in fracture mechanics, are called stress intensity
factors of the opening mode (mode I) and the sliding mode (mode II), respectively.
For simplicity, we consider the following simple crack
with only one crack tip γ = (0,0). Unfortunately, FreeFem++ cannot treat crack, so we use the modification of
the domain with U-shape channel (see Fig. 5.28) with d = 0.0001. The undeformed crack Σ is approximated by
and ΓD = Rin Fig. 5.28. In this example, we use three technique:
First example create mode I deformation by the opposed surface force on Band Tin the vertical direction of Σ, and
the displacement is fixed on R.
In a laboratory, fracture engineer use photoelasticity to make stress field visible, which shows the principal stress
difference
where σ1 and σ2 are the principal stresses. In opening mode, the photoelasticity make symmetric pattern
concentrated at γ.
Example 9.14 (Crack Opening, K2(γ) = 0)
{CrackOpen.edp}
real d = 0.0001;
int n = 5;
real cb=1, ca=1, tip=0.0;
border L1(t=0,ca-d) { x=-cb; y=-d-t; }
border L2(t=0,ca-d) { x=-cb; y=ca-t; }
border B(t=0,2) { x=cb⋆(t-1); y=-ca; }
border C1(t=0,1) { x=-ca⋆(1-t)+(tip-10⋆d)⋆t; y=d; }
border C21(t=0,1) { x=(tip-10⋆d)⋆(1-t)+tip⋆t; y=d⋆(1-t); }
border C22(t=0,1) { x=(tip-10⋆d)⋆t+tip⋆(1-t); y=-d⋆t; }
border C3(t=0,1) { x=(tip-10⋆d)⋆(1-t)-ca⋆t; y=-d; }
border C4(t=0,2⋆d) { x=-ca; y=-d+t; }
border R(t=0,2) { x=cb; y=cb⋆(t-1); }
border T(t=0,2) { x=cb⋆(1-t); y=ca; }
mesh Th = buildmesh (L1(n/2)+L2(n/2)+B(n)
+C1(n)+C21(3)+C22(3)+C3(n)+R(n)+T(n));
cb=0.1; ca=0.1;
plot(Th,wait=1);
mesh Zoom = buildmesh (L1(n/2)+L2(n/2)+B(n)+C1(n)
+C21(3)+C22(3)+C3(n)+R(n)+T(n));
plot(Zoom,wait=1);
real E = 21.5;
real sigma = 0.29;
real mu = E/(2⋆(1+sigma));
real lambda = E⋆sigma/((1+sigma)⋆(1-2⋆sigma));
fespace Vh(Th,[P2,P2]);
fespace zVh(Zoom,P2);
Vh [u,v], [w,s];
solve Problem([u,v],[w,s]) =
int2d(Th)(
2⋆mu⋆(dx(u)⋆dx(w)+ ((dx(v)+dy(u))⋆(dx(s)+dy(w)))/4 )
+ lambda⋆(dx(u)+dy(v))⋆(dx(w)+dy(s))/2
)
-int1d(Th,T)(0.1⋆(4-x)⋆s)+int1d(Th,B)(0.1⋆(4-x)⋆s)
+on(R,u=0)+on(R,v=0); //
fixed
;
zVh Sx, Sy, Sxy, N;
for (int i=1; i<=5; i++)
{
mesh Plate = movemesh(Zoom,[x+u,y+v]); // deformation near
γ
Sx = lambda⋆(dx(u)+dy(v)) + 2⋆mu⋆dx(u);
Sy = lambda⋆(dx(u)+dy(v)) + 2⋆mu⋆dy(v);
Sxy = mu⋆(dy(u) + dx(v));
N = 0.1⋆1⋆sqrt((Sx-Sy)^2+4⋆Sxy^2); // principal stress
difference
if (i==1) {
plot(Plate,ps="1stCOD.eps",bw=1); // Fig.
9.13
plot(N,ps="1stPhoto.eps",bw=1); // Fig.
9.13
} else if (i==5) {
plot(Plate,ps="LastCOD.eps",bw=1); // Fig.
9.14
plot(N,ps="LastPhoto.eps",bw=1); // Fig.
9.14
break;
}
Th=adaptmesh(Th,[u,v]);
Problem;
}
It is difficult to create mode II deformation by the opposed shear force on Band Tthat is observed in a
laboratory. So we use the body shear force along Σ, that is, the x-component f1 of the body force
is given
by
where H(t) = 1 if t > 0; = 0 if t < 0.
Example 9.15 (Crack Sliding, K2(γ) = 0)
(use the same mesh Th)
cb=0.01; ca=0.01;
mesh Zoom = buildmesh (L1(n/2)+L2(n/2)+B(n)+C1(n)
+C21(3)+C22(3)+C3(n)+R(n)+T(n));
(use same FE-space Vh and elastic modulus)
fespace Vh1(Th,P1);
Vh1 fx = ((y>0.001)⋆(y<0.1))-((y<-0.001)⋆(y>-0.1)) ;
solve Problem([u,v],[w,s]) =
int2d(Th)(
2⋆mu⋆(dx(u)⋆dx(w)+ ((dx(v)+dy(u))⋆(dx(s)+dy(w)))/4 )
+ lambda⋆(dx(u)+dy(v))⋆(dx(w)+dy(s))/2
)
-int2d(Th)(fx⋆w)
+on(R,u=0)+on(R,v=0); //
fixed
;
for (int i=1; i<=3; i++)
{
mesh Plate = movemesh(Zoom,[x+u,y+v]); // deformation near
γ
Sx = lambda⋆(dx(u)+dy(v)) + 2⋆mu⋆dx(u);
Sy = lambda⋆(dx(u)+dy(v)) + 2⋆mu⋆dy(v);
Sxy = mu⋆(dy(u) + dx(v));
N = 0.1⋆1⋆sqrt((Sx-Sy)^2+4⋆Sxy^2); // principal stress
difference
if (i==1) {
plot(Plate,ps="1stCOD2.eps",bw=1); // Fig.
9.16
plot(N,ps="1stPhoto2.eps",bw=1); // Fig.
9.15
} else if (i==3) {
plot(Plate,ps="LastCOD2.eps",bw=1); // Fig.
9.16
plot(N,ps="LastPhoto2.eps",bw=1); // Fig.
9.16
break;
}
Th=adaptmesh(Th,[u,v]);
Problem;
}
9.3 Nonlinear Static Problems
We propose how to solve the following non-linear academic problem of minimization of a functional
where
u is function of H01(Ω) and f defined by
9.3.1 Newton-Raphson algorithm
Now, we solve the Euler problem ∇J(u) = 0 with Newton-Raphson algorithm, that is,
First we introduice the two variational form vdJand vhJto compute respectively ∇J and ∇2J
// method of Newton-Raphson to solve dJ(u)=0;
//
// ---------------------------------------------
Ph dalpha ; // to store 2f′′(|∇u|2)
optimisation
// the variational form of evaluate dJ = ∇J
// --------------------------------------
// dJ = f'()⋆( dx(u)⋆dx(vh) + dy(u)⋆dy(vh)
varf vdJ(uh,vh) = int2d(Th)( alpha⋆( dx(u)⋆dx(vh) + dy(u)⋆dy(vh) ) - b⋆vh)
+ on(1,2,3,4, uh=0);
// the variational form of evaluate ddJ = ∇2J
// hJ(uh,vh) = f'()⋆( dx(uh)⋆dx(vh) + dy(uh)⋆dy(vh)
//
+ 2⋆f''()( dx(u)⋆dx(uh) + dy(u)⋆dy(uh) ) ⋆ (dx(u)⋆dx(vh) + dy(u)⋆dy(vh))
varf vhJ(uh,vh) = int2d(Th)( alpha⋆( dx(uh)⋆dx(vh) + dy(uh)⋆dy(vh) )
+ dalpha⋆( dx(u)⋆dx(vh) + dy(u)⋆dy(vh) )⋆( dx(u)⋆dx(uh) + dy(u)⋆dy(uh) ) )
+ on(1,2,3,4, uh=0);
// the Newton algorithm
Vh v,w;
u=0;
for (int i=0;i<100;i++)
{
alpha = df( dx(u)⋆dx(u) + dy(u)⋆dy(u) ) ; //
optimization
dalpha = 2⋆ddf( dx(u)⋆dx(u) + dy(u)⋆dy(u) ) ; //
optimization
v[]= vdJ(0,Vh); //
v = ∇J(u)
real res= v[]'⋆v[]; // the dot
product
cout << i << " residu^2 = " << res << endl;
if( res< 1e-12) break;
matrix H= vhJ(Vh,Vh,factorize=1,solver=LU); //
w[]=H^-1⋆v[];
u[] -= w[];
}
plot (u,wait=1,cmm="solution with Newton-Raphson");
Remark: This example is in Newton.edpfile of examples++-tutorialdirectory.
9.4 Eigenvalue Problems
This section depend on your FreeFem++compilation process (see README_arpack), of this tools. This tools is
available in FreeFem++if the word “eigenvalue” appear in line “Load:”, like:
-- FreeFem++ v1.28 (date Thu Dec 26 10:56:34 CET 2002)
file : LapEigenValue.edp
Load: lg_fem lg_mesh eigenvalue
This tools is based on the arpack++
the object-oriented version of ARPACK eigenvalue package [1].
The function EigenValue compute the generalized eigenvalue of Au = λBu where sigma =σ is the shift of the
method. The matrix OP is defined with A -σB. The return value is the number of converged eigenvalue (can be
greater than the number of eigen value nev=)
int k=EigenValue(OP,B,nev= , sigma= );
where the matrix OP = A -σB with a solver and boundary condition, and the matrix B.
Note 9.1 Boundary condition and Eigenvalue Problems
The lock (Dirichlet ) boundary condition is make with exact penalization so we put 1e30=tgv on the diagonal
term of the lock degree of freedom (see equation (6.17)). So take Dirichlet boundary condition just on A and
not on B. because we solve w = OP-1 *B *v.
If you put lock (Dirichlet ) boundary condition on B matrix (with key work on) you get small spurious modes
(10-30), du to boundary condition, but if you forget the lock boundary condition on B matrix (no key work
”on”) you get huge spurious (1030) modes associated to boundary conditon. We compute only small mode,
so we get the good one in this case.
-
sym=
- the problem is symmetric (all the eigen value are real)
-
nev=
- the number desired eigenvalues (nev) close to the shift.
-
value=
- the array to store the real part of the eigenvalues
-
ivalue=
- the array to store the imag. part of the eigenvalues
-
vector=
- the FE function array to store the eigenvectors
-
rawvector=
- an array of type real[int,int] to store eigenvectors by column. (up to version 2-17).
For real nonsymmetric problems, complex eigenvectors are given as two consecutive vectors, so if
eigenvalue k and k + 1 are complex conjugate eigenvalues, the kth vector will contain the real part
and the k + 1th vector the imaginary part of the corresponding complex conjugate eigenvectors.
-
tol=
- the relative accuracy to which eigenvalues are to be determined;
-
sigma=
- the shift value;
-
maxit=
- the maximum number of iterations allowed;
-
ncv=
- the number of Arnoldi vectors generated at each iteration of ARPACK.
Example 9.16 (lapEignenValue.edp) In the first example, we compute the eigenvalue and the eigenvector
of the Dirichlet problem on square Ω =]0,π[2.
The problem is find: λ, and ∇uλ in
×H01(Ω)
The exact eigenvalues are λn,m = (n2 + m2),(n,m)
ℕ*2 with the associated eigenvectors are um,n =
sin(nx) *sin(my).
We use the generalized inverse shift mode of the arpack++library, to find 20 eigenvalue and eigenvector
close to the shift value σ = 20.
// Computation of the eigen value and eigen vector of the
// Dirichlet problem on square ]0,π[2
// ----------------------------------------
// we use the inverse shift mode
// the shift is given with the real sigma
// -------------------------------------
// find λ and uλ
H10(Ω) such that:
// ∫
Ω ∇uλ∇v = λ∫
Ω uλv,∀v
H1
0 (Ω)
verbosity=10;
mesh Th=square(20,20,[pi⋆x,pi⋆y]);
fespace Vh(Th,P2);
Vh u1,u2;
real sigma = 20; // value of the
shift
// OP = A - sigma B ; // the shifted
matrix
varf op(u1,u2)= int2d(Th)( dx(u1)⋆dx(u2) + dy(u1)⋆dy(u2) - sigma⋆ u1⋆u2 )
+ on(1,2,3,4,u1=0) ; // Boundary
condition
varf b([u1],[u2]) = int2d(Th)( u1⋆u2 ); // no Boundary condition see note
9.1
matrix OP= op(Vh,Vh,solver=Crout,factorize=1); // crout solver because the
matrix in not positive
matrix B= b(Vh,Vh,solver=CG,eps=1e-20);
// important
remark:
// the boundary condition is make with exact
penalization:
// we put 1e30=tgv on the diagonal term of the lock degree of
freedom.
// So take Dirichlet boundary condition just on a variational
form
// and not on b variational
form.
// because we solve
w = OP-1 * B * v
int nev=20; // number of computed eigen value close to
sigma
real[int] ev(nev); // to store the nev
eigenvalue
Vh[int] eV(nev); // to store the nev eigenvector
int k=EigenValue(OP,B,sym=true,sigma=sigma,value=ev,vector=eV,
tol=1e-10,maxit=0,ncv=0);
// tol= the tolerance
// maxit= the maximum iteration see arpack doc.
// ncv see arpack doc. http://www.caam.rice.edu/software/ARPACK/
// the return value is number of converged eigen value.
for (int i=0;i<k;i++)
{
u1=eV[i];
real gg = int2d(Th)(dx(u1)⋆dx(u1) + dy(u1)⋆dy(u1));
real mm= int2d(Th)(u1⋆u1) ;
cout << " ---- " << i<< " " << ev[i]<< " err= "
<<int2d(Th)(dx(u1)⋆dx(u1) + dy(u1)⋆dy(u1) - (ev[i])⋆u1⋆u1) << " --- "<<endl;
plot(eV[i],cmm="Eigen Vector "+i+" valeur =" + ev[i] ,wait=1,value=1);
}
The output of this example is:
Nb of edges on Mortars = 0
Nb of edges on Boundary = 80, neb = 80
Nb Of Nodes = 1681
Nb of DF = 1681
Real symmetric eigenvalue problem: A⋆x - B⋆x⋆lambda
Thanks to ARPACK++ class ARrcSymGenEig
Real symmetric eigenvalue problem: A⋆x - B⋆x⋆lambda
Shift and invert mode sigma=20
Dimension of the system : 1681
Number of 'requested' eigenvalues : 20
Number of 'converged' eigenvalues : 20
Number of Arnoldi vectors generated: 41
Number of iterations taken : 2
Eigenvalues:
lambda[1]: 5.0002
lambda[2]: 8.00074
lambda[3]: 10.0011
lambda[4]: 10.0011
lambda[5]: 13.002
lambda[6]: 13.0039
lambda[7]: 17.0046
lambda[8]: 17.0048
lambda[9]: 18.0083
lambda[10]: 20.0096
lambda[11]: 20.0096
lambda[12]: 25.014
lambda[13]: 25.0283
lambda[14]: 26.0159
lambda[15]: 26.0159
lambda[16]: 29.0258
lambda[17]: 29.0273
lambda[18]: 32.0449
lambda[19]: 34.049
lambda[20]: 34.0492
---- 0 5.0002 err= -0.000225891 ---
---- 1 8.00074 err= -0.000787446 ---
---- 2 10.0011 err= -0.00134596 ---
---- 3 10.0011 err= -0.00134619 ---
---- 4 13.002 err= -0.00227747 ---
---- 5 13.0039 err= -0.004179 ---
---- 6 17.0046 err= -0.00623649 ---
---- 7 17.0048 err= -0.00639952 ---
---- 8 18.0083 err= -0.00862954 ---
---- 9 20.0096 err= -0.0110483 ---
---- 10 20.0096 err= -0.0110696 ---
---- 11 25.014 err= -0.0154412 ---
---- 12 25.0283 err= -0.0291014 ---
---- 13 26.0159 err= -0.0218532 ---
---- 14 26.0159 err= -0.0218544 ---
---- 15 29.0258 err= -0.0311961 ---
---- 16 29.0273 err= -0.0326472 ---
---- 17 32.0449 err= -0.0457328 ---
---- 18 34.049 err= -0.0530978 ---
---- 19 34.0492 err= -0.0536275 ---
9.5 Evolution Problems
FreeFem++ also solve evolution problems such as the heat problem
with a positive viscosity coefficient μ and homogeneous Neumann boundary conditions. We solve (9.26) by FEM in
space and finite differences in time. We use the definition of the partial derivative of the solution in the time
derivative,
which indicate that um(x,y) = u(x,y,mτ) imply
The time discretization of heat equation (9.27) is as follows:
which is so-called backward Euler method for (9.27). Multiplying the test function v both sides of the formula just
above, we have
By the divergence theorem, we have
By the boundary condition ∂um+1∕∂n = 0, it follows that
 | (9.28) |
Using the identity just above, we can calculate the finite element approximation uhm of um in a step-by-step manner
with respect to t.
Example 9.17 We now solve the following example with the exact solution u(x,y,t) = tx4.
// heat equation
∂tu = -μΔu = x4 - μ12tx2
mesh Th=square(16,16);
fespace Vh(Th,P1);
Vh u,v,uu,f,g;
real dt = 0.1, mu = 0.01;
problem dHeat(u,v) =
int2d(Th)( u⋆v + dt⋆mu⋆(dx(u)⋆dx(v) + dy(u)⋆dy(v)))
+ int2d(Th) (- uu⋆v - dt⋆f⋆v )
+ on(1,2,3,4,u=g);
real t = 0; // start from
t=0
uu = 0; //
u(x,y,0)=0
for (int m=0;m<=3/dt;m++)
{
t=t+dt;
f = x^4-mu⋆t⋆12⋆x^2;
g = t⋆x^4;
dHeat;
plot(u,wait=true);
uu = u;
cout <<"t="<<t<<"L^2-Error="<<sqrt( int2d(Th)((u-t⋆x^4)^2) ) << endl;
}
In the last statement, the L2-error
1∕2 is calculated at t = mτ,τ = 0.1. At t = 0.1, the error is
0.000213269. The errors increase with m and 0.00628589 at t = 3.
The iteration of the backward Euler (9.28) is made by for loop (see Section 4.9).
Note 9.2 The stiffness matrix in loop is used over and over again. FreeFem++ support reuses of stiffness
matrix.
9.5.1 Mathematical Theory on Time Difference Approximations.
In this section, we show the advantage of implicit schemes. Let V,H be separable Hilbert space and V is dense in H.
Let a be a continuous bilinear form over V ×V with coercivity and symmetry. Then
become equivalent to
the norm ∥v∥of V.
Problem Ev(f,Ω): For a given f
L2(0,T; V′),u0
H
where V′is the dual space of V. Then, there is an unique solution u
L∞(0,T; H) ∩L2(0,T; V).
Let us denote the time step by τ > 0, NT = [T∕τ]. For the discretization, we put un = u(nτ) and consider the time
difference for each θ
[0,1]
Formula (9.30) is the forward Euler scheme if θ = 0, Crank-Nicolson scheme if θ = 1∕2, the backward Euler scheme
if θ = 1.
Unknown vectors un = (uh1,
,uhM)T in
are obtained from solving the matrix
Refer [22, pp.70–75] for solvability of (9.31). The stability of (9.31) is in [22, Theorem 2.13]:
Example 9.18
mesh Th=square(12,12);
fespace Vh(Th,P1);
fespace Ph(Th,P0);
Ph h = hTriangle; // mesh sizes for each
triangle
real tau = 0.1, theta=0.;
func real f(real t) {
return x^2⋆(x-1)^2 + t⋆(-2 + 12⋆x - 11⋆x^2 - 2⋆x^3 + x^4);
}
ofstream out("err02.csv"); // file to store
calculations
out << "mesh size = "<<h[].max<<", time step = "<<tau<<endl;
for (int n=0;n<5/tau;n++) \\
out<<n⋆tau<<",";
out << endl;
Vh u,v,oldU;
Vh f1, f0;
problem aTau(u,v) =
int2d(Th)( u⋆v + theta⋆tau⋆(dx(u)⋆dx(v) + dy(u)⋆dy(v) + u⋆v))
- int2d(Th)(oldU⋆v - (1-theta)⋆tau⋆(dx(oldU)⋆dx(v)+dy(oldU)⋆dy(v)+oldU⋆v))
- int2d(Th)(tau⋆( theta⋆f1+(1-theta)⋆f0 )⋆v );
while (theta <= 1.0) {
real t = 0, T=3; // from t=0 to
T
oldU = 0; //
u(x,y,0)=0
out <<theta<<",";
for (int n=0;n<T/tau;n++) {
t = t+tau;
f0 = f(n⋆tau); f1 = f((n+1)⋆tau);
aTau;
oldU = u;
plot(u);
Vh uex = t⋆x^2⋆(1-x)^2; // exact
sol.= tx2(1 - x)2
Vh err = u - uex; // err =FE-sol -
exact
out<< abs(err[].max)/abs(uex[].max) <<","; //
∥err∥L∞(Ω)∕∥uex∥L∞(Ω)
}
out << endl;
theta = theta + 0.1;
}
We can see in Fig. 9.19 that uhn(θ) become unstable at θ = 0.4, and figures are omitted in the case
θ < 0.4.
9.5.2 Convection
The hyperbolic equation
appear frequently in scientific problems, for example, Navier-Stokes equation, Convection-Diffusion equation,
etc.
In the case of 1-dimensional space, we can easily find the general solution (x,t)
u(x,t) = u0(x -αt) of the
following equation, if α is constant,
because ∂tu+ α∂xu = -α
0 + a
0 = 0, where
0 = du0(x)∕dx. Even if α is not constant construction, the principle is
similar. One begins the ordinary differential equation (with convention which α is prolonged by zero apart from
(0,L) ×(0,T)):
In this equation τ is the variable and x,t is parameters, and we denote the solution by Xx,t(τ). Then it is noticed that
(x,t) → v(X(τ),τ) in τ = t satisfy the equation
and by the definition ∂tX = Ẋ = +α and ∂xX = ∂xx in τ = t, because if τ = t we have X(τ) = x. The general solution
of (9.35) is thus the value of the boundary condition in Xx,t(0), it is has to say u(x,t) = u0(Xx,t(0)) if Xx,t(0) is on the
x axis, u(x,t) = u0(Xx,t(0)) if Xx,t(0) is on the axis of t.
In higher dimension Ω ⊂ Rd, d = 2,3, the equation of the convection is written
where
(x,t)
Rd. FreeFem++ implements the Characteristic-Galerkin method for convection operators. Recall
that the equation (9.34) can be discretized as
where D is the total derivative operator. So a good scheme is one step of backward convection by the method of
Characteristics-Galerkin
where Xm(x) is an approximation of the solution at t = mτ of the ordinary differential equation
where
m(x) = (α1(x,mτ),α2(x,mτ)). Because, by Taylor’s expansion, we have
where Xi(t) are the i-th component of
(t), um(x) = u(x,mτ) and we used the chain rule and x =
((m + 1)τ). From
(9.37), it follows that Also we apply Taylor’s expansion for t
um(x -
m(x)t),0 ≤ t ≤ τ, then
Putting
we can get the approximation
A classical convection problem is that of the “rotating bell” (quoted from [14][p.16]). Let Ω be the unit disk
centered at 0, with its center rotating with speed α1 = y,α2 = -x We consider the problem (9.34) with f = 0 and the
initial condition u(x,0) = u0(x), that is, from (9.36)
The exact solution is u(x,t) = u(
(t)) where
equals x rotated around the origin by an angle θ = -t
(rotate in clockwise). So, if u0 in a 3D perspective looks like a bell, then u will have exactly the same
shape, but rotated by the same amount. The program consists in solving the equation until T = 2π, that
is for a full revolution and to compare the final solution with the initial one; they should be equal.
Example 9.19 (convect.edp)
border C(t=0, 2⋆pi) { x=cos(t); y=sin(t); }; // the unit
circle
mesh Th = buildmesh(C(70)); // triangulates the
disk
fespace Vh(Th,P1);
Vh u0 = exp(-10⋆((x-0.3)^2 +(y-0.3)^2)); // give
u0
real dt = 0.17,t=0; // time
step
Vh a1 = -y, a2 = x; // rotation
velocity
Vh u; //
um+1
for (int m=0; m<2⋆pi/dt ; m++) {
t += dt;
u=convect([a1,a2],-dt,u0); //
um+1 = um(Xm(x))
u0=u; //
m++
plot(u,cmm=" t="+t + ", min=" + u[].min + ", max=" + u[].max,wait=0);
};
Note 9.3 The scheme convectis unconditionally stable, then the bell become lower and lower (the maximum of
u37 is 0.406 as shown in Fig. 9.21).
9.5.3 2D Black-Scholes equation for an European Put option
In mathematical finance, an option on two assets is modeled by a Black-Scholes equations in two space variables,
(see for example Wilmott et al[39] or Achdou et al [3]).
which is to be integrated in
×
+ ×
+ subject to, in the case of a put Boundary conditions for this problem may not be so easy to device. As in the one dimensional case the PDE
contains boundary conditions on the axis x1 = 0 and on the axis x2 = 0, namely two one dimensional Black-Scholes
equations driven respectively by the data u
and u
. These will be automatically
accounted for because they are embedded in the PDE. So if we do nothing in the variational form
(i.e. if we take a Neumann boundary condition at these two axis in the strong form) there will be no
disturbance to these. At infinity in one of the variable, as in 1D, it makes sense to impose u = 0. We take
An implicit Euler scheme is used and a mesh adaptation is done every 10 time steps. To have an unconditionally
stable scheme, the first order terms are treated by the Characteristic Galerkin method, which, roughly, approximates
Example 9.20 [BlackSchol.edp]
// file
BlackScholes2D.edp
int m=30,L=80,LL=80, j=100;
real sigx=0.3, sigy=0.3, rho=0.3, r=0.05, K=40, dt=0.01;
mesh th=square(m,m,[L⋆x,LL⋆y]);
fespace Vh(th,P1);
Vh u=max(K-max(x,y),0.);
Vh xveloc, yveloc, v,uold;
for (int n=0; n⋆dt <= 1.0; n++)
{
if(j>20) { th = adaptmesh(th,u,verbosity=1,abserror=1,nbjacoby=2,
err=0.001, nbvx=5000, omega=1.8, ratio=1.8, nbsmooth=3,
splitpbedge=1, maxsubdiv=5,rescaling=1) ;
j=0;
xveloc = -x⋆r+x⋆sigx^2+x⋆rho⋆sigx⋆sigy/2;
yveloc = -y⋆r+y⋆sigy^2+y⋆rho⋆sigx⋆sigy/2;
u=u;
};
uold=u;
solve eq1(u,v,init=j,solver=LU) = int2d(th)( u⋆v⋆(r+1/dt)
+ dx(u)⋆dx(v)⋆(x⋆sigx)^2/2 + dy(u)⋆dy(v)⋆(y⋆sigy)^2/2
+ (dy(u)⋆dx(v) + dx(u)⋆dy(v))⋆rho⋆sigx⋆sigy⋆x⋆y/2)
- int2d(th)( v⋆convect([xveloc,yveloc],dt,w)/dt) + on(2,3,u=0);
j=j+1;
};
plot(u,wait=1,value=1);
Results are shown on Fig. 9.20).
9.6 Navier-Stokes Equation
9.6.1 Stokes and Navier-Stokes
The Stokes equations are: for a given
L2(Ω)2,
 | (9.43) |
where
= (u1,u2) is the velocity vector and p the pressure. For simplicity, let us choose Dirichlet boundary
conditions on the velocity,
=
Γ on Γ.
In Temam [Theorem 2.2], there ia a weak form of (9.43): Find
= (v1,v2)
(Ω)
which satisfy
Here it is used the existence p
H1(Ω) such that
= ∇p, if
Another weak form is derived as follows: We put
By multiplying the first equation in (9.43) with v
V and the second with q
W, subsequent integration over Ω, and
an application of Green’s formula, we have This yields the weak form of (9.43): Find (
,p)
×W such that for all (
,q)
V ×W, where
Now, we consider finite element spaces
h ⊂
and Wh ⊂ W, and we assume the following basis functions
The discrete weak form is: Find (
h,ph)
h ×Wh such that
 | (9.48) |
Note 9.4 Assume that:
- There is a constant αh > 0 such that
where
- There is a constant βh > 0 such that
Then we have an unique solution (
h,ph) of (9.48) satisfying
with a constant C > 0 (see e.g. [20, Theorem 10.4]).
Let us denote that
then (9.48) is written by where
Penalty method: This method consists of replacing (9.48) by a more regular problem: Find (
hϵ,phϵ)
h ×
h
satisfying
 | (9.51) |
where
h ⊂ L2(Ω). Formally, we have
and the corresponding algebraic problem
Note 9.5 We can eliminate phϵ = (1∕ϵ)BUhϵ to obtain
Since the matrix A + (1∕ϵ)B*B is symmetric, positive-definite, and sparse, (9.52) can be solved by known technique.
There is a constant C > 0 independent of ϵ such that
(see e.g. [20, 17.2])
Example 9.21 (Cavity.edp)
The driven cavity flow problem is solved first at zero Reynolds number (Stokes flow) and then at Reynolds
100. The velocity pressure formulation is used first and then the calculation is repeated with the stream
function vorticity formulation.
We solve the driven cavity problem by the penalty method (9.51) where
Γ ⋅
= 0 and
Γ ⋅
= 1 on the top
boundary and zero elsewhere (
is the unit normal to Γ, and
the unit tangent to Γ).
The mesh is constructed by
mesh Th=square(8,8);
We use a classical Taylor-Hood element technic to solve the problem:
The velocity is approximated with the P2 FE ( Xh space), and the the pressure is approximated with the P1
FE ( Mh space),
where
and
The FE spaces and functions are constructed by
fespace Xh(Th,P2); // definition of the velocity component
space
fespace Mh(Th,P1); // definition of the pressure
space
Xh u2,v2;
Xh u1,v1;
Mh p,q;
The Stokes operator is implemented as a system-solve for the velocity (u1,u2) and the pressure p. The test
function for the velocity is (v1,v2) and q for the pressure, so the variational form (9.48) in freefem language
is:
solve Stokes (u1,u2,p,v1,v2,q,solver=Crout) =
int2d(Th)( ( dx(u1)⋆dx(v1) + dy(u1)⋆dy(v1)
+ dx(u2)⋆dx(v2) + dy(u2)⋆dy(v2) )
- p⋆q⋆(0.000001)
- p⋆dx(v1) - p⋆dy(v2)
- dx(u1)⋆q - dy(u2)⋆q
)
+ on(3,u1=1,u2=0)
+ on(1,2,4,u1=0,u2=0); // see Section 5.1.1 for labels
1,2,3,4
Each unknown has its own boundary conditions.
If the streamlines are required, they can be computed by finding ψ such that rotψ = u or better,
Xh psi,phi;
solve streamlines(psi,phi) =
int2d(Th)( dx(psi)⋆dx(phi) + dy(psi)⋆dy(phi))
+ int2d(Th)( -phi⋆(dy(u1)-dx(u2)))
+ on(1,2,3,4,psi=0);
Now the Navier-Stokes equations are solved
with the same boundary conditions and with initial conditions u = 0.
This is implemented by using the convection operator convectfor the term ∂u
∂t + u ⋅∇u, giving a discretization in
time
 | (9.53) |
The term un ∘Xn(x) ≈ un(x -un(x)τ) will be computed by the operator “convect” , so we obtain
int i=0;
real nu=1./100.;
real dt=0.1;
real alpha=1/dt;
Xh up1,up2;
problem NS (u1,u2,p,v1,v2,q,solver=Crout,init=i) =
int2d(Th)(
alpha⋆( u1⋆v1 + u2⋆v2)
+ nu ⋆ ( dx(u1)⋆dx(v1) + dy(u1)⋆dy(v1)
+ dx(u2)⋆dx(v2) + dy(u2)⋆dy(v2) )
- p⋆q⋆(0.000001)
- p⋆dx(v1) - p⋆dy(v2)
- dx(u1)⋆q - dy(u2)⋆q
)
+ int2d(Th) ( -alpha⋆
convect([up1,up2],-dt,up1)⋆v1 -alpha⋆convect([up1,up2],-dt,up2)⋆v2 )
+ on(3,u1=1,u2=0)
+ on(1,2,4,u1=0,u2=0)
;
for (i=0;i<=10;i++)
{
up1=u1;
up2=u2;
NS;
if ( !(i % 10)) // plot every 10
iteration
plot(coef=0.2,cmm=" [u1,u2] and p ",p,[u1,u2]);
} ;
Notice that the stiffness matrices are reused (keyword init=i)
9.6.2 Uzawa Conjugate Gradient
We solve Stokes problem without penalty. The classical iterative method of Uzawa is described by the algorithm
(see e.g.[20, 17.3], [29, 13] or [30, 13] ):
-
Initialize:
- Let ph0 be an arbitrary chosen element of L2(Ω).
-
Calculate
h: - Once phn is known,
hn is the solution of
-
Advance ph:
- Let phn+1 be defined by
There is a constant α > 0 such that α ≤ ρn ≤ 2 for each n, then
hn converges to the solution
h, and then B
hn → 0
as n → ∞from the Advance ph. This method in general converges quite slowly.
First we define mesh, and the Taylor-Hood approximation. So Xh is the velocity space, and Mh is the pressure space.
Example 9.22 (StokesUzawa.edp)
mesh Th=square(10,10);
fespace Xh(Th,P2),Mh(Th,P1);
Xh u1,u2,v1,v2;
Mh p,q,ppp; // ppp is a working
pressure
varf bx(u1,q) = int2d(Th)( -(dx(u1)⋆q));
varf by(u1,q) = int2d(Th)( -(dy(u1)⋆q));
varf a(u1,u2)= int2d(Th)( dx(u1)⋆dx(u2) + dy(u1)⋆dy(u2) )
+ on(3,u1=1) + on(1,2,4,u1=0) ;
// remark: put the on(3,u1=1) before
on(1,2,4,u1=0)
// because we want zero on intersection
%
matrix A= a(Xh,Xh,solver=CG);
matrix Bx= bx(Xh,Mh); //
= (Bx By)
matrix By= by(Xh,Mh);
Xh bc1; bc1[] = a(0,Xh); // boundary condition contribution on
u1
Xh bc2; bc2 = O ; // no boundary condition contribution on
u2
Xh b;
phn →
A-1(-
*phn) = -div
h is realized as the function divup.
func real[int] divup(real[int] & pp)
{
// compute
u1(pp)
b[] = Bx'⋆pp; b[] ⋆=-1; b[] += bc1[] ; u1[] = A^-1⋆b[];
// compute
u2(pp)
b[] = By'⋆pp; b[] ⋆=-1; b[] += bc2[] ; u2[] = A^-1⋆b[];
//
n = A-1(BxTpn ByTpn)T
ppp[] = Bx⋆u1[]; //
ppp = Bxu1
ppp[] += By⋆u2[]; //
+Byu2
return ppp[] ;
};
Call now the conjugate gradient algorithm:
p=0;q=0; //
p0h = 0
LinearCG(divup,p[],eps=1.e-6,nbiter=50); //
pn+1h = pnh + 
n h
// if n > 50 or |pn+1h - pnh| ≤ 10-6, then the loop end.
divup(p[]); // compute the final
solution
plot([u1,u2],p,wait=1,value=true,coef=0.1);
9.6.3 NSUzawaCahouetChabart.edp
In this example we solve the Navier-Stokes equation, in the driven-cavity, with the Uzawa algorithm preconditioned
by the Cahouet-Chabart method (see [31] for all the details).
The idea of the preconditioner is that in a periodic domain, all differential operators commute and the Uzawa
algorithm comes to solving the linear operator ∇.((αId + νΔ)-1∇, where Id is the identity operator. So the
preconditioer suggested is αΔ-1 + νId.
To implement this, we reuse the previous example, by including a file. Then we define the time step Δt, viscosity,
and new variational form and matrix.
Example 9.23 (NSUzawaCahouetChabart.edp)
include "StokesUzawa.edp" // include the Stokes
part
real dt=0.05, alpha=1/dt; //
Δt
cout << " alpha = " << alpha;
real xnu=1./400; // viscosity
ν = Reynolds number-1
// the new variational form with mass term
varf at(u1,u2)= int2d(Th)( xnu⋆dx(u1)⋆dx(u2)
+ xnu⋆dy(u1)⋆dy(u2) + u1⋆u2⋆alpha )
+ on(1,2,4,u1=0) + on(3,u1=1) ;
A = at(Xh,Xh,solver=CG); // change the matrix
// set the 2 convect variational form
varf vfconv1(uu,vv) = int2d(Th,qforder=5) (convect([u1,u2],-dt,u1)⋆vv⋆alpha);
varf vfconv2(v2,v1) = int2d(Th,qforder=5) (convect([u1,u2],-dt,u2)⋆v1⋆alpha);
int idt; // index of time
set
real temps=0; // current
time
Mh pprec,prhs;
varf vfMass(p,q) = int2d(Th)(p⋆q);
matrix MassMh=vfMass(Mh,Mh,solver=CG);
varf vfLap(p,q) = int2d(Th)(dx(pprec)⋆dx(q)+dy(pprec)⋆dy(q) + pprec⋆q⋆1e-10);
matrix LapMh= vfLap(Mh,Mh,solver=Cholesky);
The function to define the preconditioner
func real[int] CahouetChabart(real[int] & xx)
{ // xx =
∫
(divu)wi
//
αLapMh-1 + νMassMh-1
pprec[]= LapMh^-1⋆ xx;
prhs[] = MassMh^-1⋆xx;
pprec[] = alpha⋆pprec[]+xnu⋆ prhs[];
return pprec[];
};
The loop in time. Warning with the stop test of the conjugate gradient, because we start from the previous
solution and the end the previous solution is close to the final solution, don’t take a relative stop test to the
first residual, take an absolute stop test ( negative here)
for (idt = 1; idt < 50; idt++)
{
temps += dt;
cout << " --------- temps " << temps << " \n ";
b1[] = vfconv1(0,Xh);
b2[] = vfconv2(0,Xh);
cout << " min b1 b2 " << b1[].min << " " << b2[].min << endl;
cout << " max b1 b2 " << b1[].max << " " << b2[].max << endl;
// call Conjugate Gradient with preconditioner
'
// warning eps < 0 => absolue stop test
LinearCG(divup,p[],eps=-1.e-6,nbiter=50,precon=CahouetChabart);
divup(p[]); // computed the
velocity
plot([u1,u2],p,wait=!(idt%10),value= 1,coef=0.1);
}
9.7 Variational inequality
We present, a classical examples of variational inequality.
Let us denote
= {u
H01(Ω),u ≤ g}
The problem is :
where
f and g are given function.
The solution is a projection on the convex
of f⋆ for the scalar product ((v,w)) = ∫
Ω∇v.∇w of H01(Ω) where f⋆ is
solution of ((f⋆,v)) = ∫
Ωfv,∀v
H01(Ω). The projection on a convex satisfy clearly ∀v
, ((u-v,u-
)) ≤ 0,
and after expanding, we get the classical inequality
We can also rewrite the problem as a saddle point problem
Find λ,u such that:
where
((u -g)+ = max(0,u -g)
This saddle point problem is equivalent to find u,λ such that:
 | (9.54) |
A algorithm to solve the previous problem is:
- k=0, and choose, λ0 belong H-1(Ω)
- loop on k = 0,.....
- set
k = {x
Ω∕λk + c *(uk+1 -g) ≤ 0}
- Vg,k+1 = {v
H01(Ω)∕v = g on Ik},
- V0,k+1 = {v
H01(Ω)∕v = 0 on Ik},
- Find uk+1
Vg,k+1 and λk+1
H-1(Ω) such that
where <,> is the duality bracket between H01(Ω) and H-1(Ω), and c is a penalty constant
(large enough).
You can find all the mathematic about this algorithm in [33].
Now how to do that in FreeFem++
The full example is:
Example 9.24 (VI.edp)
mesh Th=square(20,20);
real eps=1e-5;
fespace Vh(Th,P1); // P1 FE
space
int n = Vh.ndof; // number of Degree of
freedom
Vh uh,uhp; // solution and previous
one
Vh Ik; // to def the set where the containt is
reached.
real[int] rhs(n); // to store the right and side of the
equation
real c=1000; // the penalty parameter of the
algoritm
func f=1; // right hand side
function
func fd=0; // Dirichlet boundary condition
function
Vh g=0.05; // the discret function
g
real[int] Aii(n),Aiin(n); // to store the diagonal of the matrix 2
version
real tgv = 1e30; // a huge value for exact
penalization
// of boundary
condition
// the variatonal form of the problem:
varf a(uh,vh) = // definition of the
problem
int2d(Th)( dx(uh)⋆dx(vh) + dy(uh)⋆dy(vh) ) // bilinear
form
- int2d(Th)( f⋆vh ) // linear
form
+ on(1,2,3,4,uh=fd) ; // boundary condition
form
// two version of the matrix of the problem
matrix A=a(Vh,Vh,tgv=tgv,solver=CG); // one
changing
matrix AA=a(Vh,Vh,solver:GC); // one for computing
residual
// the mass Matrix construction:
varf vM(uh,vh) = int2d(Th)(uh⋆vh);
matrix M=vM(Vh,Vh); // to do a fast computing of L2 norm : sqrt(
u'⋆(w=M⋆u))
Aii=A.diag; // get the diagonal of the matrix (appear in version
1.46-1)
rhs = a(0,Vh,tgv=tgv);
Ik =0;
uhp=-tgv; // previous value
is
Vh lambda=0;
for(int iter=0;iter<100;++iter)
{
real[int] b(n) ; b=rhs; // get a copy of the Right hand
side
real[int] Ak(n); // the complementary of Ik ( !Ik =
(Ik-1))
// Today the operator Ik- 1. is not implement so we
do:
Ak= 1.; Ak -= Ik[]; // build Ak = !
Ik
// adding new locking condition on b and on the diagonal if (Ik ==1
)
b = Ik[] .⋆ g[]; b ⋆= tgv; b -= Ak .⋆ rhs;
Aiin = Ik[] ⋆ tgv; Aiin += Ak .⋆ Aii; // set Aii= tgv
i
Ik
A.diag = Aiin; // set the matrix diagonal (appear in version
1.46-1)
set(A,solver=CG); // important to change preconditioning for
solving
uh[] = A^-1⋆ b; // solve the problem with more locking
condition
lambda[] = AA ⋆ uh[]; // compute the residual ( fast with
matrix)
lambda[] += rhs; // remark rhs =
-∫
fv
Ik = ( lambda + c⋆( g- uh)) < 0.; // the new of locking
value
plot(Ik, wait=1,cmm=" lock set ",value=1,ps="VI-lock.eps",fill=1 );
plot(uh,wait=1,cmm="uh",ps="VI-uh.eps");
// trick to compute L2 norm of the variation (fast
method)
real[int] diff(n),Mdiff(n);
diff= uh[]-uhp[];
Mdiff = M⋆diff;
real err = sqrt(Mdiff'⋆diff);
cout << " || u_{k=1} - u_{k} ||_2 " << err << endl;
if(err< eps) break; // stop
test
uhp[]=uh[] ; // set the previous
solution
}
savemesh(Th,"mm",[x,y,uh⋆10]); // for medit
plotting
Remark, as you can see on this example, some vector , or matrix operator are not implemented so a way is to skip
the expression and we use operator +=, -=to merge the result.
9.8 Domain decomposition
We present, three classic examples, of domain decomposition technique: first, Schwarz algorithm with overlapping,
second Schwarz algorithm without overlapping (also call Shur complement), and last we show to use the conjugate
gradient to solve the boundary problem of the Shur complement.
9.8.1 Schwarz Overlap Scheme
To solve
the Schwarz algorithm runs like this where Γi is the boundary of Ωi and on the condition that Ω1 ∩Ω2 ⇔ ∅and that ui are zero at iteration 1.
Here we take Ω1 to be a quadrangle, Ω2 a disk and we apply the algorithm starting from zero.
Example 9.25 (Schwarz-overlap.edp)
int inside = 2; // inside
boundary
int outside = 1; // outside
boundary
border a(t=1,2){x=t;y=0;label=outside;};
border b(t=0,1){x=2;y=t;label=outside;};
border c(t=2,0){x=t ;y=1;label=outside;};
border d(t=1,0){x = 1-t; y = t;label=inside;};
border e(t=0, pi/2){ x= cos(t); y = sin(t);label=inside;};
border e1(t=pi/2, 2⋆pi){ x= cos(t); y = sin(t);label=outside;};
int n=4;
mesh th = buildmesh( a(5⋆n) + b(5⋆n) + c(10⋆n) + d(5⋆n));
mesh TH = buildmesh( e(5⋆n) + e1(25⋆n) );
plot(th,TH,wait=1); // to see the 2
meshes
The space and problem definition is :
fespace vh(th,P1);
fespace VH(TH,P1);
vh u=0,v; VH U,V;
int i=0;
problem PB(U,V,init=i,solver=Cholesky) =
int2d(TH)( dx(U)⋆dx(V)+dy(U)⋆dy(V) )
+ int2d(TH)( -V) + on(inside,U = u) + on(outside,U= 0 ) ;
problem pb(u,v,init=i,solver=Cholesky) =
int2d(th)( dx(u)⋆dx(v)+dy(u)⋆dy(v) )
+ int2d(th)( -v) + on(inside ,u = U) + on(outside,u = 0 ) ;
The calculation loop:
for ( i=0 ;i< 10; i++)
{
PB;
pb;
plot(U,u,wait=true);
};
9.8.2 Schwarz non Overlap Scheme
To solve
the Schwarz algorithm for domain decomposition without overlapping runs like this
Let introduce Γi is common the boundary of Ω1 and Ω2 and Γei = ∂Ωi \Γi.
The problem find λ such that (u1|Γi = u2|Γi) where ui is solution of the following Laplace problem:
To solve this problem we just make a loop with upgradingλ with
where
the sign + or -of ±is choose to have convergence.
Example 9.26 (Schwarz-no-overlap.edp)
// schwarz1 without
overlapping
int inside = 2;
int outside = 1;
border a(t=1,2){x=t;y=0;label=outside;};
border b(t=0,1){x=2;y=t;label=outside;};
border c(t=2,0){x=t ;y=1;label=outside;};
border d(t=1,0){x = 1-t; y = t;label=inside;};
border e(t=0, 1){ x= 1-t; y = t;label=inside;};
border e1(t=pi/2, 2⋆pi){ x= cos(t); y = sin(t);label=outside;};
int n=4;
mesh th = buildmesh( a(5⋆n) + b(5⋆n) + c(10⋆n) + d(5⋆n));
mesh TH = buildmesh ( e(5⋆n) + e1(25⋆n) );
plot(th,TH,wait=1,ps="schwarz-no-u.eps");
fespace vh(th,P1);
fespace VH(TH,P1);
vh u=0,v; VH U,V;
vh lambda=0;
int i=0;
problem PB(U,V,init=i,solver=Cholesky) =
int2d(TH)( dx(U)⋆dx(V)+dy(U)⋆dy(V) )
+ int2d(TH)( -V)
+ int1d(TH,inside)(-lambda⋆V) + on(outside,U= 0 ) ;
problem pb(u,v,init=i,solver=Cholesky) =
int2d(th)( dx(u)⋆dx(v)+dy(u)⋆dy(v) )
+ int2d(th)( -v)
+ int1d(th,inside)(+lambda⋆v) + on(outside,u = 0 ) ;
for ( i=0 ;i< 10; i++)
{
PB;
pb;
lambda = lambda - (u-U)/2;
plot(U,u,wait=true);
};
plot(U,u,ps="schwarz-no-u.eps");
9.8.3 Schwarz-gc.edp
To solve
the Schwarz algorithm for domain decomposition without overlapping runs like this
Let introduce Γi is common the boundary of Ω1 and Ω2 and Γei = ∂Ωi \Γi.
The problem find λ such that (u1|Γi = u2|Γi) where ui is solution of the following Laplace problem:
The version of this example for Shur componant. The border problem is solve with conjugate gradient.
First, we construct the two domain
Example 9.27 (Schwarz-gc.edp)
// Schwarz without overlapping (Shur complenement Neumann ->
Dirichet)
real cpu=clock();
int inside = 2;
int outside = 1;
border Gamma1(t=1,2){x=t;y=0;label=outside;};
border Gamma2(t=0,1){x=2;y=t;label=outside;};
border Gamma3(t=2,0){x=t ;y=1;label=outside;};
border GammaInside(t=1,0){x = 1-t; y = t;label=inside;};
border GammaArc(t=pi/2, 2⋆pi){ x= cos(t); y = sin(t);label=outside;};
int n=4;
// build the mesh of Ω1 and
Ω2
mesh Th1 = buildmesh( Gamma1(5⋆n) + Gamma2(5⋆n) + GammaInside(5⋆n) + Gamma3(5⋆n));
mesh Th2 = buildmesh ( GammaInside(-5⋆n) + GammaArc(25⋆n) );
plot(Th1,Th2);
// defined the 2 FE
space
fespace Vh1(Th1,P1), Vh2(Th2,P1);
Note 9.6 It is impossible to define a function just on a part of boundary, so the lambda function must be
defined on the all domain Ω1 such as
Vh1 lambda=0; // take
λ
Vh1
The two Poisson problem:
Vh1 u1,v1; Vh2 u2,v2;
int i=0; // for factorization
optimization
problem Pb2(u2,v2,init=i,solver=Cholesky) =
int2d(Th2)( dx(u2)⋆dx(v2)+dy(u2)⋆dy(v2) )
+ int2d(Th2)( -v2)
+ int1d(Th2,inside)(-lambda⋆v2) + on(outside,u2= 0 ) ;
problem Pb1(u1,v1,init=i,solver=Cholesky) =
int2d(Th1)( dx(u1)⋆dx(v1)+dy(u1)⋆dy(v1) )
+ int2d(Th1)( -v1)
+ int1d(Th1,inside)(+lambda⋆v1) + on(outside,u1 = 0 ) ;
or, we define a border matrix , because the lambda function is none zero inside the domain Ω1:
varf b(u2,v2,solver=CG) =int1d(Th1,inside)(u2⋆v2);
matrix B= b(Vh1,Vh1,solver=CG);
The boundary problem function,
func real[int] BoundaryProblem(real[int] &l)
{
lambda[]=l; // make FE function form
l
Pb1; Pb2;
i++; // no refactorization i
!=0
v1=-(u1-u2);
lambda[]=B⋆v1[];
return lambda[] ;
};
Note 9.7 The difference between the two notations v1and v1[]is: v1is the finite element function and
v1[]is the vector in the canonical basis of the finite element function v1.
Vh1 p=0,q=0;
// solve the problem with Conjugate
Gradient
LinearCG(BoundaryProblem,p[],eps=1.e-6,nbiter=100);
// compute the final solution, because CG works with
increment
BoundaryProblem(p[]); // solve again to have right
u1,u2
cout << " -- CPU time schwarz-gc:" << clock()-cpu << endl;
plot(u1,u2); //
plot
9.9 Fluid/Structures Coupled Problem
This problem involves the Lamé system of elasticity and the Stokes system for viscous fluids with velocity
and
pressure p:
where uΓ is the velocity of the boundaries. The force that the fluid applies to the boundaries is the normal stress
Elastic solids subject to forces deform: a point in the solid, at (x,y) goes to (X,Y) after. When the displacement
vector
= (v1,v2) = (X -x,Y -y) is small, Hooke’s law relates the stress tensor σ inside the solid to the
deformation tensor ϵ:
where
δ is the Kronecker symbol and where λ,μ are two constants describing the material mechanical properties in terms
of the modulus of elasticity, and Young’s modulus.
The equations of elasticity are naturally written in variational form for the displacement vector v(x)
V as
The
data are the gravity force
and the boundary stress
.
Example 9.28 (fluidStruct.edp) In our example the Lamé system and the Stokes system are coupled by a
common boundary on which the fluid stress creates a displacement of the boundary and hence changes the
shape of the domain where the Stokes problem is integrated. The geometry is that of a vertical driven cavity
with an elastic lid. The lid is a beam with weight so it will be deformed by its own weight and by the normal
stress due to the fluid reaction. The cavity is the 10 ×10 square and the lid is a rectangle of height l = 2.
A beam sits on a box full of fluid rotating because the left vertical side has velocity one. The beam is bent by
its own weight, but the pressure of the fluid modifies the bending.
The bending displacement of the beam is given by (uu,vv) whose solution is given as follows.
// Fluid-structure interaction for a weighting beam sitting on
a
// square cavity filled with a
fluid.
int bottombeam = 2; // label of
bottombeam
border a(t=2,0) { x=0; y=t ;label=1;}; // left
beam
border b(t=0,10) { x=t; y=0 ;label=bottombeam;}; // bottom of
beam
border c(t=0,2) { x=10; y=t ;label=1;}; // rigth
beam
border d(t=0,10) { x=10-t; y=2; label=3;}; // top
beam
real E = 21.5;
real sigma = 0.29;
real mu = E/(2⋆(1+sigma));
real lambda = E⋆sigma/((1+sigma)⋆(1-2⋆sigma));
real gravity = -0.05;
mesh th = buildmesh( b(20)+c(5)+d(20)+a(5));
fespace Vh(th,P1);
Vh uu,w,vv,s,fluidforce=0;
cout << "lambda,mu,gravity ="<<lambda<< " " << mu << " " << gravity << endl;
// deformation of a beam under its own
weight
solve bb([uu,vv],[w,s]) =
int2d(th)(
lambda⋆div(w,s)⋆div(uu,vv)
+2.⋆mu⋆( epsilon(w,s)'⋆epsilon(uu,vv) )
)
+ int2d(th) (-gravity⋆s)
+ on(1,uu=0,vv=0)
+ fluidforce[];
;
plot([uu,vv],wait=1);
mesh th1 = movemesh(th, [x+uu, y+vv]);
plot(th1,wait=1);
Then Stokes equation for fluids ast low speed are solved in the box below the beam, but the beam has
deformed the box (see border h):
// Stokes on square b,e,f,g driven cavite on left side
g
border e(t=0,10) { x=t; y=-10; label= 1; }; //
bottom
border f(t=0,10) { x=10; y=-10+t ; label= 1; }; //
right
border g(t=0,10) { x=0; y=-t ;label= 2;}; //
left
border h(t=0,10) { x=t; y=vv(t,0)⋆( t>=0.001 )⋆(t <= 9.999);
label=3;}; // top of cavity
deformed
mesh sh = buildmesh(h(-20)+f(10)+e(10)+g(10));
plot(sh,wait=1);
We use the Uzawa conjugate gradient to solve the Stokes problem like in example Section 9.6.2
fespace Xh(sh,P2),Mh(sh,P1);
Xh u1,u2,v1,v2;
Mh p,q,ppp;
varf bx(u1,q) = int2d(sh)( -(dx(u1)⋆q));
varf by(u1,q) = int2d(sh)( -(dy(u1)⋆q));
varf Lap(u1,u2)= int2d(sh)( dx(u1)⋆dx(u2) + dy(u1)⋆dy(u2) )
+ on(2,u1=1) + on(1,3,u1=0) ;
Xh bc1; bc1[] = Lap(0,Xh);
Xh brhs;
matrix A= Lap(Xh,Xh,solver=CG);
matrix Bx= bx(Xh,Mh);
matrix By= by(Xh,Mh);
Xh bcx=0,bcy=1;
func real[int] divup(real[int] & pp)
{
int verb=verbosity;
verbosity=0;
brhs[] = Bx'⋆pp; brhs[] += bc1[] .⋆bcx[];
u1[] = A^-1⋆brhs[];
brhs[] = By'⋆pp; brhs[] += bc1[] .⋆bcy[];
u2[] = A^-1⋆brhs[];
ppp[] = Bx⋆u1[];
ppp[] += By⋆u2[];
verbosity=verb;
return ppp[] ;
};
do a loop on the two problem
for(step=0;step<2;++step)
{
p=0;q=0;u1=0;v1=0;
LinearCG(divup,p[],eps=1.e-3,nbiter=50);
divup(p[]);
Now the beam will feel the stress constraint from the fluid:
Vh sigma11,sigma22,sigma12;
Vh uu1=uu,vv1=vv;
sigma11([x+uu,y+vv]) = (2⋆dx(u1)-p);
sigma22([x+uu,y+vv]) = (2⋆dy(u2)-p);
sigma12([x+uu,y+vv]) = (dx(u1)+dy(u2));
which comes as a boundary condition to the PDE of the beam:
solve bbst([uu,vv],[w,s],init=i) =
int2d(th)(
lambda⋆div(w,s)⋆div(uu,vv)
+2.⋆mu⋆( epsilon(w,s)'⋆epsilon(uu,vv) )
)
+ int2d(th) (-gravity⋆s)
+ int1d(th,bottombeam)( -coef⋆( sigma11⋆N.x⋆w + sigma22⋆N.y⋆s
+ sigma12⋆(N.y⋆w+N.x⋆s) ) )
+ on(1,uu=0,vv=0);
plot([uu,vv],wait=1);
real err = sqrt(int2d(th)( (uu-uu1)^2 + (vv-vv1)^2 ));
cout << " Erreur L2 = " << err << "----------\n";
Notice that the matrix generated by bbst is reused (see init=i). Finally we deform the beam
th1 = movemesh(th, [x+0.2⋆uu, y+0.2⋆vv]);
plot(th1,wait=1);
} // end of
loop
9.10 Transmission Problem
Consider an elastic plate whose displacement change vertically, which is made up of three plates of different
materials, welded on each other. Let Ωi,i = 1,2,3 be the domain occupied by i-th material with tension μi (see
Section 9.1.1). The computational domain Ω is the interior of Ω1 ∪Ω2 ∪Ω3. The vertical displacement u(x,y) is
obtained from
where ∂nu|Γi denotes the value of the normal derivative ∂nu on the boundary Γi of the domain Ωi.
By introducing the characteristic function χi of Ωi, that is,
 | (9.57) |
we can easily rewrite (9.55) and (9.56) to the weak form. Here we assume that u = 0 on Γ = ∂Ω.
problem Transmission: For a given function f, find u such that
where μ = μ1χ1 + μ2χ2 + μ3χ3. Here we notice that μ become the discontinuous function.
With dissipation, and at the thermal equilibrium, the temperature equation is:
This example explains the definition and manipulation of region, i.e. subdomains of the whole domain.
Consider this L-shaped domain with 3 diagonals as internal boundaries, defining 4 subdomains:
// example using region
keyword
// construct a mesh with 4 regions
(sub-domains)
border a(t=0,1){x=t;y=0;};
border b(t=0,0.5){x=1;y=t;};
border c(t=0,0.5){x=1-t;y=0.5;};
border d(t=0.5,1){x=0.5;y=t;};
border e(t=0.5,1){x=1-t;y=1;};
border f(t=0,1){x=0;y=1-t;};
// internal
boundary
border i1(t=0,0.5){x=t;y=1-t;};
border i2(t=0,0.5){x=t;y=t;};
border i3(t=0,0.5){x=1-t;y=t;};
mesh th = buildmesh (a(6) + b(4) + c(4) +d(4) + e(4) +
f(6)+i1(6)+i2(6)+i3(6));
fespace Ph(th,P0); // constant discontinuous functions /
element
fespace Vh(th,P1); // P1 continuous functions /
element
Ph reg=region; // defined the P0 function associated to region
number
plot(reg,fill=1,wait=1,value=1);
regionis a keyword of FreeFem++which is in fact a variable depending of the current position (is not a function
today, use Ph reg=region;to set a function). This variable value returned is the number of the subdomain of
the current position. This number is defined by ”buildmesh” which scans while building the mesh all
its connected component. So to get the number of a region containing a particular point one does:
int nupper=reg(0.4,0.9); // get the region number of point
(0.4,0.9)
int nlower=reg(0.9,0.1); // get the region number of point
(0.4,0.1)
cout << " nlower " << nlower << ", nupper = " << nupper<< endl;
// defined the characteristics functions of upper and lower
region
Ph nu=1+5⋆(region==nlower) + 10⋆(region==nupper);
plot(nu,fill=1,wait=1);
This is particularly useful to define discontinuous functions such as might occur when one part of the domain is
copper and the other one is iron, for example.
We this in mind we proceed to solve a Laplace equation with discontinuous coefficients (ν is 1, 6 and 11 below).
Ph nu=1+5⋆(region==nlower) + 10⋆(region==nupper);
plot(nu,fill=1,wait=1);
problem lap(u,v) = int2d(th)( nu⋆( dx(u)⋆dx(v)⋆dy(u)⋆dy(v) ))
+ int2d(-1⋆v) + on(a,b,c,d,e,f,u=0);
plot(u);
9.11 Free Boundary Problem
The domain Ω is defined with:
real L=10; // longueur du
domaine
real h=2.1; // hauteur du bord
gauche
real h1=0.35; // hauteur du bord
droite
// maillage d'un
tapeze
border a(t=0,L){x=t;y=0;}; // bottom: Γa
border b(t=0,h1){x=L;y=t;}; // right: Γb
border f(t=L,0){x=t;y=t⋆(h1-h)/L+h;}; // free surface: Γf
border d(t=h,0){x=0;y=t;}; // left: Γd
int n=4;
mesh Th=buildmesh (a(10⋆n)+b(6⋆n)+f(8⋆n)+d(3⋆n));
plot(Th,ps="dTh.eps");
The free boundary problem is:
Find u and Ω such that:
We use a fixed point method; Ω0 = Ω
in two step, fist we solve the classical following problem:
The variational formulation is:
find u on V = H1(Ωn), such than u = y on Γbn and Γf n
and secondly to construct a domain deformation
(x,y) = [x,y -v(x,y)]
where v is solution of the following problem:
The variational formulation is:
find v on V, such than v = 0 on Γan
Finally the new domain Ωn+1 =
(Ωn)
Example 9.29 (freeboundary.edp) The FreeFem++:implementation is:
real q=0.02; // flux
entrant
real K=0.5; //
permeabilité
fespace Vh(Th,P1);
int j=0;
Vh u,v,uu,vv;
problem Pu(u,uu,solver=CG) = int2d(Th)( dx(u)⋆dx(uu)+dy(u)⋆dy(uu))
+ on(b,f,u=y) ;
problem Pv(v,vv,solver=CG) = int2d(Th)( dx(v)⋆dx(vv)+dy(v)⋆dy(vv))
+ on (a, v=0) + int1d(Th,f)(vv⋆((q/K)⋆N.y- (dx(u)⋆N.x+dy(u)⋆N.y)));
real errv=1;
real erradap=0.001;
verbosity=1;
while(errv>1e-6)
{
j++;
Pu;
Pv;
plot(Th,u,v ,wait=0);
errv=int1d(Th,f)(v⋆v);
real coef=1;
//
real mintcc = checkmovemesh(Th,[x,y])/5.;
real mint = checkmovemesh(Th,[x,y-v⋆coef]);
if (mint<mintcc || j%10==0) { // mesh to bad =>
remeshing
Th=adaptmesh(Th,u,err=erradap ) ;
mintcc = checkmovemesh(Th,[x,y])/5.;
}
while (1)
{
real mint = checkmovemesh(Th,[x,y-v⋆coef]);
if (mint>mintcc) break;
cout << " min |T] " << mint << endl;
coef /= 1.5;
}
Th=movemesh(Th,[x,y-coef⋆v]); // calcul de la
deformation
cout << "\n\n"<<j <<"------------ errv = " << errv << "\n\n";
}
plot(Th,ps="d_Thf.eps");
plot(u,wait=1,ps="d_u.eps");
9.12 nolinear-elas.edp
The nonlinear elasticity problem is find the displacement (u1,u2) minimizing J
where
F2(u1,u2) = A(E[u1,u2],E[u1,u2]) and A(X,Y) is bilinear sym. positive form with respect two matrix X,Y. where f
is a given
2 function, and E[u1,u2] = (Eij)i=1,2,j=1,2 is the Green-Saint Venant deformation tensor defined
with:
denote u = (u1,u2), v = (v1,v2), w = (w1,w2).
So, the differential of J is
where
DF2(u)(v) = 2A(DE[u](v),E[u]) and DE is the first differential of E.
The second order differential is
where
and
D2E is the second differential of E.
So all notation can be define with macrolike:
macro EL(u,v) [dx(u),(dx(v)+dy(u)),dy(v)] // is
[ϵ11, 2ϵ12,ϵ22]
macro ENL(u,v) [
(dx(u)⋆dx(u)+dx(v)⋆dx(v))⋆0.5,
(dx(u)⋆dy(u)+dx(v)⋆dy(v)) ,
(dy(u)⋆dy(u)+dy(v)⋆dy(v))⋆0.5 ] // EOM
ENL
macro dENL(u,v,uu,vv) [(dx(u)⋆dx(uu)+dx(v)⋆dx(vv)),
(dx(u)⋆dy(uu)+dx(v)⋆dy(vv)+dx(uu)⋆dy(u)+dx(vv)⋆dy(v)),
(dy(u)⋆dy(uu)+dy(v)⋆dy(vv)) ] //
macro E(u,v) (EL(u,v)+ENL(u,v)) // is
[E11, 2E12,E22]
macro dE(u,v,uu,vv) (EL(uu,vv)+dENL(u,v,uu,vv)) //
macro ddE(u,v,uu,vv,uuu,vvv) dENL(uuu,vvv,uu,vv) //
macro F2(u,v) (E(u,v)'⋆A⋆E(u,v)) //
macro dF2(u,v,uu,vv) (E(u,v)'⋆A⋆dE(u,v,uu,vv)⋆2. ) //
macro ddF2(u,v,uu,vv,uuu,vvv) (
(dE(u,v,uu,vv)'⋆A⋆dE(u,v,uuu,vvv))⋆2.
+ (E(u,v)'⋆A⋆ddE(u,v,uu,vv,uuu,vvv))⋆2. ) //
EOM
The Newton Method is
choose n = 0,and uO,vO the initial displacement
- loop:
- find (du,dv) : solution of
- un = un -du, vn = vn -dv
- until (du,dv) small is enough
The way to implement this algorithm in FreeFem++ is use a macro tool to implement A and F2, f,
f′,f′′.
A macro is like in ccppreprocessor of C++ , but this begin by macroand the end of the macro definition is before
the comment ∕∕. In this case the macro is very useful because the type of parameter can be change. And it is easy to
make automatic differentiation.
// non linear elasticity
model
// for hyper elasticity
problem
//
-----------------------------
macro f(u) (u) // end of
macro
macro df(u) (1) // end of
macro
macro ddf(u) (0) // end of
macro
// -- du caouchouc --- CF cours de Herve Le
Dret.
//
-------------------------------
real mu = 0.012e5; //
kg∕cm2
real lambda = 0.4e5; //
kg∕cm2
//
// σ = 2μE + λtr(E)Id
// A(u,v) = σ(u) : E(v)
//
// ( a b )
// ( b c )
//
// tr⋆Id : (a,b,c) -> (a+c,0,a+c)
// so the associed matrix is:
// ( 1 0 1 )
// ( 0 0 0 )
// ( 1 0 1 )
//
------------------v
real a11= 2⋆mu + lambda ;
real a22= mu ; // because
[0, 2 * t12, 0]′A[0, 2 * s12, 0] =
//
= 2 * mu * (t12 * s12 + t21 * s21) = 4 * mu * t12 * s12
real a33= 2⋆mu + lambda ;
real a12= 0 ;
real a13= lambda ;
real a23= 0 ;
// symetric
part
real a21= a12 ;
real a31= a13 ;
real a32= a23 ;
// the matrix
A.
func A = [ [ a11,a12,a13],[ a21,a22,a23],[ a31,a32,a33] ];
real Pa=1e2; // a pressure of 100
Pa
// ----------------
int n=30,m=10;
mesh Th= square(n,m,[x,.3⋆y]); // label: 1 bottom, 2 right, 3 up, 4
left;
int bottom=1, right=2,upper=3,left=4;
plot(Th);
fespace Wh(Th,P1dc);
fespace Vh(Th,[P1,P1]);
fespace Sh(Th,P1);
Wh e2,fe2,dfe2,ddfe2; //
optimisation
Wh ett,ezz,err,erz; //
optimisation
Vh [uu,vv], [w,s],[un,vn];
[un,vn]=[0,0]; //
intialisation
[uu,vv]=[0,0];
varf vmass([uu,vv],[w,s],solver=CG) = int2d(Th)( uu⋆w + vv⋆s );
matrix M=vmass(Vh,Vh);
problem NonLin([uu,vv],[w,s],solver=LU)=
int2d(Th,qforder=1)( // (D2J(un))
part
dF2(un,vn,uu,vv)⋆dF2(un,vn,w,s)⋆ddfe2
+ ddF2(un,vn,w,s,uu,vv)⋆dfe2
)
- int1d(Th,3)(Pa⋆s)
- int2d(Th,qforder=1)( // (DJ(un))
part
dF2(un,vn,w,s)⋆dfe2 )
+ on(right,left,uu=0,vv=0);
;
// Newton's
method
//
---------------
Sh u1,v1;
for (int i=0;i<10;i++)
{
cout << "Loop " << i << endl;
e2 = F2(un,vn);
dfe2 = df(e2) ;
ddfe2 = ddf(e2);
cout << " e2 max " <<e2[].max << " , min" << e2[].min << endl;
cout << " de2 max "<< dfe2[].max << " , min" << dfe2[].min << endl;
cout << "dde2 max "<< ddfe2[].max << " , min" << ddfe2[].min << endl;
NonLin; // compute
[uu,vv] = (D2J(un))-1(DJ(un))
w[] = M⋆uu[];
real res = sqrt(w[]' ⋆ uu[]); // norme
L2of [uu,vv]
u1 = uu;
v1 = vv;
cout << " L^2 residual = " << res << endl;
cout << " u1 min =" <<u1[].min << ", u1 max= " << u1[].max << endl;
cout << " v1 min =" <<v1[].min << ", v2 max= " << v1[].max << endl;
plot([uu,vv],wait=1,cmm=" uu, vv " );
un[] -= uu[];
plot([un,vn],wait=1,cmm=" deplacement " );
if (res<1e-5) break;
}
plot([un,vn],wait=1);
mesh th1 = movemesh(Th, [x+un, y+vn]);
plot(th1,wait=1); // see figure
9.36
9.13 Compressible Neo-Hookean Materials: Computational Solutions
Author : Alex Sadovsky mailsashas@gmail.com
9.13.1 Notation
In what follows, the symbols
,F,B,C,σ denote, respectively, the displacement field, the deformation gradient, the
left Cauchy-Green strain tensor B = FFT, the right Cauchy-Green strain tensor C = FTF, and the
Cauchy stress tensor. We also introduce the symbols I1 := trC and J := det F. Use will be made of the
identity
 | (9.59) |
The symbol I denotes the identity tensor. The symbol Ω0 denotes the reference configuration of the body to be
deformed. The unit volume in the reference (resp., deformed) configuration is denoted dV (resp., dV0); these two are
related by
which
allows an integral over Ω involving the Cauchy stress T to be rewritten as an integral of the Kirchhoff stress κ = JT
over Ω0.
Recommended References
For an exposition of nonlinear elasticity and of the underlying linear- and tensor algebra, see [34].
For an advanced mathematical analysis of the Finite Element Method, see [35]. An explanation
of the Finite Element formulation of a nonlinear elastostatic boundary value problem, see
http://www.engin.brown.edu/courses/en222/Notes/FEMfinitestrain/FEMfinitestrain.htm.
9.13.2 A Neo-Hookean Compressible Material
Constitutive Theory and Tangent Stress Measures
The strain energy density function is given by
 | (9.60) |
(see [32], formula (12)).
The corresponding 2nd Piola-Kirchoff stress tensor is given by
 | (9.61) |
The Kirchhoff stress, then, is
 | (9.62) |
The tangent Kirchhoff stress tensor at Fn acting on δFn+1 is, consequently,
![∂-κ(F )δF = μ [F (δF )T + δF (F )T]
∂F n n+1 n n+1 n+1 n](freefem++doc714x.png) | (9.63) |
The Weak Form of the BVP in the Absence of Body (External) Forces
The Ω0 we are considering is an elliptical annulus, whose boundary consists of two concentric ellipses (each
allowed to be a circle as a special case), with the major axes parallel. Let P denote the dead stress load (traction) on
a portion ∂Ω0t (= the inner ellipse) of the boundary ∂Ω0. On the rest of the boundary, we prescribe zero
displacement.
The weak formulation of the boundary value problem is
For
brevity, in the rest of this section we assume P = 0. The provided FreeFem++ code, however, does not rely on this
assumption and allows for a general value and direction of P.
Given a Newton approximation
n of the displacement field
satisfying the BVP, we seek the correction δ
n+1 to
obtain a better approximation
by
solving the weak formulation
![∫ { } ∫ )
0 = Ω0 κ[Fn + δFn+1] : (r ⊗ w)(Fn + δFn+1)-1 - ∂Ω 0 P⋅N^0 ||||
= ∫ {κ[Fn ]+ ∂κ[Fn]δFn +1} : {(r ⊗ w)(Fn + δFn+1)-1} ||||
∫Ω0{ ∂∂Fκ } { -1 -2 } |||||
= Ω0 κ[Fn ]+ ∂F[Fn]δFn +1 : (r ⊗ w)(Fn + Fn δFn+1) |||||
∫ { } ||}
= κ[Fn] : (r ⊗ w)F-n1 |||| for all test functions w ,
∫Ω0 { -2 } |||||
- ∫Ω0 κ{[Fn] : (r ⊗ w})(F{n δFn+1) } |||||
+ Ω0 ∂∂κF[Fn]δFn+1 : (r ⊗ w)Fn-1 |||||
|)](freefem++doc720x.png) | (9.64) |
where we have taken
Note: Contrary to standard notational use, the symbol δ here bears no variational context. By δ we mean simply an
increment in the sense of Newton’s Method. The role of a variational virtual displacement here is played by
.
9.13.3 An Approach to Implementation in FreeFem++
The associated file is examples++-tutorial/nl-elast-neo-Hookean.edp.
Introducing the code-like notation, where a string in <>’s is to be read as one symbol, the individual components of
the tensor
![< TanK >:= ∂κ[F ]δF
∂F n n+1](freefem++doc723x.png) | (9.65) |
will be implemented as the macros < TanK11 >,< TanK12 >,….
The individual components of the tensor quantities
and
will
be implemented as the macros
 | (9.66) |
respectively.
In the above notation, the tangent Kirchhoff stress term becomes
 | (9.67) |
while the weak BVP formulation acquires the form
![∫ )
0 = Ω∫ 0 κ[Fn ] : D 4 |||||
- Ω∫ 0 κ{[Fn ] : D 3 } ||}
+ ∂κ[Fn ]δFn+1 : D 4|||| for all test functions w
Ω 0 ∂F |||)](freefem++doc730x.png) | (9.68) |