Actual source code: ex9.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  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";

  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.

 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: */

 33:  #include <petsc.h>

 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:     const PetscReal  r = x * x + y * y,
 38:                      r0 = 0.9,
 39:                      psi0 = PetscSqrtReal(1.0 - r0*r0),
 40:                      dpsi0 = - r0 / psi0;
 41:     if (r <= r0) {
 42:         return PetscSqrtReal(1.0 - r);
 43:     } else {
 44:         return psi0 + dpsi0 * (r - r0);
 45:     }
 46: }

 48: /*  This exact solution solves a 1D radial free-boundary problem for the
 49: Laplace equation, on the interval 0 < r < 2, with above obstacle psi(x,y).
 50: The Laplace equation applies where u(r) > psi(r),
 51:     u''(r) + r^-1 u'(r) = 0
 52: with boundary conditions including free b.c.s at an unknown location r = a:
 53:     u(a) = psi(a),  u'(a) = psi'(a),  u(2) = 0
 54: The solution is  u(r) = - A log(r) + B   on  r > a.  The boundary conditions
 55: can then be reduced to a root-finding problem for a:
 56:     a^2 (log(2) - log(a)) = 1 - a^2
 57: The solution is a = 0.697965148223374 (giving residual 1.5e-15).  Then
 58: A = a^2*(1-a^2)^(-0.5) and B = A*log(2) are as given below in the code.  */
 59: PetscReal u_exact(PetscReal x, PetscReal y) {
 60:     const PetscReal afree = 0.697965148223374,
 61:                     A     = 0.680259411891719,
 62:                     B     = 0.471519893402112;
 63:     PetscReal  r;
 64:     r = PetscSqrtReal(x * x + y * y);
 65:     return (r <= afree) ? psi(x,y)  /* active set; on the obstacle */
 66:                         : - A * PetscLogReal(r) + B; /* solves laplace eqn */
 67: }

 69: extern PetscErrorCode FormExactSolution(DMDALocalInfo*,Vec);
 70: extern PetscErrorCode FormBounds(SNES,Vec,Vec);
 71: extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscReal**,PetscReal**,void*);
 72: extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscReal**,Mat,Mat,void*);

 74: int main(int argc,char **argv) {
 75:   PetscErrorCode      ierr;
 76:   SNES                snes;
 77:   DM                  da, da_after;
 78:   Vec                 u, u_exact;
 79:   DMDALocalInfo       info;
 80:   PetscReal           error1,errorinf;

 82:   PetscInitialize(&argc,&argv,(char*)0,help);

 84:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,
 85:                       DMDA_STENCIL_STAR,5,5, /* 5x5 coarse grid; override with -da_grid_x,_y */
 86:                       PETSC_DECIDE,PETSC_DECIDE,
 87:                       1,1,  /* dof=1 and s = 1 (stencil extends out one cell) */
 88:                       NULL,NULL,&da);
 89:   DMSetFromOptions(da);
 90:   DMSetUp(da);
 91:   DMDASetUniformCoordinates(da,-2.0,2.0,-2.0,2.0,0.0,1.0);

 93:   DMCreateGlobalVector(da,&u);
 94:   VecSet(u,0.0);

 96:   SNESCreate(PETSC_COMM_WORLD,&snes);
 97:   SNESSetDM(snes,da);
 98:   SNESSetType(snes,SNESVINEWTONRSLS);
 99:   SNESVISetComputeVariableBounds(snes,&FormBounds);
100:   DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,NULL);
101:   DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,NULL);
102:   SNESSetFromOptions(snes);

104:   /* solve nonlinear system */
105:   SNESSolve(snes,NULL,u);
106:   VecDestroy(&u);
107:   DMDestroy(&da);
108:   /* DMDA after solve may be different, e.g. with -snes_grid_sequence */
109:   SNESGetDM(snes,&da_after);
110:   SNESGetSolution(snes,&u);  /* do not destroy u */
111:   DMDAGetLocalInfo(da_after,&info);
112:   VecDuplicate(u,&u_exact);
113:   FormExactSolution(&info,u_exact);
114:   VecAXPY(u,-1.0,u_exact); /* u <-- u - u_exact */
115:   VecNorm(u,NORM_1,&error1);
116:   error1 /= (PetscReal)info.mx * (PetscReal)info.my; /* average error */
117:   VecNorm(u,NORM_INFINITY,&errorinf);
118:   PetscPrintf(PETSC_COMM_WORLD,
119:              "errors on %D x %D grid:  av |u-uexact|  = %.3e,  |u-uexact|_inf = %.3e\n",
120:              info.mx,info.my,(double)error1,(double)errorinf);

122:   VecDestroy(&u_exact);
123:   SNESDestroy(&snes);
124:   DMDestroy(&da);
125:   PetscFinalize();
126:   return ierr;
127: }

129: PetscErrorCode FormExactSolution(DMDALocalInfo *info, Vec u) {
131:   PetscInt       i,j;
132:   PetscReal      **au, dx, dy, x, y;
133:   dx = 4.0 / (PetscReal)(info->mx-1);
134:   dy = 4.0 / (PetscReal)(info->my-1);
135:   DMDAVecGetArray(info->da, u, &au);
136:   for (j=info->ys; j<info->ys+info->ym; j++) {
137:     y = -2.0 + j * dy;
138:     for (i=info->xs; i<info->xs+info->xm; i++) {
139:       x = -2.0 + i * dx;
140:       au[j][i] = u_exact(x,y);
141:     }
142:   }
143:   DMDAVecRestoreArray(info->da, u, &au);
144:   return 0;
145: }

