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