-
Notifications
You must be signed in to change notification settings - Fork 10
/
Copy pathDDM-Schwarz-Lame-2d.edp
152 lines (109 loc) · 3.93 KB
/
DDM-Schwarz-Lame-2d.edp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
// NBPROC 10
// PARAM -ns -gmres 0 -k 3
// ff-mpirun -np 4 DDM-Schwarz-Lap-3d.edp -glut ffglut -n 11 -k 1 -d 1 -ns -gmres 1
/*
a first true parallele example fisrt freefem++
Ok up to 200 proc for a Poisson equation..
See the Doc for full explaiantion
F Hecht Dec. 2010.
-------------------
usage :
ff-mpirun [mpi parameter] MPIGMRES2d.edp [-glut ffglut] [-n N] [-k K] [-d D] [-ns] [-gmres [0|1]
argument:
-glut ffglut : to see graphicaly the process
-n N: set the mesh cube split NxNxN
-d D: set debug flag D must be one for mpiplot
-k K: to refined by K all elemnt
-ns: reomove script dump
-gmres 0 : use iterative schwarz algo.
1 : Algo GMRES on residu of schwarz algo.
2 : DDM GMRES
3 : DDM GMRES with coarse grid preconditionner (Good one)
*/
load "MPICG" load "medit" load "metis"
include "getARGV.idp"
include "DDM-Schwarz-macro.idp"
include "AddLayer2d.idp"
include "DDM-funcs-v2.idp"
searchMethod=1; // more safe seach algo (warning can be very expensive in case lot of ouside point)
// 0 by default
assert(version >=3.11);
real[int] ttt(10);int ittt=0;
macro settt {ttt[ittt++]=mpiWtime();}//
verbosity=getARGV("-vv",0);
int vdebug=getARGV("-d",1);
int ksplit=getARGV("-k",1);
int nloc = getARGV("-n",10);
string sff=getARGV("-p,","");
int gmres=getARGV("-gmres",3);
int nC = getARGV("-N" ,max(nloc/10,5));
int sizeoverlaps=1; // size of overlap
bool RAS=1; // select kind of of $\pi_i$
if(mpirank==0 && verbosity)
cout << " vdebug: " << vdebug << " kspilt "<< ksplit << " nloc "<< nloc << " sff "<< sff <<"."<< endl;
string sPk="P2-Lame-2gd";
func Pk=[P2,P2];
func bool plotMPIall(mesh &Th,real[int] & u,string cm)
{ PLOTMPIALLU(mesh,Pk, defPk2, Th, u, allu1, { cmm=cm,nbiso=10,fill=1,dim=3,value=1}); return 1;}
mpiComm comm(mpiCommWorld,0,0);// trick : make a no split mpiWorld
int ipart= mpiRank(comm); // current partition number
if(ipart==0) cout << " Final N=" << ksplit*nloc << " nloc =" << nloc << " split =" << ksplit << endl;
int[int] l111=[2,2,2,1];
settt
mesh Thg=square(nloc*4,nloc,[x*4,y],label=l111);
mesh ThC=square(nC*4,nC,[x*4,y],label=l111);// Coarse Mesh
fespace VhC(ThC,[P1,P1]); // of the coarse problem..
BuildPartitioning(sizeoverlaps,mesh,Thg,Thi,aThij,RAS,pii,jpart,comm,vdebug)
if(ksplit>1)
{
for(int jp=0;jp<jpart.n;++jp)
aThij[jp] = trunc(aThij[jp],1,split=ksplit);
Thi = trunc(Thi,1,split=ksplit);
}
BuildTransferMat(ipart,mesh,Pk,2,[0,1],Thi,Whi,Whij,Thij,aThij,Usend,Vrecv,jpart,vdebug)
/* the definition of the Problem .... */
// the definition of the Problem ....
real E = 21.5e4;
real sigma = 0.29;
real mu = E/(2*(1+sigma));
real lambda = E*sigma/((1+sigma)*(1-2*sigma));
real gravity = -0.05;
real sqrt2=sqrt(2.);
macro epsilon(u1,u2) [dx(u1),dy(u2),(dy(u1)+dx(u2))/sqrt2] // EOM
macro div(u1,u2) ( dx(u1)+dy(u2) ) // EOM
varf vPb([u1,u2],[v1,v2])=
int2d(Thi)(
lambda*div(u1,u2)*div(v1,v2)
+2.*mu*( epsilon(u1,u2)'*epsilon(v1,v2) ) //')
)
+ int2d(Thi) (gravity*v2)
+ on(1,10,u1=0,u2=0)
;
varf vPbC([u1,u2],[v1,v2])=
int2d(ThC)(
lambda*div(u1,u2)*div(v1,v2)
+2.*mu*( epsilon(u1,u2)'*epsilon(v1,v2) ) //')
)
+ on(1,u1=0,u2=0)
;
varf vPbon10([u1,u2],[v1,v2])=on(10,u1=1,u2=1)+on(1,u1=0,u2=0);
varf vPBC(U,V)=on(1,U=0);
real[int] onG10 = vPbon10(0,Whi); // on 1
real[int] Bi=vPb(0,Whi);
matrix Ai = vPb(Whi,Whi,solver=sparsesolver);
DMMDeffuncAndGlobals(Lame2,comm,jpart,Whi,Vhc,2,Ai,vPbC,onG10,Pii,Usend,Vrecv,[0,1])
Lame2CheckUpdate();
Whi [u,u1],[v,v1];
u[]=vPBC(0,Whi,tgv=1);
real eps=1e-10;
Lame2DDMSolver(Bi,u,v,gmres,eps,vdebug)
real errg =1,umaxg;
{
real umax = u[].max,umaxg;
real[int] aa=[umax], bb(1);
mpiAllReduce(aa,bb,comm,mpiMAX);
errg=bb[0];
if(ipart==0)
cout << " umax global = " << bb[0] << " Wtime = " << (ttt[ittt-1]-ttt[ittt-2]) << " s " << " " << Lame2kiter << endl;
}
Lame2Saveff(sff,eps,ksplit,nloc,sizeoverlaps);