Chapter 5
Mesh Generation

5.1 Commands for Mesh Generation

border, buildmesh are explained.

All the examples in this section come from the files mesh.edp and tablefunction.edp.

5.1.1 Square

For easy and simple testing, there is the command “square”. The following


mesh Th = square(4,5);

generate a 4 ×5 grid in the unit square [0,1]2 whose labels are shown in Fig. 5.1.


PIC


Figure 5.1: Boundary labels of the mesh by square(10,10)


If you want to construct a n ×m grid in the rectangle [x0,x1] ×[y0,y1], you can write


  real x0=1.2,x1=1.8;
  real y0=0,y1=1;
  int n=5,m=20;
  mesh Th=square(n,m,[x0+(x1-x0)⋆x,y0+(y1-y0)⋆y]);

Note 5.1 Adding the parameter flags=1, will produce a Union Jack flag type of mesh.


  mesh Th=square(n,m,[x0+(x1-x0)⋆x,y0+(y1-y0)⋆y],flags=1);

5.1.2 Border

A domain is defined as being on the left (resp right) of its parameterized boundary

          |
Γ j = {(x,y)||x = φx (t),y = φy (t),aj ≤ t ≤ bj}
We can easily check the orientation by drawing the curve t↦→(φx(t)y(t)),t0 t t1. If it is as in Fig. 5.2, then the domain lie on the shaded area, otherwise it lies on the opposite side

The boundaries Γj can only intersect at their end points.


PIC


Figure 5.2: Orientation of the boundary defined by (ϕx(t)y(t))


The general expression to define a triangulation with buildmesh is

mesh   Mesh_Name = buildmesh(                                     )
 Γ 1(m1) + ⋅⋅⋅ + Γ J(mj ) O ptionalParam eter;

where mj are positive or negative numbers to indicate how may point should be put on Γj,Γ = j=1JΓ J, an the optional parameter (separed with comma) can be

nbvx=<int value>
, to set the maximal number of vertices in the mesh.
fixeborder=<bool value>
, to say if the mesh generator can change the boundary mesh or not (by default the boundary mesh can change and in case of periodic boundary problem (see. 6), it can be dangerous .

We can change the orientation of boundaries by changing the sign of mj. The following example shows how to change the orientation. The example generates the unit disk with a small circular hole, and assign “1” to the unit disk (“2” to the circle inside). The boundary label must be non-zero, but it can also be omitted.


1: border a(t=0,2⋆pi){ x=cos(t); y=sin(t);label=1;}
2: border b(t=0,2⋆pi){ x=0.3+0.3⋆cos(t); y=0.3⋆sin(t);label=2;}
3: plot(a(50)+b(+30)) ;  // to see a plot of the border mesh
4: mesh Thwithouthole= buildmesh(a(50)+b(+30));
5: mesh Thwithhole   = buildmesh(a(50)+b(-30));
6: plot(Thwithouthole,wait=1,ps="Thwithouthole.eps");  // figure 5.3
7: plot(Thwithhole,wait=1,ps="Thwithhole.eps");  // figure 5.4

Note 5.2 You must notice that the orientation is changed by “b(-30)” in 5th line. In 7th line, ps="fileName"is used to generate a postscript file identification that is shown on screen.


PIC

Figure 5.3: mesh without hole

PIC

Figure 5.4: mesh with hole


Note 5.3 The border is evaluated only at the time plotor buildmeshis called so the global variable are defined at this time and the following code can not work:


   real r=1;    border a(t=0,2⋆pi){ x=r⋆cos(t); y=r⋆sin(t);label=1;}
   r=0.3    ;   border b(t=0,2⋆pi){ x=r⋆cos(t); y=r⋆sin(t);label=1;}
   mesh Thwithhole   = buildmesh(a(50)+b(-30));  // bug (a trap) because
    // the two circle have the same radius = 0.3

5.1.3 Data Structure and Read/Write Statements for a Mesh

Users who want to use a triangulation made elsewhere should see the file structure of the file generated below:


border C(t=0,2⋆pi) { x=cos(t); y=sin(t); }
mesh Th = buildmesh(C(10));
savemesh("mesh_sample.msh");

the mesh is shown on Fig. 5.5.

The informations about Th are save in the file “mesh_sample.msh”. whose structure is shown on Table 5.1.
There nv denotes the number of vertices, nt number of triangles and ns the number of edges on boundary.

For each vertex qi,i = 1,⋅⋅⋅,n v, we denote by (qxi,q yi) the x-coordinate and y-coordinate.

Each triangle Tk,k = 1,⋅⋅⋅,10 has three vertices qk1,qk2,qk3 that are oriented counterclockwise. The boundary consists of 10 lines Li,i = 1,⋅⋅⋅,10 whose end points are qi1,qi2.


PIC

Figure 5.5: mesh by buildmesh(C(10))
In the left figure, we have the following.

nv = 14,nt = 16,ns = 10

q1 = (-0.309016994375,0.951056516295)
..
. ..
. ..
.
q14 = (-0.309016994375, -0.951056516295)

The vertices of T1 are q9,q12,q10.
.
.. .
.. .
..
The vertices of T16 are q9,q10,q6.

The edge of 1st side L1 are q6,q5.
... ... ...
The edge of 10th side L10 are q10,q6.





Contents of file Explanation


14 16 10 nv nt ne
-0.309016994375 0.951056516295 1 qx1 qy1 boundary label=1
0.309016994375 0.951056516295 1 qx2 qy2 boundary label=1
⋅⋅⋅⋅⋅⋅.
..
-0.309016994375 -0.951056516295 1qx14 qy14 boundary label=1


9 12 10 0 11 12 13 region label=0
5 9 6 0 21 22 23 region label=0
⋅⋅⋅
9 10 6 0 161 162 163 region label=0


6 5 1 11 12 boundary label=1
5 2 1 21 22 boundary label=1
⋅⋅⋅
10 6 1 101 102 boundary label=1



Table 5.1: The structure of “mesh_sample.msh”

There are many mesh file formats available for communication with other tools such as emc2, modulef.. (see Section 11), The extension of a file gives the chosen type. More details can be found in the article by F. Hecht ”bamg : a bidimentional anisotropic mesh generator” available from the FreeFem web site.


A mesh file can be read back into FreeFem++ but the names of the borders are lost. So these borders have to be referenced by the number which corresponds to their order of appearance in the program, unless this number is forced by the keyword ”label”. Here are some examples:


border floor(t=0,1){ x=t; y=0; label=1;};  // the unit square
border right(t=0,1){ x=1; y=t; label=5;};
border ceiling(t=1,0){ x=t; y=1; label=5;};
border left(t=1,0){ x=0; y=t; label=5;};
int n=10;
mesh th= buildmesh(floor(n)+right(n)+ceiling(n)+left(n));
savemesh(th,"toto.am_fmt");   // "formatted Marrocco" format
savemesh(th,"toto.Th");       // "bamg"-type mesh
savemesh(th,"toto.msh");      // freefem format
savemesh(th,"toto.nopo");      // modulef format see [10]
mesh th2 = readmesh("toto.msh");  // read the mesh

Example 5.1 (Readmesh.edp)


border floor(t=0,1){ x=t; y=0; label=1;};  // the unit square
border right(t=0,1){ x=1; y=t; label=5;};
border ceiling(t=1,0){ x=t; y=1; label=5;};
border left(t=1,0){ x=0; y=t; label=5;};
int n=10;
mesh th= buildmesh(floor(n)+right(n)+ceiling(n)+left(n));
savemesh(th,"toto.am_fmt");   // format "formated Marrocco"
savemesh(th,"toto.Th");       // format database db mesh "bamg"
savemesh(th,"toto.msh");      // format freefem
savemesh(th,"toto.nopo");      // modulef format see [10]
mesh th2 = readmesh("toto.msh");
fespace femp1(th,P1);
femp1 f = sin(x)⋆cos(y),g;
{  // save solution
ofstream file("f.txt");
file << f[] << endl;
}   // close the file (end block)
{   // read
ifstream file("f.txt");
file >> g[] ;
}  // close reading file (end block)
fespace Vh2(th2,P1);
Vh2 u,v;
plot(g);
// find u such that
// u + Δu = g in Ω ,
// u = 0 on Γ1 and ∂∂un = g on Γ2
solve pb(u,v) =
    int2d(th)( u⋆v - dx(u)⋆dx(v)-dy(u)⋆dy(v) )
  + int2d(th)(-g⋆v)
  + int1d(th,5)( g⋆v)  // ∂∂un = g on Γ2
  + on(1,u=0) ;
plot (th2,u);

5.1.4 Mesh Connectivity

The following example explains methods to obtain mesh information.


{  // get mesh information (version 1.37)
  mesh Th=square(2,2);
   // get data of the mesh
  int nbtriangles=Th.nt;
  cout << " nb of Triangles = " << nbtriangles << endl;
  for (int i=0;i<nbtriangles;i++)
    for (int j=0; j <3; j++)
      cout << i << " " << j << " Th[i][j] = "
           << Th[i][j] << "  x = "<< Th[i][j].x  << " , y= "<< Th[i][j].y
           << ",  label=" << Th[i][j].label << endl;
// Th(i) return the vextex i of Th
// Th[k] return the triangle k of Th
  fespace femp1(Th,P1);
  femp1 Thx=x,Thy=y;  // hack of get vertex coordinates
   // get vertices information :
  int nbvertices=Th.nv;
  cout << " nb of vertices = " << nbvertices << endl;
  for (int i=0;i<nbvertices;i++)
        cout << "Th(" <<i  << ") : "    // << endl;
             << Th(i).x << " " << Th(i).y  << " " << Th(i).label  // v 2.19
             << "       old method: " << Thx[][i] << " " << Thy[][i] << endl;
// method to find information of point (0.55,0.6)
  int it00 = Th(0.55,0.6).nuTriangle; // then triangle number
  int nr00 = Th(0.55,0.6).region;  //
// info of a triangle
  real area00 = Th[it00].area;  // new in version 2.19
  real nrr00 = Th[it00].region;  // new in version 2.19
  real nll00 = Th[it00].label;  // same as region in this case.
// Hack to get a triangle contening point x,y
// or region number (old method)
// -------------------------------------------------------
  fespace femp0(Th,P0);
  femp0 nuT;  // a P0 function to get triangle numbering
    for (int i=0;i<Th.nt;i++)
     nuT[][i]=i;
  femp0 nuReg=region;  // a P0 function to get the region number
   // inquire
  int it0=nuT(0.55,0.6); // number of triangle Th's contening (0.55,0,6);
  int nr0=nuReg(0.55,0.6); // number of region of Th's contening (0.55,0,6);
  // dump
  // -------------------------------------------------------
  cout << "  point (0.55,0,6) :triangle number " << it00 << " " << it00
       << ", region = " << nr0 << " == " << nr00 << ",  area K " << area00 << endl;
}

the output is:


 -- square mesh : nb vertices  =9 ,  nb triangles = 8 ,  nb boundary edges 8
    Nb of Vertices 9 ,  Nb of Triangles 8
    Nb of edge on user boundary  8 ,  Nb of edges on true boundary  8
 number of real boundary edges 8
 nb of Triangles = 8
0 0 Th[i][j] = 0  x = 0 , y= 0,  label=4
0 1 Th[i][j] = 1  x = 0.5 , y= 0,  label=1
0 2 Th[i][j] = 4  x = 0.5 , y= 0.5,  label=0
...
6 0 Th[i][j] = 4  x = 0.5 , y= 0.5,  label=0
6 1 Th[i][j] = 5  x = 1 , y= 0.5,  label=2
6 2 Th[i][j] = 8  x = 1 , y= 1,  label=3
7 0 Th[i][j] = 4  x = 0.5 , y= 0.5,  label=0
7 1 Th[i][j] = 8  x = 1 , y= 1,  label=3
7 2 Th[i][j] = 7  x = 0.5 , y= 1,  label=3
 Nb Of Nodes = 9
 Nb of DF = 9
 -- vector function's bound  0 1
 -- vector function's bound  0 1
 nb of vertices = 9
Th(0) : 0 0 4       old method: 0 0
Th(1) : 0.5 0 1       old method: 0.5 0
...
Th(7) : 0.5 1 3       old method: 0.5 1
Th(8) : 1 1 3       old method: 1 1
 Nb Of Nodes = 8
 Nb of DF = 8

5.1.5 The keyword ”triangulate”

FreeFem++ is able to build a triangulation from a set of points. This triangulation is a Delaunay mesh of the convex hull of the set of points. It can be useful to build a mesh form a table function.

The coordinates of the points and the value of the table function are defined separately with rows of the form: x y f(x,y) in a file such as:


0.51387 0.175741 0.636237
0.308652 0.534534 0.746765
0.947628 0.171736 0.899823
0.702231 0.226431 0.800819
0.494773 0.12472 0.580623
0.0838988 0.389647 0.456045
...............


PIC

Figure 5.6: Delaunay mesh of the convex hull of point set in file xyf

PIC

Figure 5.7: Isovalue of table function


The third column of each line is left untouched by the triangulate command. But you can use this third value to define a table function with rows of the form: x y f(x,y).

The following example shows how to make a mesh from the file “xyf” with the format stated just above. The command triangulate command use only use 1st and 2nd rows.


mesh Thxy=triangulate("xyf");  // build the Delaunay mesh of the convex hull
// points are defined by the first 2 columns of file xyf
plot(Thxy,ps="Thxyf.ps");  // (see figure 5.6)
fespace Vhxy(Thxy,P1);  // create a P1 interpolation
Vhxy fxy;  // the function
// reading the 3rd row to define the function
{ ifstream file("xyf");
   real xx,yy;
   for(int i=0;i<fxy.n;i++)
   file >> xx >>yy >> fxy[][i];   // to read third row only.
    // xx and yy are just skipped
}
plot(fxu,ps="xyf.eps");  // plot the function (see figure 5.7)

On new way to bluid a mesh from two array one the x values and the other for the y values (version 2.23-2):


Vhxy xx=x,yy=y; // to set two arrays for the x's and y's
mesh Th=triangulate(xx[],yy[]);

5.2 Boundary FEM Spaces Built as Empty Meshes

To define a Finite Element space on boundary, we came up with the idea of a mesh with no internal points (call empty mesh). It can be useful when you have a Lagrange multiplier definied on the border.

So the function emptymesh remove all the internal points of a mesh except points is on internal boundaries.


{   // new stuff 2004 emptymesh (version 1.40)
  // -- useful to build Multiplicator space
  // build a mesh without internal point
  // with the same boundary
  // -----
  assert(version>=1.40);
  border a(t=0,2⋆pi){ x=cos(t); y=sin(t);label=1;}
  mesh Th=buildmesh(a(20));
   Th=emptymesh(Th);
  plot(Th,wait=1,ps="emptymesh-1.eps"); // see figure 5.8
}

It is also possible to build an empty mesh of a pseudo subregion with emptymesh(Th,ssd) with the set of edges of the mesh Th; a edge e is in this set if with the two adjacent triangles e = t1 t2 and ssd[T1] ssd[T2] where ssd refers to the pseudo region numbering of triangles, when they are stored in an int[int] array of size the number of triangles.


{   // new stuff 2004 emptymesh (version 1.40)
  // -- useful to build Multiplicator space
  // build a mesh without internal point
  // of peusdo sub domain
  // -----
  assert(version>=1.40);
  mesh Th=square(10,10);
  int[int] ssd(Th.nt);
  for(int i=0;i<ssd.n;i++)  // build the pseudo region numbering
   {  int iq=i/2;    // because 2 triangle per quad
      int ix=iq%10;  //
      int iy=iq/10;  //
    ssd[i]= 1 + (ix>=5) +  (iy>=5)⋆2;
   }
  Th=emptymesh(Th,ssd);  // build emtpy with
   // all edge e = T1 T2 and ssd[T1] ssd[T2]
  plot(Th,wait=1,ps="emptymesh-2.eps"); // see figure 5.9
  savemesh(Th,"emptymesh-2.msh");
}


PIC

Figure 5.8: The empty mesh with boundary

PIC

Figure 5.9: An empty mesh defined from a pseudo region numbering of triangle


5.3 Remeshing

5.3.1 Movemesh

Meshes can be translated, rotated and deformed by movemesh; this is useful for elasticity to watch the deformation due to the displacement Φ(x,y) = (Φ1(x,y),Φ2(x,y)) of shape. It is also useful to handle free boundary problems or optimal shape problems.

If Ω is triangulated as Th(Ω), and Φ is a displacement vector then Φ(Th) is obtained by


mesh  Th=movemesh(Th,[Φ1,Φ2]);

Sometimes the moved mesh is invalid because some triangle becomes reversed (with a negative area). This is why we should check the minimum triangle area in the transformed mesh with checkmovemesh before any real transformation.

Example 5.2 Φ1(x,y) = x + k *sin(y *π)10), Φ2(x,y) = y + k *cos()10) for a big number k > 1.


