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 const PetscReal r = x * x + y * y, r0 = 0.9, psi0 = PetscSqrtReal(1.0 - r0 * r0), dpsi0 = -r0 / psi0; 38 if (r <= r0) { 39 return PetscSqrtReal(1.0 - r); 40 } else { 41 return psi0 + dpsi0 * (r - r0); 42 } 43 } 44 45 /* This exact solution solves a 1D radial free-boundary problem for the 46 Laplace equation, on the interval 0 < r < 2, with above obstacle psi(x,y). 47 The Laplace equation applies where u(r) > psi(r), 48 u''(r) + r^-1 u'(r) = 0 49 with boundary conditions including free b.c.s at an unknown location r = a: 50 u(a) = psi(a), u'(a) = psi'(a), u(2) = 0 51 The solution is u(r) = - A log(r) + B on r > a. The boundary conditions 52 can then be reduced to a root-finding problem for a: 53 a^2 (log(2) - log(a)) = 1 - a^2 54 The solution is a = 0.697965148223374 (giving residual 1.5e-15). Then 55 A = a^2*(1-a^2)^(-0.5) and B = A*log(2) are as given below in the code. */ 56 PetscReal u_exact(PetscReal x, PetscReal y) { 57 const PetscReal afree = 0.697965148223374, A = 0.680259411891719, B = 0.471519893402112; 58 PetscReal r; 59 r = PetscSqrtReal(x * x + y * y); 60 return (r <= afree) ? psi(x, y) /* active set; on the obstacle */ 61 : -A * PetscLogReal(r) + B; /* solves laplace eqn */ 62 } 63 64 extern PetscErrorCode FormExactSolution(DMDALocalInfo *, Vec); 65 extern PetscErrorCode FormBounds(SNES, Vec, Vec); 66 extern PetscErrorCode FormFunctionLocal(DMDALocalInfo *, PetscReal **, PetscReal **, void *); 67 extern PetscErrorCode FormJacobianLocal(DMDALocalInfo *, PetscReal **, Mat, Mat, void *); 68 69 int main(int argc, char **argv) { 70 SNES snes; 71 DM da, da_after; 72 Vec u, u_exact; 73 DMDALocalInfo info; 74 PetscReal error1, errorinf; 75 76 PetscFunctionBeginUser; 77 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 78 79 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 */ 80 PETSC_DECIDE, PETSC_DECIDE, 1, 1, /* dof=1 and s = 1 (stencil extends out one cell) */ 81 NULL, NULL, &da)); 82 PetscCall(DMSetFromOptions(da)); 83 PetscCall(DMSetUp(da)); 84 PetscCall(DMDASetUniformCoordinates(da, -2.0, 2.0, -2.0, 2.0, 0.0, 1.0)); 85 86 PetscCall(DMCreateGlobalVector(da, &u)); 87 PetscCall(VecSet(u, 0.0)); 88 89 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 90 PetscCall(SNESSetDM(snes, da)); 91 PetscCall(SNESSetType(snes, SNESVINEWTONRSLS)); 92 PetscCall(SNESVISetComputeVariableBounds(snes, &FormBounds)); 93 PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunction)FormFunctionLocal, NULL)); 94 PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobian)FormJacobianLocal, NULL)); 95 PetscCall(SNESSetFromOptions(snes)); 96 97 /* solve nonlinear system */ 98 PetscCall(SNESSolve(snes, NULL, u)); 99 PetscCall(VecDestroy(&u)); 100 PetscCall(DMDestroy(&da)); 101 /* DMDA after solve may be different, e.g. with -snes_grid_sequence */ 102 PetscCall(SNESGetDM(snes, &da_after)); 103 PetscCall(SNESGetSolution(snes, &u)); /* do not destroy u */ 104 PetscCall(DMDAGetLocalInfo(da_after, &info)); 105 PetscCall(VecDuplicate(u, &u_exact)); 106 PetscCall(FormExactSolution(&info, u_exact)); 107 PetscCall(VecAXPY(u, -1.0, u_exact)); /* u <-- u - u_exact */ 108 PetscCall(VecNorm(u, NORM_1, &error1)); 109 error1 /= (PetscReal)info.mx * (PetscReal)info.my; /* average error */ 110 PetscCall(VecNorm(u, NORM_INFINITY, &errorinf)); 111 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)); 112 PetscCall(VecDestroy(&u_exact)); 113 PetscCall(SNESDestroy(&snes)); 114 PetscCall(DMDestroy(&da)); 115 PetscCall(PetscFinalize()); 116 return 0; 117 } 118 119 PetscErrorCode FormExactSolution(DMDALocalInfo *info, Vec u) { 120 PetscInt i, j; 121 PetscReal **au, dx, dy, x, y; 122 dx = 4.0 / (PetscReal)(info->mx - 1); 123 dy = 4.0 / (PetscReal)(info->my - 1); 124 PetscCall(DMDAVecGetArray(info->da, u, &au)); 125 for (j = info->ys; j < info->ys + info->ym; j++) { 126 y = -2.0 + j * dy; 127 for (i = info->xs; i < info->xs + info->xm; i++) { 128 x = -2.0 + i * dx; 129 au[j][i] = u_exact(x, y); 130 } 131 } 132 PetscCall(DMDAVecRestoreArray(info->da, u, &au)); 133 return 0; 134 } 135 136 PetscErrorCode FormBounds(SNES snes, Vec Xl, Vec Xu) { 137 DM da; 138 DMDALocalInfo info; 139 PetscInt i, j; 140 PetscReal **aXl, dx, dy, x, y; 141 142 PetscCall(SNESGetDM(snes, &da)); 143 PetscCall(DMDAGetLocalInfo(da, &info)); 144 dx = 4.0 / (PetscReal)(info.mx - 1); 145 dy = 4.0 / (PetscReal)(info.my - 1); 146 PetscCall(DMDAVecGetArray(da, Xl, &aXl)); 147 for (j = info.ys; j < info.ys + info.ym; j++) { 148 y = -2.0 + j * dy; 149 for (i = info.xs; i < info.xs + info.xm; i++) { 150 x = -2.0 + i * dx; 151 aXl[j][i] = psi(x, y); 152 } 153 } 154 PetscCall(DMDAVecRestoreArray(da, Xl, &aXl)); 155 PetscCall(VecSet(Xu, PETSC_INFINITY)); 156 return 0; 157 } 158 159 PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **au, PetscScalar **af, void *user) { 160 PetscInt i, j; 161 PetscReal dx, dy, x, y, ue, un, us, uw; 162 163 PetscFunctionBeginUser; 164 dx = 4.0 / (PetscReal)(info->mx - 1); 165 dy = 4.0 / (PetscReal)(info->my - 1); 166 for (j = info->ys; j < info->ys + info->ym; j++) { 167 y = -2.0 + j * dy; 168 for (i = info->xs; i < info->xs + info->xm; i++) { 169 x = -2.0 + i * dx; 170 if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 171 af[j][i] = 4.0 * (au[j][i] - u_exact(x, y)); 172 } else { 173 uw = (i - 1 == 0) ? u_exact(x - dx, y) : au[j][i - 1]; 174 ue = (i + 1 == info->mx - 1) ? u_exact(x + dx, y) : au[j][i + 1]; 175 us = (j - 1 == 0) ? u_exact(x, y - dy) : au[j - 1][i]; 176 un = (j + 1 == info->my - 1) ? u_exact(x, y + dy) : au[j + 1][i]; 177 af[j][i] = -(dy / dx) * (uw - 2.0 * au[j][i] + ue) - (dx / dy) * (us - 2.0 * au[j][i] + un); 178 } 179 } 180 } 181 PetscCall(PetscLogFlops(12.0 * info->ym * info->xm)); 182 PetscFunctionReturn(0); 183 } 184 185 PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **au, Mat A, Mat jac, void *user) { 186 PetscInt i, j, n; 187 MatStencil col[5], row; 188 PetscReal v[5], dx, dy, oxx, oyy; 189 190 PetscFunctionBeginUser; 191 dx = 4.0 / (PetscReal)(info->mx - 1); 192 dy = 4.0 / (PetscReal)(info->my - 1); 193 oxx = dy / dx; 194 oyy = dx / dy; 195 for (j = info->ys; j < info->ys + info->ym; j++) { 196 for (i = info->xs; i < info->xs + info->xm; i++) { 197 row.j = j; 198 row.i = i; 199 if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { /* boundary */ 200 v[0] = 4.0; 201 PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES)); 202 } else { /* interior grid points */ 203 v[0] = 2.0 * (oxx + oyy); 204 col[0].j = j; 205 col[0].i = i; 206 n = 1; 207 if (i - 1 > 0) { 208 v[n] = -oxx; 209 col[n].j = j; 210 col[n++].i = i - 1; 211 } 212 if (i + 1 < info->mx - 1) { 213 v[n] = -oxx; 214 col[n].j = j; 215 col[n++].i = i + 1; 216 } 217 if (j - 1 > 0) { 218 v[n] = -oyy; 219 col[n].j = j - 1; 220 col[n++].i = i; 221 } 222 if (j + 1 < info->my - 1) { 223 v[n] = -oyy; 224 col[n].j = j + 1; 225 col[n++].i = i; 226 } 227 PetscCall(MatSetValuesStencil(jac, 1, &row, n, col, v, INSERT_VALUES)); 228 } 229 } 230 } 231 232 /* Assemble matrix, using the 2-step process: */ 233 PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 234 PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 235 if (A != jac) { 236 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 237 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 238 } 239 PetscCall(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