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