verbosity=4;
border a(t=0,1){x=t;y=0;label=1;};
border b(t=0,0.5){x=1;y=t;label=1;};
border c(t=0,0.5){x=1-t;y=0.5;label=1;};
border d(t=0.5,1){x=0.5;y=t;label=1;};
border e(t=0.5,1){x=1-t;y=1;label=1;};
border f(t=0,1){x=0;y=1-t;label=1;};
func uu= sin(y⋆pi)/10;
func vv= cos(x⋆pi)/10;
mesh Th = buildmesh ( a(6) + b(4) + c(4) +d(4) + e(4) + f(6));
plot(Th,wait=1,fill=1,ps="Lshape.eps"); // see figure 5.10
real coef=1;
real minT0= checkmovemesh(Th,[x,y]);  // the min triangle area
while(1)  // find a correct move mesh
{
  real minT=checkmovemesh(Th,[x+coef⋆uu,y+coef⋆vv]); // the min triangle area
  if (minT > minT0/5) break ;  // if big enough
  coef=/1.5;
}
Th=movemesh(Th,[x+coef⋆uu,y+coef⋆vv]);
plot(Th,wait=1,fill=1,ps="movemesh.eps"); // see figure 5.11


PIC

Figure 5.10: L-shape

PIC

Figure 5.11: moved L-shape


Note 5.4 Consider a function u defined on a mesh Th. A statement like Th=movemesh(Th...) does not change u and so the old mesh still exists. It will be destroyed when no function use it. A statement like u = u redefines u on the new mesh Thwith interpolation and therefore destroys the old Thif u was the only function using it.

Example 5.3 (movemesh.edp) Now, we given an example of moving mesh with a lagrangian function u defined on the moving mesh.


// simple movemesh example
mesh Th=square(10,10);
fespace Vh(Th,P1);
real t=0;
// ---
// the problem is how to build data without interpolation
// so the data u is moving with the mesh as you can see in the plot
// ---
Vh u=y;
for (int i=0;i<4;i++)
{
 t=i⋆0.1;
 Vh f= x⋆t;
 real minarea=checkmovemesh(Th,[x,y+f]);
 if (minarea >0 )  // movemesh will be ok
   Th=movemesh(Th,[x,y+f]);
 cout << " Min area  " << minarea << endl;
 real[int] tmp(u[].n);
 tmp=u[];   // save the value
 u=0;         // to change the FEspace and mesh associated with u
 u[]=tmp;   // set the value of u without any mesh update
 plot(Th,u,wait=1);
};
// In this program, since u is only defined on the last mesh, all the
// previous meshes are deleted from memory.
// --------

5.4 Regular Triangulation: hTriangle

For a set S, we define the diameter of S by

diam (S ) = sup{|x - y |;x,y ∈ S}

The sequence {Th}h0 of Ω is called regular if they satisfy the following:

  1. lim m ax{diam (Tk)|Tk ∈ Th } = 0
h↓0

  2. There is a number σ > 0 independent of h such that
    --ρ(Tk)--≥ σ     for all T ∈ T
diam (Tk)                 k    h

    where ρ(Tk) are the diameter of the inscribed circle of Tk.

