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 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 82 83 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 */ 84 PETSC_DECIDE,PETSC_DECIDE, 1,1, /* dof=1 and s = 1 (stencil extends out one cell) */ 85 NULL,NULL,&da)); 86 PetscCall(DMSetFromOptions(da)); 87 PetscCall(DMSetUp(da)); 88 PetscCall(DMDASetUniformCoordinates(da,-2.0,2.0,-2.0,2.0,0.0,1.0)); 89 90 PetscCall(DMCreateGlobalVector(da,&u)); 91 PetscCall(VecSet(u,0.0)); 92 93 PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 94 PetscCall(SNESSetDM(snes,da)); 95 PetscCall(SNESSetType(snes,SNESVINEWTONRSLS)); 96 PetscCall(SNESVISetComputeVariableBounds(snes,&FormBounds)); 97 PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,NULL)); 98 PetscCall(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,NULL)); 99 PetscCall(SNESSetFromOptions(snes)); 100 101 /* solve nonlinear system */ 102 PetscCall(SNESSolve(snes,NULL,u)); 103 PetscCall(VecDestroy(&u)); 104 PetscCall(DMDestroy(&da)); 105 /* DMDA after solve may be different, e.g. with -snes_grid_sequence */ 106 PetscCall(SNESGetDM(snes,&da_after)); 107 PetscCall(SNESGetSolution(snes,&u)); /* do not destroy u */ 108 PetscCall(DMDAGetLocalInfo(da_after,&info)); 109 PetscCall(VecDuplicate(u,&u_exact)); 110 PetscCall(FormExactSolution(&info,u_exact)); 111 PetscCall(VecAXPY(u,-1.0,u_exact)); /* u <-- u - u_exact */ 112 PetscCall(VecNorm(u,NORM_1,&error1)); 113 error1 /= (PetscReal)info.mx * (PetscReal)info.my; /* average error */ 114 PetscCall(VecNorm(u,NORM_INFINITY,&errorinf)); 115 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)); 116 PetscCall(VecDestroy(&u_exact)); 117 PetscCall(SNESDestroy(&snes)); 118 PetscCall(DMDestroy(&da)); 119 PetscCall(PetscFinalize()); 120 return 0; 121 } 122 123 PetscErrorCode FormExactSolution(DMDALocalInfo *info, Vec u) 124 { 125 PetscInt i,j; 126 PetscReal **au, dx, dy, x, y; 127 dx = 4.0 / (PetscReal)(info->mx-1); 128 dy = 4.0 / (PetscReal)(info->my-1); 129 PetscCall(DMDAVecGetArray(info->da, u, &au)); 130 for (j=info->ys; j<info->ys+info->ym; j++) { 131 y = -2.0 + j * dy; 132 for (i=info->xs; i<info->xs+info->xm; i++) { 133 x = -2.0 + i * dx; 134 au[j][i] = u_exact(x,y); 135 } 136 } 137 PetscCall(DMDAVecRestoreArray(info->da, u, &au)); 138 return 0; 139 } 140 141 PetscErrorCode FormBounds(SNES snes, Vec Xl, Vec Xu) 142 { 143 DM da; 144 DMDALocalInfo info; 145 PetscInt i, j; 146 PetscReal **aXl, dx, dy, x, y; 147 148 PetscCall(SNESGetDM(snes,&da)); 149 PetscCall(DMDAGetLocalInfo(da,&info)); 150 dx = 4.0 / (PetscReal)(info.mx-1); 151 dy = 4.0 / (PetscReal)(info.my-1); 152 PetscCall(DMDAVecGetArray(da, Xl, &aXl)); 153 for (j=info.ys; j<info.ys+info.ym; j++) { 154 y = -2.0 + j * dy; 155 for (i=info.xs; i<info.xs+info.xm; i++) { 156 x = -2.0 + i * dx; 157 aXl[j][i] = psi(x,y); 158 } 159 } 160 PetscCall(DMDAVecRestoreArray(da, Xl, &aXl)); 161 PetscCall(VecSet(Xu,PETSC_INFINITY)); 162 return 0; 163 } 164 165 PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **au, PetscScalar **af, void *user) 166 { 167 PetscInt i,j; 168 PetscReal dx,dy,x,y,ue,un,us,uw; 169 170 PetscFunctionBeginUser; 171 dx = 4.0 / (PetscReal)(info->mx-1); 172 dy = 4.0 / (PetscReal)(info->my-1); 173 for (j=info->ys; j<info->ys+info->ym; j++) { 174 y = -2.0 + j * dy; 175 for (i=info->xs; i<info->xs+info->xm; i++) { 176 x = -2.0 + i * dx; 177 if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 178 af[j][i] = 4.0 * (au[j][i] - u_exact(x,y)); 179 } else { 180 uw = (i-1 == 0) ? u_exact(x-dx,y) : au[j][i-1]; 181 ue = (i+1 == info->mx-1) ? u_exact(x+dx,y) : au[j][i+1]; 182 us = (j-1 == 0) ? u_exact(x,y-dy) : au[j-1][i]; 183 un = (j+1 == info->my-1) ? u_exact(x,y+dy) : au[j+1][i]; 184 af[j][i] = - (dy/dx) * (uw - 2.0 * au[j][i] + ue) - (dx/dy) * (us - 2.0 * au[j][i] + un); 185 } 186 } 187 } 188 PetscCall(PetscLogFlops(12.0*info->ym*info->xm)); 189 PetscFunctionReturn(0); 190 } 191 192 PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **au, Mat A, Mat jac, void *user) 193 { 194 PetscInt i,j,n; 195 MatStencil col[5],row; 196 PetscReal v[5],dx,dy,oxx,oyy; 197 198 PetscFunctionBeginUser; 199 dx = 4.0 / (PetscReal)(info->mx-1); 200 dy = 4.0 / (PetscReal)(info->my-1); 201 oxx = dy / dx; 202 oyy = dx / dy; 203 for (j=info->ys; j<info->ys+info->ym; j++) { 204 for (i=info->xs; i<info->xs+info->xm; i++) { 205 row.j = j; row.i = i; 206 if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { /* boundary */ 207 v[0] = 4.0; 208 PetscCall(MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES)); 209 } else { /* interior grid points */ 210 v[0] = 2.0 * (oxx + oyy); col[0].j = j; col[0].i = i; 211 n = 1; 212 if (i-1 > 0) { 213 v[n] = -oxx; col[n].j = j; col[n++].i = i-1; 214 } 215 if (i+1 < info->mx-1) { 216 v[n] = -oxx; col[n].j = j; col[n++].i = i+1; 217 } 218 if (j-1 > 0) { 219 v[n] = -oyy; col[n].j = j-1; col[n++].i = i; 220 } 221 if (j+1 < info->my-1) { 222 v[n] = -oyy; col[n].j = j+1; col[n++].i = i; 223 } 224 PetscCall(MatSetValuesStencil(jac,1,&row,n,col,v,INSERT_VALUES)); 225 } 226 } 227 } 228 229 /* Assemble matrix, using the 2-step process: */ 230 PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 231 PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 232 if (A != jac) { 233 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 234 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 235 } 236 PetscCall(PetscLogFlops(2.0*info->ym*info->xm)); 237 PetscFunctionReturn(0); 238 } 239 240 /*TEST 241 242 build: 243 requires: !complex 244 245 test: 246 suffix: 1 247 requires: !single 248 nsize: 1 249 args: -da_refine 1 -snes_monitor_short -snes_type vinewtonrsls 250 251 test: 252 suffix: 2 253 requires: !single 254 nsize: 2 255 args: -da_refine 1 -snes_monitor_short -snes_type vinewtonssls 256 257 test: 258 suffix: 3 259 requires: !single 260 nsize: 2 261 args: -snes_grid_sequence 2 -snes_vi_monitor -snes_type vinewtonrsls 262 263 test: 264 suffix: mg 265 requires: !single 266 nsize: 4 267 args: -snes_grid_sequence 3 -snes_converged_reason -pc_type mg 268 269 test: 270 suffix: 4 271 nsize: 1 272 args: -mat_is_symmetric 273 274 test: 275 suffix: 5 276 nsize: 1 277 args: -ksp_converged_reason -snes_fd_color 278 279 test: 280 suffix: 6 281 requires: !single 282 nsize: 2 283 args: -snes_grid_sequence 2 -pc_type mg -snes_monitor_short -ksp_converged_reason 284 285 test: 286 suffix: 7 287 nsize: 2 288 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 289 TODO: fix nasty memory leak in SNESCOMPOSITE 290 291 test: 292 suffix: 8 293 nsize: 2 294 args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additive -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor 295 TODO: fix nasty memory leak in SNESCOMPOSITE 296 297 test: 298 suffix: 9 299 nsize: 2 300 args: -da_refine 1 -snes_monitor_short -snes_type composite -snes_composite_type additiveoptimal -snes_composite_sneses vinewtonrsls -sub_0_snes_vi_monitor 301 TODO: fix nasty memory leak in SNESCOMPOSITE 302 303 TEST*/ 304