xref: /petsc/src/snes/tutorials/ex9.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static const char help[] = "Solves obstacle problem in 2D as a variational inequality\n\
2 or nonlinear complementarity problem.  This is a form of the Laplace equation in\n\
3 which the solution u is constrained to be above a given function psi.  In the\n\
4 problem here an exact solution is known.\n";
5 
6 /*  On a square S = {-2<x<2,-2<y<2}, the PDE
7     u_{xx} + u_{yy} = 0
8 is solved on the set where membrane is above obstacle (u(x,y) >= psi(x,y)).
9 Here psi is the upper hemisphere of the unit ball.  On the boundary of S
10 we have Dirichlet boundary conditions from the exact solution.  Uses centered
11 FD scheme.  This example contributed by Ed Bueler.
12 
13 Example usage:
14   * get help:
15     ./ex9 -help
16   * monitor run:
17     ./ex9 -da_refine 2 -snes_vi_monitor
18   * use other SNESVI type (default is SNESVINEWTONRSLS):
19     ./ex9 -da_refine 2 -snes_vi_monitor -snes_type vinewtonssls
20   * use FD evaluation of Jacobian by coloring, instead of analytical:
21     ./ex9 -da_refine 2 -snes_fd_color
22   * X windows visualizations:
23     ./ex9 -snes_monitor_solution draw -draw_pause 1 -da_refine 4
24     ./ex9 -snes_vi_monitor_residual -draw_pause 1 -da_refine 4
25   * full-cycle multigrid:
26     ./ex9 -snes_converged_reason -snes_grid_sequence 4 -pc_type mg
27   * serial convergence evidence:
28     for M in 3 4 5 6 7; do ./ex9 -snes_grid_sequence $M -pc_type mg; done
29   * FIXME sporadic parallel bug:
30     mpiexec -n 4 ./ex9 -snes_converged_reason -snes_grid_sequence 4 -pc_type mg
31 */
32 
33 #include <petsc.h>
34 
35 /* z = psi(x,y) is the hemispherical obstacle, but made C^1 with "skirt" at r=r0 */
36 PetscReal psi(PetscReal x, PetscReal y)
37 {
38     const PetscReal  r = x * x + y * y,r0 = 0.9,psi0 = PetscSqrtReal(1.0 - r0*r0),dpsi0 = - r0 / psi0;
39     if (r <= r0) {
40       return PetscSqrtReal(1.0 - r);
41     } else {
42       return psi0 + dpsi0 * (r - r0);
43     }
44 }
45 
46 /*  This exact solution solves a 1D radial free-boundary problem for the
47 Laplace equation, on the interval 0 < r < 2, with above obstacle psi(x,y).
48 The Laplace equation applies where u(r) > psi(r),
49     u''(r) + r^-1 u'(r) = 0
50 with boundary conditions including free b.c.s at an unknown location r = a:
51     u(a) = psi(a),  u'(a) = psi'(a),  u(2) = 0
52 The solution is  u(r) = - A log(r) + B   on  r > a.  The boundary conditions
53 can then be reduced to a root-finding problem for a:
54     a^2 (log(2) - log(a)) = 1 - a^2
55 The solution is a = 0.697965148223374 (giving residual 1.5e-15).  Then
56 A = a^2*(1-a^2)^(-0.5) and B = A*log(2) are as given below in the code.  */
57 PetscReal u_exact(PetscReal x, PetscReal y)
58 {
59     const PetscReal afree = 0.697965148223374,
60                     A     = 0.680259411891719,
61                     B     = 0.471519893402112;
62     PetscReal  r;
63     r = PetscSqrtReal(x * x + y * y);
64     return (r <= afree) ? psi(x,y)  /* active set; on the obstacle */
65                         : - A * PetscLogReal(r) + B; /* solves laplace eqn */
66 }
67 
68 extern PetscErrorCode FormExactSolution(DMDALocalInfo*,Vec);
69 extern PetscErrorCode FormBounds(SNES,Vec,Vec);
70 extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscReal**,PetscReal**,void*);
71 extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscReal**,Mat,Mat,void*);
72 
73 int main(int argc,char **argv)
74 {
75   SNES                snes;
76   DM                  da, da_after;
77   Vec                 u, u_exact;
78   DMDALocalInfo       info;
79   PetscReal           error1,errorinf;
80 
81   PetscFunctionBeginUser;
82   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
83 
84   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,5,5, /* 5x5 coarse grid; override with -da_grid_x,_y */
85                          PETSC_DECIDE,PETSC_DECIDE, 1,1,  /* dof=1 and s = 1 (stencil extends out one cell) */
86                          NULL,NULL,&da));
87   PetscCall(DMSetFromOptions(da));
88   PetscCall(DMSetUp(da));
89   PetscCall(DMDASetUniformCoordinates(da,-2.0,2.0,-2.0,2.0,0.0,1.0));
90 
91   PetscCall(DMCreateGlobalVector(da,&u));
92   PetscCall(VecSet(u,0.0));
93 
94   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
95   PetscCall(SNESSetDM(snes,da));
96   PetscCall(SNESSetType(snes,SNESVINEWTONRSLS));
97   PetscCall(SNESVISetComputeVariableBounds(snes,&FormBounds));
98   PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,NULL));
99   PetscCall(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,NULL));
100   PetscCall(SNESSetFromOptions(snes));
101 
102   /* solve nonlinear system */
103   PetscCall(SNESSolve(snes,NULL,u));
104   PetscCall(VecDestroy(&u));
105   PetscCall(DMDestroy(&da));
106   /* DMDA after solve may be different, e.g. with -snes_grid_sequence */
107   PetscCall(SNESGetDM(snes,&da_after));
108   PetscCall(SNESGetSolution(snes,&u)); /* do not destroy u */
109   PetscCall(DMDAGetLocalInfo(da_after,&info));
110   PetscCall(VecDuplicate(u,&u_exact));
111   PetscCall(FormExactSolution(&info,u_exact));
112   PetscCall(VecAXPY(u,-1.0,u_exact)); /* u <-- u - u_exact */
113   PetscCall(VecNorm(u,NORM_1,&error1));
114   error1 /= (PetscReal)info.mx * (PetscReal)info.my; /* average error */
115   PetscCall(VecNorm(u,NORM_INFINITY,&errorinf));
116   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"errors on %" PetscInt_FMT " x %" PetscInt_FMT " grid:  av |u-uexact|  = %.3e,  |u-uexact|_inf = %.3e\n",info.mx,info.my,(double)error1,(double)errorinf));
117   PetscCall(VecDestroy(&u_exact));
118   PetscCall(SNESDestroy(&snes));
119   PetscCall(DMDestroy(&da));
120   PetscCall(PetscFinalize());
121   return 0;
122 }
123 
124 PetscErrorCode FormExactSolution(DMDALocalInfo *info, Vec u)
125 {
126   PetscInt       i,j;
127   PetscReal      **au, dx, dy, x, y;
128   dx = 4.0 / (PetscReal)(info->mx-1);
129   dy = 4.0 / (PetscReal)(info->my-1);
130   PetscCall(DMDAVecGetArray(info->da, u, &au));
131   for (j=info->ys; j<info->ys+info->ym; j++) {
132     y = -2.0 + j * dy;
133     for (i=info->xs; i<info->xs+info->xm; i++) {
134       x = -2.0 + i * dx;
135       au[j][i] = u_exact(x,y);
136     }
137   }
138   PetscCall(DMDAVecRestoreArray(info->da, u, &au));
139   return 0;
140 }
141 
142 PetscErrorCode FormBounds(SNES snes, Vec Xl, Vec Xu)
143 {
144   DM             da;
145   DMDALocalInfo  info;
146   PetscInt       i, j;
147   PetscReal      **aXl, dx, dy, x, y;
148 
149   PetscCall(SNESGetDM(snes,&da));
150   PetscCall(DMDAGetLocalInfo(da,&info));
151   dx = 4.0 / (PetscReal)(info.mx-1);
152   dy = 4.0 / (PetscReal)(info.my-1);
153   PetscCall(DMDAVecGetArray(da, Xl, &aXl));
154   for (j=info.ys; j<info.ys+info.ym; j++) {
155     y = -2.0 + j * dy;
156     for (i=info.xs; i<info.xs+info.xm; i++) {
157       x = -2.0 + i * dx;
158       aXl[j][i] = psi(x,y);
159     }
160   }
161   PetscCall(DMDAVecRestoreArray(da, Xl, &aXl));
162   PetscCall(VecSet(Xu,PETSC_INFINITY));
163   return 0;
164 }
165 
166 PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **au, PetscScalar **af, void *user)
167 {
168   PetscInt       i,j;
169   PetscReal      dx,dy,x,y,ue,un,us,uw;
170 
171   PetscFunctionBeginUser;
172   dx = 4.0 / (PetscReal)(info->mx-1);
173   dy = 4.0 / (PetscReal)(info->my-1);
174   for (j=info->ys; j<info->ys+info->ym; j++) {
175     y = -2.0 + j * dy;
176     for (i=info->xs; i<info->xs+info->xm; i++) {
177       x = -2.0 + i * dx;
178       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
179         af[j][i] = 4.0 * (au[j][i] - u_exact(x,y));
180       } else {
181         uw = (i-1 == 0)          ? u_exact(x-dx,y) : au[j][i-1];
182         ue = (i+1 == info->mx-1) ? u_exact(x+dx,y) : au[j][i+1];
183         us = (j-1 == 0)          ? u_exact(x,y-dy) : au[j-1][i];
184         un = (j+1 == info->my-1) ? u_exact(x,y+dy) : au[j+1][i];
185         af[j][i] = - (dy/dx) * (uw - 2.0 * au[j][i] + ue) - (dx/dy) * (us - 2.0 * au[j][i] + un);
186       }
187     }
188   }
189   PetscCall(PetscLogFlops(12.0*info->ym*info->xm));
190   PetscFunctionReturn(0);
191 }
192 
193 PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **au, Mat A, Mat jac, void *user)
194 {
195   PetscInt       i,j,n;
196   MatStencil     col[5],row;
197   PetscReal      v[5],dx,dy,oxx,oyy;
198 
199   PetscFunctionBeginUser;
200   dx  = 4.0 / (PetscReal)(info->mx-1);
201   dy  = 4.0 / (PetscReal)(info->my-1);
202   oxx = dy / dx;
203   oyy = dx / dy;
204   for (j=info->ys; j<info->ys+info->ym; j++) {
205     for (i=info->xs; i<info->xs+info->xm; i++) {
206       row.j = j; row.i = i;
207       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { /* boundary */
208         v[0] = 4.0;
209         PetscCall(MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES));
210       } else { /* interior grid points */
211         v[0] = 2.0 * (oxx + oyy);  col[0].j = j;  col[0].i = i;
212         n = 1;
213         if (i-1 > 0) {
214           v[n] = -oxx;  col[n].j = j;  col[n++].i = i-1;
215         }
216         if (i+1 < info->mx-1) {
217           v[n] = -oxx;  col[n].j = j;  col[n++].i = i+1;
218         }
219         if (j-1 > 0) {
220           v[n] = -oyy;  col[n].j = j-1;  col[n++].i = i;
221         }
222         if (j+1 < info->my-1) {
223           v[n] = -oyy;  col[n].j = j+1;  col[n++].i = i;
224         }
225         PetscCall(MatSetValuesStencil(jac,1,&row,n,col,v,INSERT_VALUES));
226       }
227     }
228   }
229 
230   /* Assemble matrix, using the 2-step process: */
231   PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
232   PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
233   if (A != jac) {
234     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
235     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
236   }
237   PetscCall(PetscLogFlops(2.0*info->ym*info->xm));
238   PetscFunctionReturn(0);
239 }
240 
241 /*TEST
242 
243    build:
244       requires: !complex
245 
246    test:
247       suffix: 1
248       requires: !single
249       nsize: 1
250       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonrsls
251 
252    test:
253       suffix: 2
254       requires: !single
255       nsize: 2
256       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonssls
257 
258    test:
259       suffix: 3
260       requires: !single
261       nsize: 2
262       args: -snes_grid_sequence 2 -snes_vi_monitor -snes_type vinewtonrsls
263 
264    test:
265       suffix: mg
266       requires: !single
267       nsize: 4
268       args: -snes_grid_sequence 3 -snes_converged_reason -pc_type mg
269 
270    test:
271       suffix: 4
272       nsize: 1
273       args: -mat_is_symmetric
274 
275    test:
276       suffix: 5
277       nsize: 1
278       args: -ksp_converged_reason -snes_fd_color
279 
280    test:
281       suffix: 6
282       requires: !single
283       nsize: 2
284       args: -snes_grid_sequence 2 -pc_type mg -snes_monitor_short -ksp_converged_reason
285 
286    test:
287       suffix: 7
288       nsize: 2
289       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type multiplicative -snes_composite_sneses vinewtonrsls,vinewtonssls -sub_0_snes_vi_monitor -sub_1_snes_vi_monitor
290       TODO: fix nasty memory leak in SNESCOMPOSITE
291 
292    test:
293       suffix: 8
294       nsize: 2
295       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additive -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
296       TODO: fix nasty memory leak in SNESCOMPOSITE
297 
298    test:
299       suffix: 9
300       nsize: 2
301       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additiveoptimal -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
302       TODO: fix nasty memory leak in SNESCOMPOSITE
303 
304 TEST*/
305