We put h(Th) = max{diam(Tk)|Tk ∈ Th}, which is obtained by


mesh Th = ......;
fespace Ph(Th,P0);
Ph h = hTriangle;
cout << "size of mesh = " << h[].max << endl;

5.5 Adaptmesh

The function

              3   3      -1
f (x,y) = 10.0x + y  + tan  [ε∕(sin(5.0y) - 2.0x )]     ε = 0.000 1

sharply varies in value. However, the initial mesh given by the command in Section 5.1 cannot reflect its sharp variations.

Example 5.4  


real eps =  0.0001;
real h=1;
real hmin=0.05;
func f = 10.0⋆x^3+y^3+h⋆atan2(eps,sin(5.0⋆y)-2.0⋆x);
mesh Th=square(5,5,[-1+2⋆x,-1+2⋆y]);
fespace Vh(Th,P1);
Vh fh=f;
plot(fh);
for (int i=0;i<2;i++)
 {
   Th=adaptmesh(Th,fh);
   fh=f;   // old mesh is deleted
   plot(Th,fh,wait=1);
 }


PIC


Figure 5.12: 3D graphs for the initial mesh and 1st and 2nd mesh adaptation


FreeFem++ uses a variable metric/Delaunay automatic meshing algorithm. The command


mesh ATh = adaptmesh(Th, f);

create the new mesh ATh by the Hessian

D2f = (∂2f∕∂x 2,∂ 2f∕∂x∂y,∂2f ∕∂y2)
of a function (formula or FE-function). Mesh adaptation is a very powerful tool when the solution of a problem vary locally and sharply.

Here we solve the problem (2.1)-(2.2), when f = 1 and Ω is a L-shape domain.


PIC

Figure 5.13: L-shape domain and its boundary name

PIC

Figure 5.14: Final solution after 4-times adaptation


Example 5.5 (Adapt.edp) The solution has the singularity r32,r = |x -γ|at the point γ of the intersection of two lines bc and bd (see Fig. 5.13).


border ba(t=0,1.0){x=t;   y=0;  label=1;};
border bb(t=0,0.5){x=1;   y=t;  label=1;};
border bc(t=0,0.5){x=1-t; y=0.5;label=1;};
border bd(t=0.5,1){x=0.5; y=t;  label=1;};
border be(t=0.5,1){x=1-t; y=1;  label=1;};
border bf(t=0.0,1){x=0;   y=1-t;label=1;};
mesh Th = buildmesh ( ba(6)+bb(4)+bc(4)+bd(4)+be(4)+bf(6) );
fespace Vh(Th,P1);   // set FE space
Vh u,v;              // set unknown and test function
func f = 1;
real error=0.1;         // level of error
problem Poisson(u,v,solver=CG,eps=1.0e-6) =
    int2d(Th)(  dx(u)⋆dx(v) + dy(u)⋆dy(v))
  - int2d(Th) ( f⋆v )
  + on(1,u=0)  ;
for (int i=0;i< 4;i++)
{
  Poisson;
  Th=adaptmesh(Th,u,err=error);
  error = error/2;
} ;
plot(u);

To speed up the adaptation we change by hand a default parameter err of adaptmesh, which specifies the required precision, so as to make the new mesh finer. The problem is coercive and symmetric, so the linear system can be solved with the conjugate gradient method (parameter solver=CG with the stopping criteria on the residual, here eps=1.0e-6). By adaptmesh, we get good slope of the final solution near the point of intersection of bc and bd as in Fig. 5.14.

This method is described in detail in [9]. It has a number of default parameters which can be modified :

Si f1,f2 sont des functions et thold, Thnew des maillages.


    Thnew = adaptmesh(Thold, f1  ...  );
    Thnew = adaptmesh(Thold, f1,f2  ...  ]);
    Thnew = adaptmesh(Thold, [f1,f2]  ...  );

Les paramters additionnels de adaptmesh représente par les ”...”

hmin=
Minimum edge size. (val is a real. Its default is related to the size of the domain to be meshed and the precision of the mesh generator).
hmax=
Maximum edge size. (val is a real. It defaults to the diameter of the domain to be meshed)
err=
P1 interpolation error level (0.01 is the default).
errg=
Relative geometrical error. By default this error is 0.01, and in any case it must be lower than 1√ --
  2. Meshes created with this option may have some edges smaller than the -hmin due to geometrical constraints.
nbvx=
Maximum number of vertices generated by the mesh generator (9000 is the default).
nbsmooth=
number of iterations of the smoothing procedure (5 is the default).
nbjacoby=
number of iterations in a smoothing procedure during the metric construction, 0 means no smoothing (6 is the default).
ratio=
ratio for a prescribed smoothing on the metric. If the value is 0 or less than 1.1 no smoothing is done on the metric (1.8 is the default).

If ratio > 1.1, the speed of mesh size variations is bounded by log(ratio). Note: As ratio gets closer to 1, the number of generated vertices increases. This may be useful to control the thickness of refined regions near shocks or boundary layers .

omega=
relaxation parameter for the smoothing procedure (1.0 is the default).
iso=
If true, forces the metric to be isotropic (false is the default).
abserror=
If false, the metric is evaluated using the criterium of equi-repartion of relative error (false is the default). In this case the metric is defined by
     (                            )
          1              |H |       p
M  =  ---------2  max-(Cut-Off,|η|)
      er rcoef
(5.1)

otherwise, the metric is evaluated using the criterium of equi-distribution of errors. In this case the metric is defined by

     (                          )p
M  =  ----1-----  -----|H-|-----  .
      err coef2   sup(η) - inf(η)
(5.2)

cutoff=
lower limit for the relative error evaluation (1.0e-6 is the default).
verbosity=
informational messages level (can be chosen between 0 and ). Also changes the value of the global variable verbosity (obsolete).
inquire=
To inquire graphically about the mesh (false is the default).
splitpbedge=
If true, splits all internal edges in half with two boundary vertices (true is the default).
maxsubdiv=
Changes the metric such that the maximum subdivision of a background edge is bound by val (always limited by 10, and 10 is also the default).
rescaling=
if true, the function with respect to which the mesh is adapted is rescaled to be between 0 and 1 (true is the default).
keepbackvertices=
if true, tries to keep as many vertices from the original mesh as possible (true is the default).
isMetric=
if true, the metric is defined explicitly (false is the default). If the 3 functions m11,m12,m22 are given, they directly define a symmetric matrix field whose Hessian is computed to define a metric. If only one function is given, then it represents the isotropic mesh size at every point.

For example, if the partial derivatives fxx (= 2f∕∂x2), fxy (= 2f∕∂x∂y), fyy (= 2f∕∂y2) are given, we can set

T h=a daptmes h(Th ,f xx, fxy ,f yy, IsM et ric =1, nb vx= 100 00, hm �
power=
exponent power of the Hessian used to compute the metric (1 is the default).
thetamax=
minimum corner angle in degrees (default is 0).
splitin2=
boolean value. If true, splits all triangles of the final mesh into 4 sub-triangles.
metric=
an array of 3 real arrays to set or get metric data information. The size of these three arrays must be the number of vertices. So if m11,m12,m22 are three P1 finite elements related to the mesh to adapt, you can write: metric=[m11[],m12[],m22[]] (see file convect-apt.edp for a full example)
nomeshgeneration=
If true, no adapted mesh is generated (useful to compute only a metric).
periodic=
As writing periodic=[[4,y],[2,y],[1,x],[3,x]]; it builds an adapted periodic mesh. The sample build a biperiodic mesh of a square. (see periodic finite element spaces 6, and see sphere.edp for a full example)

Example 5.6 uniformmesh.edp

We can use the command adaptmeshto build uniform mesh with a contant mesh size.

So to buid a mesh with a constant mesh size equal to 1-
30 do


mesh Th=square(2,2);  // to have initial mesh
plot(Th,wait=1,ps="square-0.eps");
Th= adaptmesh(Th,1./30.,IsMetric=1,nbvx=10000); //
plot(Th,wait=1,ps="square-1.eps");
Th= adaptmesh(Th,1./30.,IsMetric=1,nbvx=10000); // more the one time du to
Th= adaptmesh(Th,1./30.,IsMetric=1,nbvx=10000); // adaptation bound maxsubdiv=
plot(Th,wait=1,ps="square-2.eps");


PIC

Figure 5.15: Initial mesh

PIC

Figure 5.16: first iteration

PIC

Figure 5.17: last iteration


5.6 Trunc

Two operators have been introduce to remove triangles from a mesh or to divide them. Operator trunc has two parameters

label=
sets the label number of new boundary item (one by default)
split=
sets the level n of triangle splitting. each triangle is splitted in n ×n ( one by default).

To create the mesh Th3 where alls triangles of a mesh Th are splitted in 3×3 , just write:


  mesh Th3 = trunc(Th,1,split=3);

The truncmesh.edp example construct all ”trunc” mesh to the support of the basic function of the space Vh (cf. abs(u)>0), split all the triangles in 5×5, and put a label number to 2 on new boundary.


mesh Th=square(3,3);
fespace Vh(Th,P1);
Vh u;
int i,n=u.n;
u=0;
for (i=0;i<n;i++)   // all degree of freedom
 {
  u[][i]=1;         // the basic function i
  plot(u,wait=1);
  mesh Sh1=trunc(Th,abs(u)>1.e-10,split=5,label=2);
  plot(Th,Sh1,wait=1,ps="trunc"+i+".eps"); // plot the mesh of
   // the function's support
  u[][i]=0;       // reset
 }


PIC

Figure 5.18: mesh of support the function P1 number 0, splitted in 5×5

PIC

Figure 5.19: mesh of support the function P1 number 6, splitted in 5×5


5.7 Splitmesh

Another way to split mesh triangles is to use splitmesh, for example:


{   // new stuff 2004 splitmesh (version 1.37)
  assert(version>=1.37);
  border a(t=0,2⋆pi){ x=cos(t); y=sin(t);label=1;}
  plot(Th,wait=1,ps="nosplitmesh.eps");  // see figure 5.20
  mesh Th=buildmesh(a(20));
  plot(Th,wait=1);
  Th=splitmesh(Th,1+5⋆(square(x-0.5)+y⋆y));
  plot(Th,wait=1,ps="splitmesh.eps");  // see figure 5.21
}


PIC

Figure 5.20: initial mesh

PIC

Figure 5.21: all left mesh triangle is split conformaly in int(1+5⋆(square(x-0.5)+y⋆y)^2  triangles.


5.8 Meshing Examples

Example 5.7 (Two rectangles touching by a side)  


border a(t=0,1){x=t;y=0;};
border b(t=0,1){x=1;y=t;};
border c(t=1,0){x=t ;y=1;};
border d(t=1,0){x = 0; y=t;};
border c1(t=0,1){x=t ;y=1;};
border e(t=0,0.2){x=1;y=1+t;};
border f(t=1,0){x=t ;y=1.2;};
border g(t=0.2,0){x=0;y=1+t;};
int n=1;
mesh th = buildmesh(a(10⋆n)+b(10⋆n)+c(10⋆n)+d(10⋆n));
mesh TH = buildmesh ( c1(10⋆n) + e(5⋆n) + f(10⋆n) + g(5⋆n) );
plot(th,TH,ps="TouchSide.esp");  // Fig. 5.22

Example 5.8 (NACA0012 Airfoil)  


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); }
mesh Th = buildmesh(c(30)+upper(35)+lower(35));
plot(Th,ps="NACA0012.eps",bw=1);   // Fig. 5.23