147: PetscErrorCode FormBounds(SNES snes, Vec Xl, Vec Xu) {
149:   DM             da;
150:   DMDALocalInfo  info;
151:   PetscInt       i, j;
152:   PetscReal      **aXl, dx, dy, x, y;
153:   SNESGetDM(snes,&da);
154:   DMDAGetLocalInfo(da,&info);
155:   dx = 4.0 / (PetscReal)(info.mx-1);
156:   dy = 4.0 / (PetscReal)(info.my-1);
157:   DMDAVecGetArray(da, Xl, &aXl);
158:   for (j=info.ys; j<info.ys+info.ym; j++) {
159:     y = -2.0 + j * dy;
160:     for (i=info.xs; i<info.xs+info.xm; i++) {
161:       x = -2.0 + i * dx;
162:       aXl[j][i] = psi(x,y);
163:     }
164:   }
165:   DMDAVecRestoreArray(da, Xl, &aXl);
166:   VecSet(Xu,PETSC_INFINITY);
167:   return 0;
168: }

170: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **au, PetscScalar **af, void *user) {
172:   PetscInt       i,j;
173:   PetscReal      dx,dy,x,y,ue,un,us,uw;

176:   dx = 4.0 / (PetscReal)(info->mx-1);
177:   dy = 4.0 / (PetscReal)(info->my-1);
178:   for (j=info->ys; j<info->ys+info->ym; j++) {
179:     y = -2.0 + j * dy;
180:     for (i=info->xs; i<info->xs+info->xm; i++) {
181:       x = -2.0 + i * dx;
182:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
183:         af[j][i] = 4.0 * (au[j][i] - u_exact(x,y));
184:       } else {
185:         uw = (i-1 == 0)          ? u_exact(x-dx,y) : au[j][i-1];
186:         ue = (i+1 == info->mx-1) ? u_exact(x+dx,y) : au[j][i+1];
187:         us = (j-1 == 0)          ? u_exact(x,y-dy) : au[j-1][i];
188:         un = (j+1 == info->my-1) ? u_exact(x,y+dy) : au[j+1][i];
189:         af[j][i] = - (dy/dx) * (uw - 2.0 * au[j][i] + ue)
190:                    - (dx/dy) * (us - 2.0 * au[j][i] + un);
191:       }
192:     }
193:   }
194:   PetscLogFlops(12.0*info->ym*info->xm);
195:   return(0);
196: }

198: PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **au, Mat A, Mat jac, void *user) {
200:   PetscInt       i,j,n;
201:   MatStencil     col[5],row;
202:   PetscReal      v[5],dx,dy,oxx,oyy;

205:   dx  = 4.0 / (PetscReal)(info->mx-1);
206:   dy  = 4.0 / (PetscReal)(info->my-1);
207:   oxx = dy / dx;
208:   oyy = dx / dy;
209:   for (j=info->ys; j<info->ys+info->ym; j++) {
210:     for (i=info->xs; i<info->xs+info->xm; i++) {
211:       row.j = j; row.i = i;
212:       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { /* boundary */
213:         v[0] = 4.0;
214:         MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
215:       } else { /* interior grid points */
216:         v[0] = 2.0 * (oxx + oyy);  col[0].j = j;  col[0].i = i;
217:         n = 1;
218:         if (i-1 > 0) {
219:           v[n] = -oxx;  col[n].j = j;  col[n++].i = i-1;
220:         }
221:         if (i+1 < info->mx-1) {
222:           v[n] = -oxx;  col[n].j = j;  col[n++].i = i+1;
223:         }
224:         if (j-1 > 0) {
225:           v[n] = -oyy;  col[n].j = j-1;  col[n++].i = i;
226:         }
227:         if (j+1 < info->my-1) {
228:           v[n] = -oyy;  col[n].j = j+1;  col[n++].i = i;
229:         }
230:         MatSetValuesStencil(jac,1,&row,n,col,v,INSERT_VALUES);
231:       }
232:     }
233:   }

235:   /* Assemble matrix, using the 2-step process: */
236:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
237:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
238:   if (A != jac) {
239:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
240:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
241:   }
242:   PetscLogFlops(2.0*info->ym*info->xm);
243:   return(0);
244: }



248: /*TEST

250:    build:
251:       requires: !complex

253:    test:
254:       suffix: 1
255:       requires: !single
256:       nsize: 1
257:       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonrsls

259:    test:
260:       suffix: 2
261:       requires: !single
262:       nsize: 2
263:       args: -da_refine 1 -snes_monitor_short -snes_type vinewtonssls

265:    test:
266:       suffix: 3
267:       requires: !single
268:       nsize: 2
269:       args: -snes_grid_sequence 2 -snes_vi_monitor -snes_type vinewtonrsls

271:    test:
272:       suffix: 4
273:       nsize: 1
274:       args: -mat_is_symmetric

276:    test:
277:       suffix: 5
278:       nsize: 1
279:       args: -ksp_converged_reason -snes_fd_color

281:    test:
282:       suffix: 6
283:       requires: !single
284:       nsize: 2
285:       args: -snes_grid_sequence 2 -pc_type mg -snes_monitor_short -ksp_converged_reason

287:    test:
288:       suffix: 7
289:       nsize: 2
290:       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
291:       TODO: fix nasty memory leak in SNESCOMPOSITE

293:    test:
294:       suffix: 8
295:       nsize: 2
296:       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additive -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
297:       TODO: fix nasty memory leak in SNESCOMPOSITE

299:    test:
300:       suffix: 9
301:       nsize: 2
302:       args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additiveoptimal -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor
303:       TODO: fix nasty memory leak in SNESCOMPOSITE

305: TEST*/