PIC

Figure 5.22: Two rectangles touching by a side

PIC

Figure 5.23: NACA0012 Airfoil


Example 5.9 (Cardioid)  


real b = 1, a = b;
border C(t=0,2⋆pi) { x=(a+b)⋆cos(t)-b⋆cos((a+b)⋆t/b);
                        y=(a+b)⋆sin(t)-b⋆sin((a+b)⋆t/b); }
mesh Th = buildmesh(C(50));
plot(Th,ps="Cardioid.eps",bw=1);  // Fig. 5.24

Example 5.10 (Cassini Egg)  


border C(t=0,2⋆pi) { x=(2⋆cos(2⋆t)+3)⋆cos(t);
                      y=(2⋆cos(2⋆t)+3)⋆sin(t); }
mesh Th = buildmesh(C(50));
plot(Th,ps="Cassini.eps",bw=1);  // Fig. 5.25


PIC

Figure 5.24: Domain with Cardioid curve boundary

PIC

Figure 5.25: Domain with Cassini Egg curve boundary


Example 5.11 (By cubic Bezier curve)  


// A cubic Bezier curve connecting two points with two control points
func real bzi(real p0,real p1,real q1,real q2,real t)
{
  return p0⋆(1-t)^3+q1⋆3⋆(1-t)^2⋆t+q2⋆3⋆(1-t)⋆t^2+p1⋆t^3;
}
real[int] p00=[0,1], p01=[0,-1], q00=[-2,0.1], q01=[-2,-0.5];
real[int] p11=[1,-0.9], q10=[0.1,-0.95], q11=[0.5,-1];
real[int] p21=[2,0.7], q20=[3,-0.4], q21=[4,0.5];
real[int] q30=[0.5,1.1], q31=[1.5,1.2];
border G1(t=0,1) { x=bzi(p00[0],p01[0],q00[0],q01[0],t);
                   y=bzi(p00[1],p01[1],q00[1],q01[1],t); }
border G2(t=0,1) { x=bzi(p01[0],p11[0],q10[0],q11[0],t);
                   y=bzi(p01[1],p11[1],q10[1],q11[1],t); }
border G3(t=0,1) { x=bzi(p11[0],p21[0],q20[0],q21[0],t);
                   y=bzi(p11[1],p21[1],q20[1],q21[1],t); }
border G4(t=0,1) { x=bzi(p21[0],p00[0],q30[0],q31[0],t);
                   y=bzi(p21[1],p00[1],q30[1],q31[1],t); }
int m=5;
mesh Th = buildmesh(G1(2⋆m)+G2(m)+G3(3⋆m)+G4(m));
plot(Th,ps="Bezier.eps",bw=1);   // Fig 5.26

Example 5.12 (Section of Engine)  


real a= 6., b= 1., c=0.5;
border L1(t=0,1) { x= -a; y= 1+b - 2⋆(1+b)⋆t; }
border L2(t=0,1) { x= -a+2⋆a⋆t; y= -1-b⋆(x/a)⋆(x/a)⋆(3-2⋆abs(x)/a );}
border L3(t=0,1) { x= a; y=-1-b + (1+ b )⋆t; }
border L4(t=0,1) { x= a - a⋆t;   y=0; }
border L5(t=0,pi) { x= -c⋆sin(t)/2; y=c/2-c⋆cos(t)/2; }
border L6(t=0,1) { x= a⋆t;  y=c; }
border L7(t=0,1) { x= a;  y=c + (1+ b-c )⋆t; }
border L8(t=0,1) { x= a-2⋆a⋆t; y= 1+b⋆(x/a)⋆(x/a)⋆(3-2⋆abs(x)/a); }
mesh Th = buildmesh(L1(8)+L2(26)+L3(8)+L4(20)+L5(8)+L6(30)+L7(8)+L8(30));
plot(Th,ps="Engine.eps",bw=1);  // Fig. 5.27


PIC

Figure 5.26: Boundary drawed by Bezier curves

  
PIC

Figure 5.27: Section of Engine


Example 5.13 (Domain with U-shape channel)  


real d = 0.1;  // width of U-shape
border L1(t=0,1-d) { x=-1; y=-d-t; }
border L2(t=0,1-d) { x=-1; y=1-t; }
border B(t=0,2) { x=-1+t; y=-1; }
border C1(t=0,1) { x=t-1; y=d; }
border C2(t=0,2⋆d) { x=0; y=d-t; }
border C3(t=0,1) { x=-t; y=-d; }
border R(t=0,2) { x=1; y=-1+t; }
border T(t=0,2) { x=1-t; y=1; }
int n = 5;
mesh Th = buildmesh (L1(n/2)+L2(n/2)+B(n)+C1(n)+C2(3)+C3(n)+R(n)+T(n));
plot(Th,ps="U-shape.eps",bw=1);  // Fig 5.28

Example 5.14 (Domain with V-shape cut)  


real dAg = 0.01;  // angle of V-shape
border C(t=dAg,2⋆pi-dAg) { x=cos(t); y=sin(t); };
real[int] pa(2), pb(2), pc(2);
pa[0] = cos(dAg); pa[1] = sin(dAg);
pb[0] = cos(2⋆pi-dAg); pb[1] = sin(2⋆pi-dAg);
pc[0] = 0; pc[1] = 0;
border seg1(t=0,1) { x=(1-t)⋆pb[0]+t⋆pc[0]; y=(1-t)⋆pb[1]+t⋆pc[1]; };
border seg2(t=0,1) { x=(1-t)⋆pc[0]+t⋆pa[0]; y=(1-t)⋆pc[1]+t⋆pa[1]; };
mesh Th = buildmesh(seg1(20)+C(40)+seg2(20));
plot(Th,ps="V-shape.eps",bw=1);   // Fig. 5.29


PIC

Figure 5.28: Domain with U-shape channel changed by d

PIC

Figure 5.29: Domain with V-shape cut changed by dAg


Example 5.15 (Smiling face)  


real d=0.1;
int m=5;
real a=1.5, b=2, c=0.7, e=0.01;
border F(t=0,2⋆pi) { x=a⋆cos(t); y=b⋆sin(t); }
border E1(t=0,2⋆pi) { x=0.2⋆cos(t)-0.5; y=0.2⋆sin(t)+0.5; }
border E2(t=0,2⋆pi) { x=0.2⋆cos(t)+0.5; y=0.2⋆sin(t)+0.5; }
func real st(real t) {
   return sin(pi⋆t)-pi/2;
}
border C1(t=-0.5,0.5) { x=(1-d)⋆c⋆cos(st(t)); y=(1-d)⋆c⋆sin(st(t)); }
border C2(t=0,1){x=((1-d)+d⋆t)⋆c⋆cos(st(0.5));y=((1-d)+d⋆t)⋆c⋆sin(st(0.5));}
border C3(t=0.5,-0.5) { x=c⋆cos(st(t)); y=c⋆sin(st(t)); }
border C4(t=0,1) { x=(1-d⋆t)⋆c⋆cos(st(-0.5)); y=(1-d⋆t)⋆c⋆sin(st(-0.5));}
border C0(t=0,2⋆pi) { x=0.1⋆cos(t); y=0.1⋆sin(t); }
mesh Th=buildmesh(F(10⋆m)+C1(2⋆m)+C2(3)+C3(2⋆m)+C4(3)
                  +C0(m)+E1(-2⋆m)+E2(-2⋆m));
plot(Th,ps="SmileFace.eps",bw=1);   // see Fig. 5.30
}

Example 5.16 (3point bending)  


// Square for Three-Point Bend Specimens fixed on Fix1, Fix2
// It will be loaded on Load.
real a=1, b=5, c=0.1;
int n=5, m=b⋆n;
border Left(t=0,2⋆a) { x=-b; y=a-t; }
border Bot1(t=0,b/2-c) { x=-b+t; y=-a; }
border Fix1(t=0,2⋆c) { x=-b/2-c+t; y=-a; }
border Bot2(t=0,b-2⋆c) { x=-b/2+c+t; y=-a; }
border Fix2(t=0,2⋆c) { x=b/2-c+t; y=-a; }
border Bot3(t=0,b/2-c) { x=b/2+c+t; y=-a; }
border Right(t=0,2⋆a) { x=b; y=-a+t; }
border Top1(t=0,b-c) { x=b-t; y=a; }
border Load(t=0,2⋆c) { x=c-t; y=a; }
border Top2(t=0,b-c) { x=-c-t; y=a; }
mesh Th = buildmesh(Left(n)+Bot1(m/4)+Fix1(5)+Bot2(m/2)+Fix2(5)+Bot3(m/4)
                    +Right(n)+Top1(m/2)+Load(10)+Top2(m/2));
plot(Th,ps="ThreePoint.eps",bw=1);  // Fig. 5.31


PIC

Figure 5.30: Smiling face (Mouth is changeable)

  
PIC

Figure 5.31: Domain for three-point bending test


5.9 How to change the label of elements and border elements of a mesh in FreeFem++ ?

Changing the label of elements and border elements will be done using the keyword change. The parameters for this command line are for a two dimensional case:

refe =
is a vector of integer that contains the old labels and the new labels of edges.
reft =
is a vector of integer that contains the old labels and the new labels of triangles.

and for a three dimensional case:

reftet =
is a vector of integer that contains the old labels and the new labels of tetrahedrons.
refface =
is a vector of integer that contains the old labels and the new labels of triangles.

These vectors are composed of nl succesive sub vectors of size two where nl is the number of label that we want to change. The first element and second element of sub vector are respectively the old label and the new label. For example, we have

refe = [V 1,...,Vn ].
                l
(5.3)

where the sub vector Vi is defined by Vi = [old label of set i, new label of set i].

An example of using this function is given in ”glumesh2D.edp”:

Example 5.17 (glumesh2D.edp)  


 1:
 2: mesh Th1=square(10,10);
 3: mesh Th2=square(20,10,[x+1,y]);
 4: verbosity=3;
 5: int[int] r1=[2,0],  r2=[4,0];
 6: plot(Th1,wait=1);
 7: Th1=change(Th1,refe=r1); // Change the label of Edges of Th1 with label 2 into label 0.
 8: plot(Th1,wait=1);
 9: Th2=change(Th2,refe=r2); // Change the label of Edges of Th2 with label 4 into label 0.
10: mesh Th=Th1+Th2;          // ‘‘gluing together'' of meshes Th1 and Th2
11: cout << " nb ref = " << int1d(Th1,1,3,4)(1./lenEdge)+int1d(Th2,1,2,3)(1./lenEdge)
12:          << " == " << int1d(Th,1,2,3,4)(1./lenEdge) <<" == " << ((10+20)+10)⋆2 << endl;
13: plot(Th,wait=1);
14: fespace Vh(Th,P1);
15: macro Grad(u) [dx(u),dy(u)];  // definition of a macro
16: Vh u,v;
17: solve P(u,v)=int2d(Th)(Grad(u)'⋆Grad(v))-int2d(Th)(v)+on(1,3,u=0);
18: plot(u,wait=1);

“gluing” different mesh In line 10 of previous file, the method to “gluing” different mesh of the same dimension in FreeFem++ is using. This function is just call using the symbol addition ”+” between meshes. The method implemented need that the point in adjacent mesh are the same.

5.10 Mesh in three dimension

5.10.1 Read/Write Statements for a Mesh in 3D

In three dimension, the file mesh format supported for input and output files by FreeFem++ are the extension .msh and .mesh. These format are described in Chapter Mesh Files in two dimension.

extension file .msh The structure of the extension .msh in three dimensional are given in Table 5.2. In this structure, nv denotes the number of vertices, ntet the number of tetrahedra and ntri the number of triangles For each vertex qi,i = 1,⋅⋅⋅,n v, we denote by (qxi,q yi,q zi) the x-coordinate, the y-coordinate and the z-coordinate. Each tetrahedra Tk,k = 1,⋅⋅⋅,ntet has four vertices qk1,qk2,qk3,qk4. The boundary consists of an union of triangles. Each triangle bej,j = 1,⋅⋅⋅,ntri has three vertices qj1,qj2,qj3.


     





nv ntet ntri
qx1 qy1 qz1 Vertex label
qx2 qy2 qz2 Vertex label
.
.. .
.. .
.. .
..
qxnv q ynv q znv Vertex label
11 12 13 14 region label
21 22 23 24 region label
.
.. .
.. .
.. .
.. .
..
(ntet)1 (ntet)2(ntet)3 (ntet)4 region label
11 12 13 boundary label
21 22 23 boundary label
.
.. .
.. .
.. .
..
(ntri)1(ntri)2(ntri)3boundary label






Table 5.2: The structure of mesh file format “.msh” in three dimension.

extension file .mesh The data structure for a three dimensional mesh is composed of the data structure presented in Section ?? and a data structure for tetrahedrons. The tetrahedrons of a three dimensional mesh are refereed using the following field:

This field is express with the notation of Section ??.

5.10.2 TeGen: A tetrahedral mesh generator

TetGen TetGen is a software developed by Hang Si of Weierstrass Institute for Applied Analysis and Stochastics of Berlin in Germany [36]. TetGen is a free for research and non-commercial uses. For any commercial licence utilization, a commercial licence is available upon request of Hang Si.

This software is a tetrahedral mesh generator of a three dimensional domain defined by its boundary. The input domain take into account a polyhedral or a piecewise linear complex. This tetrahelization is a constrained Delaunay tetrahelization.

The method used in TetGen to control the quality of the mesh is a Delaunay refinement from Shewchuk [37]. The quality measure of this algorithm is the Radius-Edge Ratio (see Section 1.3.1 [36] for more details). A theorical bounds of this ratio of the algorithm of Shewchuk is obtained for a given complex of vertices, constrained segments and facets of surface mesh, with no input angle less than 90 degree. This theorical bounds is 2.0.

The call of Tetgen is done with the keyword tetg. The parameters of this command line is:

refface =
is a vector of integer that contains the old labels and the new labels of Triangles. This parameters is initialized as refface for the keyword change (5.3).
switch =
A string expression. This string corresponds to the command line switch of Tetgen see Section 3.2 of [36].
nbofholes=
Number of holes.
holelist =
This array correspond to holelist of tetgenio data structure [36]. A real vector of size 3 ×nbofholes. In TetGen, each hole is associated with a point inside this domain. This vector is x1h,y1h,z1h,x2h,y2h,z2h,⋅⋅⋅, where xih,yih,zih is the associated point with the ith hole.
nbofregions =
Number of regions.
regionlist =
This array correspond to regionlist of tetgenio data structure [36]. The attribute and the volume constraint of region are given in this real vector of size 5 ×nbofregions. The ith region is described by five elements: x-coordinate, y-coordinate and z-coordinate of a point inside this domain (xi,yi,zi); the attribute (ati) and the maximum volume for tetrahedrons (mvoli) for this region. The regionlistvector is: x1,y1,z1,at1,mvol1,x2,y2,z2,at2,mvol2,⋅⋅⋅.
nboffacetcl=
Number of facets constraints.
facetcl=
This array correspond to facetconstraintlist of tetgenio data structure [36]. The ith facet constraint is defined by the facet marker Refifc and the maximum area for faces mareaifc. The facetclarray is: Ref1fc,marea1fc,Ref2fc,marea2fc,⋅⋅⋅. This parameters has no effect if switch qis not selected.

Principal switch parameters in TetGen:

To obtain a tetrahedral mesh generator with tetgen, we need the surface mesh of three dimensional domain. We give now the command line in FreeFem++ to construct these meshes.

keyword: “movemesh23” A simple method to construct a surface is to place a two dimensional domain in a three dimensional space. This corresponding to move the domain by a displacement vector of this form Φ(x,y) = (Φ1(x,y),Φ2(x,y),Φ3(x,y)). The result of moving a two dimensional mesh Th2 by this three dimensional displacement is obtained using:


mesh3 Th3 = movemesh23(Th2,transfo=[Φ1,Φ2,Φ3]);

The parameters of this command line are:

transfo =
[Φ1, Φ2, Φ3] set the displacement vector of transformation Φ(x,y) = [Φ1(x,y),Φ2(x,y),Φ3(x,y)].
refface =
set integer label of triangles
orientation=
set integer orientation of mesh.
ptmerge =
A real expression. When you transform a mesh, some points can be merged. This parameters is the criteria to define two merging points. By default, we use
ptmerge = 1e - 7Vol(B),
where B is the smallest axis parallel boxes containing the discretize domain of Ω and Vol(B) is the volume of this box.

We can “gluing” surface meshes using the process given in Section 5.9. An example of obtaining a three dimensional mesh using the command line tetgand movemesh23is given in the file tetgencube.edp.

Example 5.18 (tetgencube.edp)  


// file tetgencube.edp
load "msh3"
load "tetgen"
real x0,x1,y0,y1;
x0=1.; x1=2.; y0=0.; y1=2⋆pi;
mesh Thsq1 = square(5,35,[x0+(x1-x0)⋆x,y0+(y1-y0)⋆y]);
func ZZ1min = 0;
func ZZ1max = 1.5;
func XX1 = x;
func YY1 = y;
mesh3 Th31h = movemesh23(Thsq1,transfo=[XX1,YY1,ZZ1max]);
mesh3 Th31b = movemesh23(Thsq1,transfo=[XX1,YY1,ZZ1min]);
// ///////////////////////////////
x0=1.; x1=2.; y0=0.; y1=1.5;
mesh Thsq2 = square(5,8,[x0+(x1-x0)⋆x,y0+(y1-y0)⋆y]);
func ZZ2 = y;
func XX2 = x;
func YY2min = 0.;
func YY2max = 2⋆pi;
mesh3 Th32h = movemesh23(Thsq2,transfo=[XX2,YY2max,ZZ2]);
mesh3 Th32b = movemesh23(Thsq2,transfo=[XX2,YY2min,ZZ2]);
// ///////////////////////////////
x0=0.; x1=2⋆pi; y0=0.; y1=1.5;
mesh Thsq3 = square(35,8,[x0+(x1-x0)⋆x,y0+(y1-y0)⋆y]);
func XX3min = 1.;
func XX3max = 2.;
func YY3 = x;
func ZZ3 = y;
mesh3 Th33h = movemesh23(Thsq3,transfo=[XX3max,YY3,ZZ3]);
mesh3 Th33b = movemesh23(Thsq3,transfo=[XX3min,YY3,ZZ3]);
// //////////////////////////////
mesh3 Th33 = Th31h+Th31b+Th32h+Th32b+Th33h+Th33b; // "gluing" surface meshs to obtain the surface of cube
savemesh(Th33,"Th33.mesh");
// build a mesh of a axis parallel box with TetGen
real[int] domain =[1.5,pi,0.75,145,0.0025];
mesh3 Thfinal = tetg(Th33,switch="paAAQY",nbofregions=1,regionlist=domain);     // Tetrahelize the interior of the cube with tetgen
savemesh(Thfinal,"Thfinal.mesh");
// build a mesh of a half cylindrical shell of interior radius 1. and exterior radius 2 and heigh 1.5
func mv2x = x⋆cos(y);
func mv2y = x⋆sin(y);
func mv2z = z;
mesh3 Thmv2 = movemesh3(Thfinal, transfo=[mv2x,mv2y,mv2z]);
savemesh(Thmv2,"halfcylindricalshell.mesh")

The command movemeshis describe in the following section.

The keyword “tetgtransfo” This keyword correspond to a composition of command line tetgand movemesh23:


tetgtransfo( Th2, transfo= [Φ1, Φ2, Φ3] ), ⋅⋅⋅ ) = tetg( Th3surf, ⋅⋅⋅ ),

where Th3surf = movemesh23( Th2,tranfo=[Φ1, Φ2, Φ3] ) and Th2 is the input two dimensional mesh of tetgtransfo.

The parameters of this command line are on the one hand the parameters:
      refface, switch, regionlistnboffacetclfacetcl
of keyword tetgand on the other hand the parameter ptmergeof keyword movemesh23.

Remark: To use tetgtransfo, the result’s mesh of movemesh23must be an enclosed surface and defined one region. Therefore, the parameter regionlistis defined for one region.

An example of this keyword can be found in line of file “buildlayers.edp”

The keyword ”tetgconvexhull” FreeFem++ , using tetgen, is able to build a tetrahelization from a set of points. This tetrahelization is a Delaunay mesh of the convex hull of the set of points.

The coordinates of the points can be initializes in two way. The first is a file that contains the coordinate of points Xi = (xi,yi,zi). This files is organized as follow:

nv
x1   y1  z1
x2   y2  z2
 ..   ..    ..
 .   .    .
xnv ynv  znv
The second way is to give three arrays that corresponding respectively to x-coordinates, y-coordinates and z-coordinates.

The parameters of this command line are

switch =
A string expression. This string corresponds to the command line switch of TetGen see Section 3.2 of [36].
reftet =
An integer expression. set the label of tetrahedrons.
refface =
An integer expression. set the label of triangles.

In the string switch, we can’t used the option ’p’ and ’q’ of tetgen.

5.10.3 Reconstruct/Refine a three dimensional mesh with TetGen

Meshes in three dimension can be refined using TetGen with the command line tetgreconstruction.

The parameter of this keyword are

reftet=
an integer array that allow to change the label of tetrahedrons. This array is defined as the parameter reftetin the keyword change.
refface=
an integer array that allow to change the label of triangles. This array is defined as the parameter reffacein the keyword change.
sizevolume=
a function. This function allows to constraint volume size of tetrahedrons in the domain.

The parameter switchnbofregions, regionlist, nboffacetcland facetclof the command line which call TetGen (tetg) is used for tetgrefine.

In the parameter switch=, the character ’r’ should be used without the character ’p’. For instance, see the manual of TetGen [36] for effect of ’r’ to other character.

The parameter regionlistallow to define a new volume contraint in the region. The label in the regionlist will be the previous label of region. This parameter and nbofregionscan’t be used with parameter sizevolume.

Example:

Example 5.19 (refinesphere.edp)


// file refinesphere.edp
load "msh3"
load "tetgen"
load "medit"
mesh Th=square(10,20,[x⋆pi-pi/2,2⋆y⋆pi]);   // ]-p2i,frac-pi2[×]0, 2π[
// a parametrization of a sphere
func f1 =cos(x)⋆cos(y);
func f2 =cos(x)⋆sin(y);
func f3 = sin(x);
// partiel derivative of the parametrization DF
func f1x=sin(x)⋆cos(y);
func f1y=-cos(x)⋆sin(y);
func f2x=-sin(x)⋆sin(y);
func f2y=cos(x)⋆cos(y);
func f3x=cos(x);
func f3y=0;
// M = DFtDF
func m11=f1x^2+f2x^2+f3x^2;
func m21=f1x⋆f1y+f2x⋆f2y+f3x⋆f3y;
func m22=f1y^2+f2y^2+f3y^2;
func perio=[[4,y],[2,y],[1,x],[3,x]];
real hh=0.1;
real vv= 1/square(hh);
verbosity=2;
Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,periodic=perio);
Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,periodic=perio);
plot(Th,wait=1);
verbosity=2;
// construction of the surface of spheres
real Rmin  = 1.;
func f1min = Rmin⋆f1;
func f2min = Rmin⋆f2;
func f3min = Rmin⋆f3;
mesh3 Th3=movemesh23(Th,transfo=[f1min,f2min,f3min]);
real[int] domain = [0.,0.,0.,145,0.01];
mesh3 Th3sph=tetg(Th3,switch="paAAQYY",nbofregions=1,regionlist=domain);
int[int] newlabel = [145,18];
real[int] domainrefine = [0.,0.,0.,145,0.0001];
mesh3 Th3sphrefine=tetgreconstruction(Th3sph,switch="raAQ",reftet=newlabel,
nbofregions=1,regionlist=domainrefine);
int[int] newlabel2 = [145,53];
func fsize = 0.01/(( 1 + 5⋆sqrt( (x-0.5)^2+(y-0.5)^2+(z-0.5)^2) )^3);
mesh3 Th3sphrefine2=tetgreconstruction(Th3sph,switch="raAQ",reftet=newlabel2,
sizeofvolume=fsize);
medit(‘‘sphere'',Th3sph);
medit(‘‘isotroperefine'' ,Th3sphrefine);
medit(‘‘anisotroperefine'',Th3sphrefine2);

5.10.4 Moving mesh in three dimension

Meshes in three dimension can be translated rotated and deformed using the command line movemesh as in the 2D case (see section movemesh in chapiter 5). If Ω is tetrahelize as Th(Ω), and Φ(x,y) = (Φ1(x,y,z),Φ1(x,y,z),Φ3(x,y,z)) is a displacement vector then Φ(Th) is obtained by


mesh3 Th = movemesh( Th, transfo=[Φ1, Φ2, Φ3], ... );

The parameters of movemesh in three dimension are

transfo =
[Φ1,Φ2, Φ3] set the displacement vector of 3D transformation [Φ1(x,y,z),Φ2(x,y,z),Φ3(x,y,z)].
reftet =
set integer label of tetrahedrons. 0 by default.
refface =
set the label of faces of border. This parameters is initialized as refface for the keyword change (5.3).
facemerge =
An integer expression. When you transform a mesh, some faces can be merged. This parameters equals to one if merge’s faces is considered. Otherwise equals to zero. By default, this parameter is equals to 1.
ptmerge =
A real expression. When you transform a mesh, some points can be merged. This parameters is the criteria to define two merging points. By default, we use
ptmerge = 1e - 7Vol(B),
where B is the smallest axis parallel boxes containing the discretion domain of Ω and Vol(B) is the volume of this box.

An example of this command can be found in the file ”Poisson3d.edp” situed in the directory examples++-3d.

5.10.5 Layer mesh

In this section, we present the command line to obtain a Layer mesh: buildlayermesh. This mesh is obtaining in extending a two dimensional mesh in the z-axis.

The domain Ω3d defined by the layer mesh is equals to Ω3d = Ω2d ×[zmin,zmax] where Ω2d is the domain define by the two dimensional mesh, zmin and zmax are function of Ω2d in R that defines respectively the lower surface and upper surface of Ω3d.


      PIC

Figure 5.32: Example of Layer mesh in three dimension.


For a vertex of two dimensional mesh Vi2d = (xi,yi), we introduce the number of associated vertices in the z-axis Mi + 1. We denote by M the maximum of Mi over the vertices of the two dimensional mesh. This value are called the number of layers (if i,Mi = M then there are M layers in the mesh of Ω3d). Vi2d generated M + 1 vertices which are defined by

                   3d
∀j = 0,...,M,     Vi,j = (xi,yi,θi(zi,j)),
where (zi,j)j=0,,M are the M + 1 equidistant points on the interval [zmin(Vi2d),zmax(Vi2d)]:
                                   2d         2d
zi,j = jδα +zmin(V2d),   δα = zmax(Vi-)--zmin(V-i-).
                i                     M
The function θi, defined on [zmin(Vi2d),zmax(Vi2d)], is given by
      (
      ||| θi,0 ifz = zmin(V2id),
θi(z) = {|| θi,j ifz ∈]θi,j-1,θi,j],
      |(
with (θi,j)j=0,,Mi are the Mi + 1 equidistant points on the interval [zmin(Vi2d),zmax(Vi2d)].

Set a triangle K = (Vi12d, Vi22d, Vi32d) of the two dimensional mesh. K is associated with a triangle on the upper surface (resp. on the lower surface) of layer mesh: (Vi1,M3d,Vi2,M3d,Vi3,M3d) (resp. (Vi1,03d,Vi2,03d,Vi3,03d)).

Also K is associated with M volume primatic elements which are defined by

∀j = 0,...,M,   Hj = (V 3i1d,j,V3id2,j,V 3i3d,j,V3id1,j+1,V3id2,j+1,Vi3d3,j+1).

Theses volume elements can be have some merged point:

The volume with merged points are called degenerate elements. To obtain a meshing with tetrahedrons, we decompose the pyramid into two tetrahedrons and the prism into three tetrahedrons. These tetrahedons are obtaining in cut the quadrilateral face of pyramid and prism with the diagonal with have the vertex with the maximum index (see [8] for explication of this choice).

The triangles on the middle surface obtained with the decomposition of the volume prismatic elements are the triangles generated by the edges on the border of the two dimensional mesh. The label of triangles on the border elements and tetrahedrons are defined with the label of these associated elements.

The arguments of buildlayermeshis a two dimensional mesh and the number of layers M.

The parameters of this command are:

zbound =
[zmin,zmax] where zmin and zmax are functions expression. Theses functions define the lower surface mesh and upper mesh of surface mesh.
coef =
A function expression between [0,1]. This parameter is used to introduce degenerate element in mesh. The number of associated points of vertex Vi2d is the integer part of coef(Vi2d)M.
reftet =
This vector is used to initialized the labels of tetrahedrons. This vector contains these labels and labels of triangles of the 2D mesh.
reffacemid =
This vector is used to initialized the labels of triangles of the layermesh of middle surface mesh. This vector contains these labels and labels of edges of the 2D mesh.
reffaceup =
This vector is used to initialized the label of triangles of the layermesh of upper surface mesh. This vector contains these labels and labels of triangles of the 2D mesh.
reffacelow =
This vector is used to initialized the label of triangles of the layermesh of lower surface mesh. This vector contains these labels and labels of triangles of the 2D mesh.

Moreover, we also add post processing parameters that allow to moving the mesh. These parameters correspond to parameters transfo, facemergeand ptmergeof the command line movemesh.

The vector reftet, reffacemid, reffaceupand reffaceloware composed of nl succesive sub vectors of size two where nl is the number of label that we want to initialized. The first element and second element of sub vector are respectively a label of elements of the two dimensional mesh and the label of elements of layer mesh which is associated with these elements.

An example of this command line is given in “buildlayermesh.edp”.

Example 5.20 (buildlayermesh.edp)  


// file buildlayermesh.edp
load "msh3"
load "tetgen"
// Test 1
int C1=99, C2=98;  // could be anything
border C01(t=0,pi){ x=t;  y=0;      label=1;}
border C02(t=0,2⋆pi){ x=pi; y=t;  label=1;}
border C03(t=0,pi){ x=pi-t;  y=2⋆pi;    label=1;}
border C04(t=0,2⋆pi){ x=0;    y=2⋆pi-t; label=1;}
border C11(t=0,0.7){ x=0.5+t;  y=2.5;      label=C1;}
border C12(t=0,2){ x=1.2;    y=2.5+t;  label=C1;}
border C13(t=0,0.7){ x=1.2-t;  y=4.5;     label=C1;}
border C14(t=0,2){ x=0.5;    y=4.5-t; label=C1;}
border C21(t=0,0.7){ x= 2.3+t;     y=2.5;  label=C2;}
border C22(t=0,2){        x=3;   y=2.5+t;  label=C2;}
border C23(t=0,0.7){   x=3-t;     y=4.5;  label=C2;}
border C24(t=0,2){       x=2.3;   y=4.5-t; label=C2;}
mesh Th=buildmesh(    C01(10)+C02(10)+ C03(10)+C04(10)
                    + C11(5)+C12(5)+C13(5)+C14(5)
                    + C21(-5)+C22(-5)+C23(-5)+C24(-5));
mesh Ths=buildmesh(    C01(10)+C02(10)+ C03(10)+C04(10)
                    + C11(5)+C12(5)+C13(5)+C14(5) );
// construction of a box with one hole and two regions
func zmin=0.;
func zmax=1.;
int MaxLayer=10;
func XX = x⋆cos(y);
func YY = x⋆sin(y);
func ZZ = z;
int[int] r1=[0,41], r2=[98,98,99,99,1,56];
int[int] r3=[4,12];    // The triangles of uppper surface mesh generated by the triangle in the 2D region of mesh Th of label 4 as label 12.
int[int] r4=[4,45];     // The triangles of lower surface mesh generated by the triangle in the 2D region of mesh Th of label 4 as label 45.
mesh3 Th3=buildlayers( Th, MaxLayer, zbound=[zmin,zmax], reftet=r1, reffacemid=r2, reffaceup = r3, reffacelow = r4 );
savemesh(Th3,"box2region1hole.mesh");
// construction of a sphere with TetGen
func XX1 = cos(y)⋆sin(x);
func YY1 = sin(y)⋆sin(x);
func ZZ1 = cos(x);
string test="paACQ";
cout << "test=" << test << endl;
mesh3 Th3sph=tetgtransfo(Ths,transfo=[XX1,YY1,ZZ1],switch=test,nbofregions=1,regionlist=domain);
savemesh(Th3sph,"sphere2region.mesh");

5.11 Meshing examples

Example 5.21 (lac.edp) // file ”lac.edp”


load ‘‘msh3''
int nn=5;
border cc(t=0,2⋆pi){x=cos(t);y=sin(t);label=1;}
mesh Th2 = buildmesh(cc(100));
fespace Vh2(Th2,P2);
Vh2 ux,uy,p2;
int[int] rup=[0,2], rdlow=[0,1], rmid=[1,1,2,1,3,1,4,1];
func zmin = 2-sqrt(4-(x⋆x+y⋆y));
func zmax = 2-sqrt(3.);
mesh3 Th = buildlayers(Th2,nn,
  coeff = max((zmax-zmin)/zmax, 1./nn),
  zbound=[zmin,zmax],
  reffacemid=rmid;
  reffaceup=rup;
  reffacelow=rlow);
savemesh(Th,''Th.meshb'');
exec(‘‘medit Th; Th.meshb'');

Example 5.22 (tetgenholeregion.edp)


// file ‘‘tetgenholeregion.edp''
load "msh3''
load "tetgen"
mesh Th=square(10,20,[x⋆pi-pi/2,2⋆y⋆pi]);   // ]-pi
 2,frac-pi2[×]0, 2π[
// a parametrization of a sphere
func f1 =cos(x)⋆cos(y);
func f2 =cos(x)⋆sin(y);
func f3 = sin(x);
// partiel derivative of the parametrization DF
func f1x=sin(x)⋆cos(y);
func f1y=-cos(x)⋆sin(y);
func f2x=-sin(x)⋆sin(y);
func f2y=cos(x)⋆cos(y);
func f3x=cos(x);
func f3y=0;
// M = DFtDF
func m11=f1x^2+f2x^2+f3x^2;
func m21=f1x⋆f1y+f2x⋆f2y+f3x⋆f3y;
func m22=f1y^2+f2y^2+f3y^2;
func perio=[[4,y],[2,y],[1,x],[3,x]];
real hh=0.1;
real vv= 1/square(hh);
verbosity=2;
Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,periodic=perio);
Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,periodic=perio);
plot(Th,wait=1);
verbosity=2;
// construction of the surface of spheres
real Rmin  = 1.;
func f1min = Rmin⋆f1;
func f2min = Rmin⋆f2;
func f3min = Rmin⋆f3;
mesh3 Th3sph = movemesh23(Th,transfo=[f1min,f2min,f3min]);
real Rmax  = 2.;
func f1max = Rmax⋆f1;
func f2max = Rmax⋆f2;
func f3max = Rmax⋆f3;
mesh3 Th3sph2 = movemesh23(Th,transfo=[f1max,f2max,f3max]);
cout << "addition" << endl;
mesh3 Th3 = Th3sph+Th3sph2;
real[int] domain2 = [1.5,0.,0.,145,0.001,0.5,0.,0.,18,0.001];
cout << "==============================" << endl;
cout << " tetgen call without hole " << endl;
cout << "==============================" << endl;
mesh3 Th3fin = tetg(Th3,switch="paAAQYY",nbofregions=2,regionlist=domain2);
cout << "=============================" << endl;
cout << "finish tetgen call without hole" << endl;
cout << "=============================" << endl;
savemesh(Th3fin,"spherewithtworegion.mesh");
real[int] hole = [0.,0.,0.];
real[int] domain = [1.5,0.,0.,53,0.001];
cout << "=============================" << endl;
cout << "  tetgen call with hole   " << endl;
cout << "=============================" << endl;
mesh3 Th3finhole=tetg(Th3,switch="paAAQYY",nbofholes=1,holelist=hole,
nbofregions=1,regionlist=domain);
cout << "=============================" << endl;
cout << "finish tetgen call with hole   " << endl;
cout << "=============================" << endl;
savemesh(Th3finhole,"spherewithahole.mesh");

5.11.1 Build a 3d mesh of a cube with a ballon incrustation


 mesh Th;  // to store the final mesh
{  // to clean all intermedial variable
  mesh3 ThHex;
  real volumetet;   // use in tetg.
  {
     // first build the 6 faces of the hex.
    real x0=-1,x1=1;
    real y0=-1.1,y1=1.1;
    real z0=-1.2,z1=1.2;
    int nx=19,ny=20,nz=21;
     // a volume of on tet.
    volumetet= (x1-x0)⋆(y1-y0)⋆(z1-z0)/ (nx⋆ny⋆ny) /6.;
    mesh Thx = square(ny,nz,[y0+(y1-y0)⋆x,z0+(z1-z0)⋆y]);
    mesh Thy = square(nx,nz,[x0+(x1-x0)⋆x,z0+(z1-z0)⋆y]);
    mesh Thz = square(nx,ny,[x0+(x1-x0)⋆x,y0+(y1-y0)⋆y]);
    int[int] refz=[0,5];   // bas
    int[int] refZ=[0,6];    // haut
    int[int] refy=[0,3];   // devant
    int[int] refY=[0,4];    // derriere
    int[int] refx=[0,1];   // gauche
    int[int] refX=[0,2];    // droite
     // buil the mesh of the 6 faces ..
    mesh3 Thx0 = movemesh23(Thx,transfo=[x0,x,y],orientation=-1,refface=refx);
    mesh3 Thx1 = movemesh23(Thx,transfo=[x1,x,y],orientation=1,refface=refX);
    mesh3 Thy0 = movemesh23(Thy,transfo=[x,y0,y],orientation=+1,refface=refy);
    mesh3 Thy1 = movemesh23(Thy,transfo=[x,y1,y],orientation=-1,refface=refY);
    mesh3 Thz0 = movemesh23(Thz,transfo=[x,y,z0],orientation=-1,refface=refz);
    mesh3 Thz1 = movemesh23(Thz,transfo=[x,y,z1],orientation=+1,refface=refZ);
     // medit(" --- ", Thx0,Thx1,Thy0,Thy1,Thz0,Thz1);
    ThHex = Thx0+Thx1+Thy0+Thy1+Thz0+Thz1;
  }
  mesh3 Thsph;  //
  {
    mesh  Th=square(10,20,[x⋆pi-pi/2,2⋆y⋆pi]);   // ]-pi
 2,frac-pi2[×]0, 2π[
     // a paratrization of a sphere
    func f1 =cos(x)⋆cos(y);
    func f2 =cos(x)⋆sin(y);
    func f3 = sin(x);
     // de partiel derivatrive of the parametrization DF
    func f1x=sin(x)⋆cos(y);
    func f1y=-cos(x)⋆sin(y);
    func f2x=-sin(x)⋆sin(y);
    func f2y=cos(x)⋆cos(y);
    func f3x=cos(x);
    func f3y=0;
     // M = DFtDF
    func m11=f1x^2+f2x^2+f3x^2;
    func m21=f1x⋆f1y+f2x⋆f2y+f3x⋆f3y;
    func m22=f1y^2+f2y^2+f3y^2;
    func perio=[[4,y],[2,y],[1,x],[3,x]];   // to store the periodic condition
     // the intial mesh
    savemesh(Th,"sphere",[f1,f2,f3]);
    real R=0.5,hh=0.1/R; // hh taille du maille sur la shere unite.
    real vv= 1/square(hh);
    verbosity=2;
    Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,inquire=1,periodic=perio);
    plot(Th,wait=1);
    Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,periodic=perio);
    plot(Th,wait=1);
    Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,periodic=perio);
    plot(Th,wait=1);
    Th=adaptmesh(Th,m11⋆vv,m21⋆vv,m22⋆vv,IsMetric=1,periodic=perio);
    Thsph = movemesh23(Th,transfo=[f1⋆R,f2⋆R,f3⋆R],orientation=-1);
  }
   // //////////////////////////////
  mesh3 ThS = ThHex+Thsph;  // "gluing" surface meshs to total boundary meshes
  medit("Bounday mesh",ThS,wait=1);
   // build a mesh of a axis parallel box with TetGen
  real[int] domaine = [0,0,0,1,volumetet,0,0,0.7,2,volumetet];
  Th = tetg(ThS,switch="pqaAAYYQ",nbofregions=2,regionlist=domaine);
   // Tetrahelize the interior of the cube with tetgen
  medit("tetg",Th);
  savemesh(Th,"Th-hex-sph.mesh");
}  // to clean all intermedial variable

5.12 Write solution at the format .sol and .solb

With the keyword savesol, we can store a scalar functions, a scalar FE functions, a vector fields, a vector FE fields, a symmetric tensor and a symmetric FE tensor.. Such format is used in medit.

extension file .sol The first two lines of the file are

The following fields are begin with one of the following keyword: SolAtVertices, SolAtEdges, SolAtTriangles, SolAtQuadrilaterals, SolAtTetrahedra, SolAtPentahedra, SolAtHexahedra.

In each field, we give then in the next line the number of elements in solutions (SolAtVertices: number of vertices, SolAtTriangles: number of triangles, ...). In other lines, we give the number of solutions , the type of solution (1: scalar, 2: vector, 3: symmetric tensor). And finaly, we give the solutions at the differents elements.

The file must be finished with the keyword End.

The real element of symmetric tensor

       (   3d    3d     3d)
       |||ST xx  STxy  ST xz |||                ( ST2d  ST 2d)
ST3d = |||||ST 3ydx  ST3ydy  ST 3yzd|||||         ST 2d =     x2xd     xy2d               (5.4 )
       (ST 3zdx  ST3zdy  ST 3zzd)                  STyx  ST yy
stored in the extension .sol are respectively STxx3d,STyx3d,STyy3d,STzx3d,STzy3d,STzz3d and STxx2d,STyx2d,STyy2d

An example of field with the keyword SolAtTetrahedra:

where

This field is express with the notation of Section ??. The format .solb is the format .sol writing in binary.

A real scalar functions f1, a vector fields Φ = [Φ1,Φ2,Φ3] and a symmetric tensor ST3d (5.4) at the vertices of the three dimensional mesh Th3 is stored in the file ”f1PhiTh3.sol” using


savesol("f1PhiST3dTh3.sol",Th3, f 1, [Φ1, Φ2, Φ3], VV3, order=1);

where VV3 = [STxx3d,STyx3d,STyy3d,STzx3d,STzy3d,STzz3d]. For a two dimension mesh Th, A real scalar functions f2, a vector fields Ψ = [Ψ1,Ψ2] and a symmetric tensor ST2d (5.4) at triangles is stored in the file ”f2PsiST2dTh3.solb” using


savesol("f2PsiST2dTh3.solb",Th, f 2, [Ψ1, Ψ2], VV2, order=0);

where VV2 = [STxx2d,STyx2d,STyy2d] The arguments of savesolfunctions are the name of a file, a mesh and solutions. These arguments must be given in this order.

The parmameters of this keyword are

order =
0 is the solution is given at the center of gravity of elements. 1 is the solution is given at the vertices of elements.

In the file, solutions are stored in this order : scalar solutions, vector solutions and finaly symmetric tensor solutions.

5.13 Call medit with the keyword medit

The keyword medit allows to dipslay a mesh or mesh and solution using medit. Medit is a freeware display package by Pascal Frey using OpenGL. To use this command we need to install medit.

A vizualisation with medit of scalar solutions f1 and f2 at vertices of the mesh Th is obtained using


medit("sol1 sol2",Th, f 1, f 2, order=1);

The first plot named “sol1” display f1. The second plot names “sol2” display f2.

The arguments of function meditare the name of the differents scenes (separated by a space) of medit, a mesh and solutions. Each solution is associated with one scene. The scalar, vector and symmetric tensor solutions are given as in the keyword savesol.

The parmameters of this command line are

order =
0 is the solution is given at the center of gravity of elements. 1 is the solution is given at the vertices of elements.
meditff =
set the name of execute command of medit. By default, this string is medit.
save =
set the name of a file .sol or .solb to save solutions.

This commandline allows also to represent two differents meshes with solutions in the same windows. The nature of solutions must be the same. Hence, we can vizualize in the same window the different domains in a domain decomposition method. A vizualisation with medit of scalar solutions h1 and h2 at vertices of the mesh Th1 and Th2 respectively are obtained using


medit("sol2domain",Th1, h1, Th2, h2, order=1);

Example 5.23 (meditddm.edp)


// meditddm.edp
load "medit"
// Initial Problem:
// Resolution of the following EDP:
// -Δu = f on Ω = {(x,y)|1 sqrt(x2 + y2) 2}
// -Δu = f 1 on Ω1 = {(x,y)|0.5 sqrt(x2 + y2) 1.}
// u = 1 on Γ + Null Neumman condition on Γ1 and on Γ2
// We find the solution u in solving two EDP defined on domain Ω and Ω1
// This solution is vizualize with medit
verbosity=3;
border Gamma(t=0,2⋆pi){x=cos(t); y=sin(t); label=1;};
border Gamma1(t=0,2⋆pi){x=2⋆cos(t); y=2⋆sin(t); label=2;};
border Gamma2(t=0,2⋆pi){x=0.5⋆cos(t); y=0.5⋆sin(t); label=3;};
// construction of mesh of domain Ω
mesh Th=buildmesh(Gamma1(40)+Gamma(-40));
fespace Vh(Th,P2);
func f=sqrt(x⋆x+y⋆y);
Vh us,v;
macro Grad2(us) [dx(us),dy(us)]   // EOM
problem Lap2dOmega(us,v,init=false)=int2d(Th)(Grad2(v)' ⋆Grad2(us)) - int2d(Th)(f⋆v)+on(Gamma,us=1) ;
// Resolution of EDP defined on the domain Ω
// -Δu = f on Ω
// u = 1 on Γ
// + Null Neumann condition on Γ1
Lap2dOmega;
// construction of mesh of domain Ω1
mesh Th1=buildmesh(Gamma(40)+Gamma2(-40));
fespace Vh1(Th1,P2);
func f1=10⋆sqrt(x⋆x+y⋆y);
Vh1 u1,v1;
macro Grad21(u1) [dx(u1),dy(u1)]   // EOM
problem Lap2dOmega1(u1,v1,init=false)=int2d(Th1)(Grad21(v1)' ⋆Grad21(u1)) - int2d(Th1)(f1⋆v1)+on(Gamma,u1=1) ;
// Resolution of EDP defined on the domain Ω1
// -Δu = f 1 on Ω1
// u = 1 on Γ
// + Null Neumann condition on Γ2
Lap2dOmega1;
// vizualisation of solution of the initial problem
medit("solution",Th,us,Th1,u1,order=1,save="testsavemedit.solb");

Example 5.24 (StockesUzawa.edp)


// signe of pressure is correct
assert(version>1.18);
real s0=clock();
mesh Th=square(10,10);
fespace Xh(Th,P2),Mh(Th,P1);
Xh u1,u2,v1,v2;
Mh p,q,ppp;
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(1,2,4,u1=0)  +  on(3,u1=1) ;
Xh bc1; bc1[] = a(0,Xh);
Xh b;
matrix A= a(Xh,Xh,solver=CG);
matrix Bx= bx(Xh,Mh);
matrix By= by(Xh,Mh);
Xh bcx=1,bcy=0;
func real[int] divup(real[int] & pp)
{
  int verb=verbosity;
   verbosity=0;
   b[]  = Bx'⋆pp; b[] += bc1[] .⋆bcx[];
   u1[] = A^-1⋆b[];
   b[]  = By'⋆pp; b[] += bc1[] .⋆bcy[];
   u2[] = A^-1⋆b[];
   ppp[] =   Bx⋆u1[];
   ppp[] +=  By⋆u2[];
   verbosity=verb;
   return ppp[] ;
};
p=0;q=0;u1=0;v1=0;
LinearCG(divup,p[],q[],eps=1.e-6,nbiter=50);
divup(p[]);
plot([u1,u2],p,wait=1,value=true,coef=0.1);
medit(‘‘velocity pressure'',Th,[u1,u2],p,order=1);