1c4762a1bSJed Brown 2c4762a1bSJed Brown #include <petscsnes.h> 3c4762a1bSJed Brown #include <petscdm.h> 4c4762a1bSJed Brown #include <petscdmda.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown static const char help[] = "Parallel version of the minimum surface area problem in 2D using DMDA.\n\ 7c4762a1bSJed Brown It solves a system of nonlinear equations in mixed\n\ 8c4762a1bSJed Brown complementarity form.This example is based on a\n\ 9c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and\n\ 10c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\ 11c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\ 12c4762a1bSJed Brown This application solves this problem using complimentarity -- We are actually\n\ 13c4762a1bSJed Brown solving the system (grad f)_i >= 0, if x_i == l_i \n\ 14c4762a1bSJed Brown (grad f)_i = 0, if l_i < x_i < u_i \n\ 15c4762a1bSJed Brown (grad f)_i <= 0, if x_i == u_i \n\ 16c4762a1bSJed Brown where f is the function to be minimized. \n\ 17c4762a1bSJed Brown \n\ 18c4762a1bSJed Brown The command line options are:\n\ 19c4762a1bSJed Brown -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\ 20c4762a1bSJed Brown -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\ 21c4762a1bSJed Brown -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise\n\ 22c4762a1bSJed Brown -lb <value>, lower bound on the variables\n\ 23c4762a1bSJed Brown -ub <value>, upper bound on the variables\n\n"; 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* 26c4762a1bSJed Brown User-defined application context - contains data needed by the 27c4762a1bSJed Brown application-provided call-back routines, FormJacobian() and 28c4762a1bSJed Brown FormFunction(). 29c4762a1bSJed Brown */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* 32c4762a1bSJed Brown This is a new version of the ../tests/ex8.c code 33c4762a1bSJed Brown 34c4762a1bSJed Brown Run, for example, with the options ./ex58 -snes_vi_monitor -ksp_monitor -mg_levels_ksp_monitor -pc_type mg -pc_mg_levels 2 -pc_mg_galerkin pmat -ksp_type fgmres 35c4762a1bSJed Brown 36c4762a1bSJed Brown Or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of 37c4762a1bSJed Brown multigrid levels, it will be determined automatically based on the number of refinements done) 38c4762a1bSJed Brown 39c4762a1bSJed Brown ./ex58 -pc_type mg -ksp_monitor -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 40c4762a1bSJed Brown -mg_levels_ksp_monitor -snes_vi_monitor -mg_levels_pc_type sor -pc_mg_type full 41c4762a1bSJed Brown 42c4762a1bSJed Brown */ 43c4762a1bSJed Brown 44c4762a1bSJed Brown typedef struct { 45c4762a1bSJed Brown PetscScalar *bottom, *top, *left, *right; 46c4762a1bSJed Brown PetscScalar lb,ub; 47c4762a1bSJed Brown } AppCtx; 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 50c4762a1bSJed Brown 51c4762a1bSJed Brown extern PetscErrorCode FormBoundaryConditions(SNES,AppCtx**); 52c4762a1bSJed Brown extern PetscErrorCode DestroyBoundaryConditions(AppCtx**); 53c4762a1bSJed Brown extern PetscErrorCode ComputeInitialGuess(SNES,Vec,void*); 54c4762a1bSJed Brown extern PetscErrorCode FormGradient(SNES,Vec,Vec,void*); 55c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 56c4762a1bSJed Brown extern PetscErrorCode FormBounds(SNES,Vec,Vec); 57c4762a1bSJed Brown 58c4762a1bSJed Brown int main(int argc, char **argv) 59c4762a1bSJed Brown { 60c4762a1bSJed Brown Vec x,r; /* solution and residual vectors */ 61c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 62c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 63c4762a1bSJed Brown DM da; 64c4762a1bSJed Brown 65*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, (char*)0, help)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Create distributed array to manage the 2d grid */ 685f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Extract global vectors from DMDA; */ 735f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &r)); 75c4762a1bSJed Brown 765f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATAIJ)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&J)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Create nonlinear solver context */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes,da)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Set function evaluation and Jacobian evaluation routines */ 845f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormGradient,NULL)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,NULL)); 86c4762a1bSJed Brown 875f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetComputeApplicationContext(snes,(PetscErrorCode (*)(SNES,void**))FormBoundaryConditions,(PetscErrorCode (*)(void**))DestroyBoundaryConditions)); 88c4762a1bSJed Brown 895f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetComputeInitialGuess(snes,ComputeInitialGuess,NULL)); 90c4762a1bSJed Brown 915f80ce2aSJacob Faibussowitsch CHKERRQ(SNESVISetComputeVariableBounds(snes,FormBounds)); 92c4762a1bSJed Brown 935f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* Solve the application */ 965f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* Free memory */ 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Free user-created data structures */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 106c4762a1bSJed Brown 107*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 108*b122ec5aSJacob Faibussowitsch return 0; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* FormBounds - sets the upper and lower bounds 114c4762a1bSJed Brown 115c4762a1bSJed Brown Input Parameters: 116c4762a1bSJed Brown . snes - the SNES context 117c4762a1bSJed Brown 118c4762a1bSJed Brown Output Parameters: 119c4762a1bSJed Brown . xl - lower bounds 120c4762a1bSJed Brown . xu - upper bounds 121c4762a1bSJed Brown */ 122c4762a1bSJed Brown PetscErrorCode FormBounds(SNES snes, Vec xl, Vec xu) 123c4762a1bSJed Brown { 124c4762a1bSJed Brown AppCtx *ctx; 125c4762a1bSJed Brown 126c4762a1bSJed Brown PetscFunctionBeginUser; 1275f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetApplicationContext(snes,&ctx)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xl,ctx->lb)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xu,ctx->ub)); 130c4762a1bSJed Brown PetscFunctionReturn(0); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* FormGradient - Evaluates gradient of f. 136c4762a1bSJed Brown 137c4762a1bSJed Brown Input Parameters: 138c4762a1bSJed Brown . snes - the SNES context 139c4762a1bSJed Brown . X - input vector 140c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 141c4762a1bSJed Brown 142c4762a1bSJed Brown Output Parameters: 143c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 144c4762a1bSJed Brown */ 145c4762a1bSJed Brown PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr) 146c4762a1bSJed Brown { 147c4762a1bSJed Brown AppCtx *user; 148c4762a1bSJed Brown PetscInt i,j; 149c4762a1bSJed Brown PetscInt mx, my; 150c4762a1bSJed Brown PetscScalar hx,hy, hydhx, hxdhy; 151c4762a1bSJed Brown PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 152c4762a1bSJed Brown PetscScalar df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 153c4762a1bSJed Brown PetscScalar **g, **x; 154c4762a1bSJed Brown PetscInt xs,xm,ys,ym; 155c4762a1bSJed Brown Vec localX; 156c4762a1bSJed Brown DM da; 157c4762a1bSJed Brown 158c4762a1bSJed Brown PetscFunctionBeginUser; 1595f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&da)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetApplicationContext(snes,&user)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 162c4762a1bSJed Brown hx = 1.0/(mx+1);hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy; 163c4762a1bSJed Brown 1645f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(G,0.0)); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Get local vector */ 1675f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localX)); 168c4762a1bSJed Brown /* Get ghost points */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 171c4762a1bSJed Brown /* Get pointer to local vector data */ 1725f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,localX, &x)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,G, &g)); 174c4762a1bSJed Brown 1755f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 176c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */ 177c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 178c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 179c4762a1bSJed Brown 180c4762a1bSJed Brown xc = x[j][i]; 181c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 182c4762a1bSJed Brown 183c4762a1bSJed Brown if (i==0) { /* left side */ 184c4762a1bSJed Brown xl = user->left[j+1]; 185c4762a1bSJed Brown xlt = user->left[j+2]; 186c4762a1bSJed Brown } else xl = x[j][i-1]; 187c4762a1bSJed Brown 188c4762a1bSJed Brown if (j==0) { /* bottom side */ 189c4762a1bSJed Brown xb = user->bottom[i+1]; 190c4762a1bSJed Brown xrb = user->bottom[i+2]; 191c4762a1bSJed Brown } else xb = x[j-1][i]; 192c4762a1bSJed Brown 193c4762a1bSJed Brown if (i+1 == mx) { /* right side */ 194c4762a1bSJed Brown xr = user->right[j+1]; 195c4762a1bSJed Brown xrb = user->right[j]; 196c4762a1bSJed Brown } else xr = x[j][i+1]; 197c4762a1bSJed Brown 198c4762a1bSJed Brown if (j+1==0+my) { /* top side */ 199c4762a1bSJed Brown xt = user->top[i+1]; 200c4762a1bSJed Brown xlt = user->top[i]; 201c4762a1bSJed Brown } else xt = x[j+1][i]; 202c4762a1bSJed Brown 203c4762a1bSJed Brown if (i>0 && j+1<my) xlt = x[j+1][i-1]; /* left top side */ 204c4762a1bSJed Brown if (j>0 && i+1<mx) xrb = x[j-1][i+1]; /* right bottom */ 205c4762a1bSJed Brown 206c4762a1bSJed Brown d1 = (xc-xl); 207c4762a1bSJed Brown d2 = (xc-xr); 208c4762a1bSJed Brown d3 = (xc-xt); 209c4762a1bSJed Brown d4 = (xc-xb); 210c4762a1bSJed Brown d5 = (xr-xrb); 211c4762a1bSJed Brown d6 = (xrb-xb); 212c4762a1bSJed Brown d7 = (xlt-xl); 213c4762a1bSJed Brown d8 = (xt-xlt); 214c4762a1bSJed Brown 215c4762a1bSJed Brown df1dxc = d1*hydhx; 216c4762a1bSJed Brown df2dxc = (d1*hydhx + d4*hxdhy); 217c4762a1bSJed Brown df3dxc = d3*hxdhy; 218c4762a1bSJed Brown df4dxc = (d2*hydhx + d3*hxdhy); 219c4762a1bSJed Brown df5dxc = d2*hydhx; 220c4762a1bSJed Brown df6dxc = d4*hxdhy; 221c4762a1bSJed Brown 222c4762a1bSJed Brown d1 /= hx; 223c4762a1bSJed Brown d2 /= hx; 224c4762a1bSJed Brown d3 /= hy; 225c4762a1bSJed Brown d4 /= hy; 226c4762a1bSJed Brown d5 /= hy; 227c4762a1bSJed Brown d6 /= hx; 228c4762a1bSJed Brown d7 /= hy; 229c4762a1bSJed Brown d8 /= hx; 230c4762a1bSJed Brown 231c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 232c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 233c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 234c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 235c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 236c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 237c4762a1bSJed Brown 238c4762a1bSJed Brown df1dxc /= f1; 239c4762a1bSJed Brown df2dxc /= f2; 240c4762a1bSJed Brown df3dxc /= f3; 241c4762a1bSJed Brown df4dxc /= f4; 242c4762a1bSJed Brown df5dxc /= f5; 243c4762a1bSJed Brown df6dxc /= f6; 244c4762a1bSJed Brown 245c4762a1bSJed Brown g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0; 246c4762a1bSJed Brown 247c4762a1bSJed Brown } 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown /* Restore vectors */ 2515f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,localX, &x)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,G, &g)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localX)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(67.0*mx*my)); 255c4762a1bSJed Brown PetscFunctionReturn(0); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 259c4762a1bSJed Brown /* 260c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 261c4762a1bSJed Brown 262c4762a1bSJed Brown Input Parameters: 263c4762a1bSJed Brown . snes - SNES context 264c4762a1bSJed Brown . X - input vector 265c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian() 266c4762a1bSJed Brown 267c4762a1bSJed Brown Output Parameters: 268c4762a1bSJed Brown . tH - Jacobian matrix 269c4762a1bSJed Brown 270c4762a1bSJed Brown */ 271c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *ptr) 272c4762a1bSJed Brown { 273c4762a1bSJed Brown AppCtx *user; 274c4762a1bSJed Brown PetscInt i,j,k; 275c4762a1bSJed Brown PetscInt mx, my; 276c4762a1bSJed Brown MatStencil row,col[7]; 277c4762a1bSJed Brown PetscScalar hx, hy, hydhx, hxdhy; 278c4762a1bSJed Brown PetscScalar f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 279c4762a1bSJed Brown PetscScalar hl,hr,ht,hb,hc,htl,hbr; 280c4762a1bSJed Brown PetscScalar **x, v[7]; 281c4762a1bSJed Brown PetscBool assembled; 282c4762a1bSJed Brown PetscInt xs,xm,ys,ym; 283c4762a1bSJed Brown Vec localX; 284c4762a1bSJed Brown DM da; 285c4762a1bSJed Brown 286c4762a1bSJed Brown PetscFunctionBeginUser; 2875f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&da)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetApplicationContext(snes,&user)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 290c4762a1bSJed Brown hx = 1.0/(mx+1); hy=1.0/(my+1); hydhx=hy/hx; hxdhy=hx/hy; 291c4762a1bSJed Brown 292c4762a1bSJed Brown /* Set various matrix options */ 2935f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssembled(H,&assembled)); 2945f80ce2aSJacob Faibussowitsch if (assembled) CHKERRQ(MatZeroEntries(H)); 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* Get local vector */ 2975f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localX)); 298c4762a1bSJed Brown /* Get ghost points */ 2995f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 3005f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* Get pointers to vector data */ 3035f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,localX, &x)); 304c4762a1bSJed Brown 3055f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 306c4762a1bSJed Brown /* Compute Jacobian over the locally owned part of the mesh */ 307c4762a1bSJed Brown for (j=ys; j< ys+ym; j++) { 308c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 309c4762a1bSJed Brown xc = x[j][i]; 310c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 311c4762a1bSJed Brown 312c4762a1bSJed Brown /* Left */ 313c4762a1bSJed Brown if (i==0) { 314c4762a1bSJed Brown xl = user->left[j+1]; 315c4762a1bSJed Brown xlt = user->left[j+2]; 316c4762a1bSJed Brown } else xl = x[j][i-1]; 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* Bottom */ 319c4762a1bSJed Brown if (j==0) { 320c4762a1bSJed Brown xb =user->bottom[i+1]; 321c4762a1bSJed Brown xrb = user->bottom[i+2]; 322c4762a1bSJed Brown } else xb = x[j-1][i]; 323c4762a1bSJed Brown 324c4762a1bSJed Brown /* Right */ 325c4762a1bSJed Brown if (i+1 == mx) { 326c4762a1bSJed Brown xr =user->right[j+1]; 327c4762a1bSJed Brown xrb = user->right[j]; 328c4762a1bSJed Brown } else xr = x[j][i+1]; 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* Top */ 331c4762a1bSJed Brown if (j+1==my) { 332c4762a1bSJed Brown xt =user->top[i+1]; 333c4762a1bSJed Brown xlt = user->top[i]; 334c4762a1bSJed Brown } else xt = x[j+1][i]; 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* Top left */ 337c4762a1bSJed Brown if (i>0 && j+1<my) xlt = x[j+1][i-1]; 338c4762a1bSJed Brown 339c4762a1bSJed Brown /* Bottom right */ 340c4762a1bSJed Brown if (j>0 && i+1<mx) xrb = x[j-1][i+1]; 341c4762a1bSJed Brown 342c4762a1bSJed Brown d1 = (xc-xl)/hx; 343c4762a1bSJed Brown d2 = (xc-xr)/hx; 344c4762a1bSJed Brown d3 = (xc-xt)/hy; 345c4762a1bSJed Brown d4 = (xc-xb)/hy; 346c4762a1bSJed Brown d5 = (xrb-xr)/hy; 347c4762a1bSJed Brown d6 = (xrb-xb)/hx; 348c4762a1bSJed Brown d7 = (xlt-xl)/hy; 349c4762a1bSJed Brown d8 = (xlt-xt)/hx; 350c4762a1bSJed Brown 351c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 352c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 353c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 354c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 355c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 356c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 357c4762a1bSJed Brown 358c4762a1bSJed Brown hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+ 359c4762a1bSJed Brown (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 360c4762a1bSJed Brown hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+ 361c4762a1bSJed Brown (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 362c4762a1bSJed Brown ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+ 363c4762a1bSJed Brown (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 364c4762a1bSJed Brown hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+ 365c4762a1bSJed Brown (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 366c4762a1bSJed Brown 367c4762a1bSJed Brown hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 368c4762a1bSJed Brown htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 369c4762a1bSJed Brown 370c4762a1bSJed Brown hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + 371c4762a1bSJed Brown hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) + 372c4762a1bSJed Brown (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2.0*d1*d4)/(f2*f2*f2) + 373c4762a1bSJed Brown (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2.0*d2*d3)/(f4*f4*f4); 374c4762a1bSJed Brown 375c4762a1bSJed Brown hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0; 376c4762a1bSJed Brown 377c4762a1bSJed Brown k =0; 378c4762a1bSJed Brown row.i = i;row.j= j; 379c4762a1bSJed Brown /* Bottom */ 380c4762a1bSJed Brown if (j>0) { 381c4762a1bSJed Brown v[k] =hb; 382c4762a1bSJed Brown col[k].i = i; col[k].j=j-1; k++; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown /* Bottom right */ 386c4762a1bSJed Brown if (j>0 && i < mx -1) { 387c4762a1bSJed Brown v[k] =hbr; 388c4762a1bSJed Brown col[k].i = i+1; col[k].j = j-1; k++; 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 391c4762a1bSJed Brown /* left */ 392c4762a1bSJed Brown if (i>0) { 393c4762a1bSJed Brown v[k] = hl; 394c4762a1bSJed Brown col[k].i = i-1; col[k].j = j; k++; 395c4762a1bSJed Brown } 396c4762a1bSJed Brown 397c4762a1bSJed Brown /* Centre */ 398c4762a1bSJed Brown v[k]= hc; col[k].i= row.i; col[k].j = row.j; k++; 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* Right */ 401c4762a1bSJed Brown if (i < mx-1) { 402c4762a1bSJed Brown v[k] = hr; 403c4762a1bSJed Brown col[k].i= i+1; col[k].j = j;k++; 404c4762a1bSJed Brown } 405c4762a1bSJed Brown 406c4762a1bSJed Brown /* Top left */ 407c4762a1bSJed Brown if (i>0 && j < my-1) { 408c4762a1bSJed Brown v[k] = htl; 409c4762a1bSJed Brown col[k].i = i-1;col[k].j = j+1; k++; 410c4762a1bSJed Brown } 411c4762a1bSJed Brown 412c4762a1bSJed Brown /* Top */ 413c4762a1bSJed Brown if (j < my-1) { 414c4762a1bSJed Brown v[k] = ht; 415c4762a1bSJed Brown col[k].i = i; col[k].j = j+1; k++; 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 4185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesStencil(H,1,&row,k,col,v,INSERT_VALUES)); 419c4762a1bSJed Brown } 420c4762a1bSJed Brown } 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* Assemble the matrix */ 4235f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 4245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,localX,&x)); 4255f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 4265f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localX)); 427c4762a1bSJed Brown 4285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(199.0*mx*my)); 429c4762a1bSJed Brown PetscFunctionReturn(0); 430c4762a1bSJed Brown } 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 433c4762a1bSJed Brown /* 434c4762a1bSJed Brown FormBoundaryConditions - Calculates the boundary conditions for 435c4762a1bSJed Brown the region. 436c4762a1bSJed Brown 437c4762a1bSJed Brown Input Parameter: 438c4762a1bSJed Brown . user - user-defined application context 439c4762a1bSJed Brown 440c4762a1bSJed Brown Output Parameter: 441c4762a1bSJed Brown . user - user-defined application context 442c4762a1bSJed Brown */ 443c4762a1bSJed Brown PetscErrorCode FormBoundaryConditions(SNES snes,AppCtx **ouser) 444c4762a1bSJed Brown { 445c4762a1bSJed Brown PetscInt i,j,k,limit=0,maxits=5; 446c4762a1bSJed Brown PetscInt mx,my; 447c4762a1bSJed Brown PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 448c4762a1bSJed Brown PetscScalar one =1.0, two=2.0, three=3.0; 449c4762a1bSJed Brown PetscScalar det,hx,hy,xt=0,yt=0; 450c4762a1bSJed Brown PetscReal fnorm, tol=1e-10; 451c4762a1bSJed Brown PetscScalar u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 452c4762a1bSJed Brown PetscScalar b=-0.5, t=0.5, l=-0.5, r=0.5; 453c4762a1bSJed Brown PetscScalar *boundary; 454c4762a1bSJed Brown AppCtx *user; 455c4762a1bSJed Brown DM da; 456c4762a1bSJed Brown 457c4762a1bSJed Brown PetscFunctionBeginUser; 4585f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&da)); 4595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 460c4762a1bSJed Brown *ouser = user; 461c4762a1bSJed Brown user->lb = .05; 462c4762a1bSJed Brown user->ub = PETSC_INFINITY; 4635f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 464c4762a1bSJed Brown 465c4762a1bSJed Brown /* Check if lower and upper bounds are set */ 4665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL, "-lb", &user->lb, 0)); 4675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL, "-ub", &user->ub, 0)); 468c4762a1bSJed Brown bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2; 469c4762a1bSJed Brown 4705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bsize, &user->bottom)); 4715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(tsize, &user->top)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(lsize, &user->left)); 4735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(rsize, &user->right)); 474c4762a1bSJed Brown 475c4762a1bSJed Brown hx= (r-l)/(mx+1.0); hy=(t-b)/(my+1.0); 476c4762a1bSJed Brown 477c4762a1bSJed Brown for (j=0; j<4; j++) { 478c4762a1bSJed Brown if (j==0) { 479c4762a1bSJed Brown yt = b; 480c4762a1bSJed Brown xt = l; 481c4762a1bSJed Brown limit = bsize; 482c4762a1bSJed Brown boundary = user->bottom; 483c4762a1bSJed Brown } else if (j==1) { 484c4762a1bSJed Brown yt = t; 485c4762a1bSJed Brown xt = l; 486c4762a1bSJed Brown limit = tsize; 487c4762a1bSJed Brown boundary = user->top; 488c4762a1bSJed Brown } else if (j==2) { 489c4762a1bSJed Brown yt = b; 490c4762a1bSJed Brown xt = l; 491c4762a1bSJed Brown limit = lsize; 492c4762a1bSJed Brown boundary = user->left; 493c4762a1bSJed Brown } else { /* if (j==3) */ 494c4762a1bSJed Brown yt = b; 495c4762a1bSJed Brown xt = r; 496c4762a1bSJed Brown limit = rsize; 497c4762a1bSJed Brown boundary = user->right; 498c4762a1bSJed Brown } 499c4762a1bSJed Brown 500c4762a1bSJed Brown for (i=0; i<limit; i++) { 501c4762a1bSJed Brown u1=xt; 502c4762a1bSJed Brown u2=-yt; 503c4762a1bSJed Brown for (k=0; k<maxits; k++) { 504c4762a1bSJed Brown nf1 = u1 + u1*u2*u2 - u1*u1*u1/three-xt; 505c4762a1bSJed Brown nf2 = -u2 - u1*u1*u2 + u2*u2*u2/three-yt; 506c4762a1bSJed Brown fnorm = PetscRealPart(PetscSqrtScalar(nf1*nf1+nf2*nf2)); 507c4762a1bSJed Brown if (fnorm <= tol) break; 508c4762a1bSJed Brown njac11=one+u2*u2-u1*u1; 509c4762a1bSJed Brown njac12=two*u1*u2; 510c4762a1bSJed Brown njac21=-two*u1*u2; 511c4762a1bSJed Brown njac22=-one - u1*u1 + u2*u2; 512c4762a1bSJed Brown det = njac11*njac22-njac21*njac12; 513c4762a1bSJed Brown u1 = u1-(njac22*nf1-njac12*nf2)/det; 514c4762a1bSJed Brown u2 = u2-(njac11*nf2-njac21*nf1)/det; 515c4762a1bSJed Brown } 516c4762a1bSJed Brown 517c4762a1bSJed Brown boundary[i]=u1*u1-u2*u2; 518c4762a1bSJed Brown if (j==0 || j==1) xt=xt+hx; 519c4762a1bSJed Brown else yt=yt+hy; /* if (j==2 || j==3) */ 520c4762a1bSJed Brown } 521c4762a1bSJed Brown } 522c4762a1bSJed Brown PetscFunctionReturn(0); 523c4762a1bSJed Brown } 524c4762a1bSJed Brown 525c4762a1bSJed Brown PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser) 526c4762a1bSJed Brown { 527c4762a1bSJed Brown AppCtx *user = *ouser; 528c4762a1bSJed Brown 529c4762a1bSJed Brown PetscFunctionBeginUser; 5305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user->bottom)); 5315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user->top)); 5325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user->left)); 5335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user->right)); 5345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(*ouser)); 535c4762a1bSJed Brown PetscFunctionReturn(0); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown 538c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 539c4762a1bSJed Brown /* 540c4762a1bSJed Brown ComputeInitialGuess - Calculates the initial guess 541c4762a1bSJed Brown 542c4762a1bSJed Brown Input Parameters: 543c4762a1bSJed Brown . user - user-defined application context 544c4762a1bSJed Brown . X - vector for initial guess 545c4762a1bSJed Brown 546c4762a1bSJed Brown Output Parameters: 547c4762a1bSJed Brown . X - newly computed initial guess 548c4762a1bSJed Brown */ 549c4762a1bSJed Brown PetscErrorCode ComputeInitialGuess(SNES snes, Vec X,void *dummy) 550c4762a1bSJed Brown { 551c4762a1bSJed Brown PetscInt i,j,mx,my; 552c4762a1bSJed Brown DM da; 553c4762a1bSJed Brown AppCtx *user; 554c4762a1bSJed Brown PetscScalar **x; 555c4762a1bSJed Brown PetscInt xs,xm,ys,ym; 556c4762a1bSJed Brown 557c4762a1bSJed Brown PetscFunctionBeginUser; 5585f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&da)); 5595f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetApplicationContext(snes,&user)); 560c4762a1bSJed Brown 5615f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 5625f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 563c4762a1bSJed Brown 564c4762a1bSJed Brown /* Get pointers to vector data */ 5655f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,X,&x)); 566c4762a1bSJed Brown /* Perform local computations */ 567c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 568c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 569c4762a1bSJed Brown x[j][i] = (((j+1.0)*user->bottom[i+1]+(my-j+1.0)*user->top[i+1])/(my+2.0)+((i+1.0)*user->left[j+1]+(mx-i+1.0)*user->right[j+1])/(mx+2.0))/2.0; 570c4762a1bSJed Brown } 571c4762a1bSJed Brown } 572c4762a1bSJed Brown /* Restore vectors */ 5735f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,X,&x)); 574c4762a1bSJed Brown PetscFunctionReturn(0); 575c4762a1bSJed Brown } 576c4762a1bSJed Brown 577c4762a1bSJed Brown /*TEST 578c4762a1bSJed Brown 579c4762a1bSJed Brown test: 580c4762a1bSJed Brown args: -snes_type vinewtonrsls -pc_type mg -ksp_monitor_short -pc_mg_galerkin pmat -da_refine 5 -snes_vi_monitor -pc_mg_type full -snes_max_it 100 -snes_converged_reason 581c4762a1bSJed Brown requires: !single 582c4762a1bSJed Brown 583c4762a1bSJed Brown test: 584c4762a1bSJed Brown suffix: 2 585c4762a1bSJed Brown args: -snes_type vinewtonssls -pc_type mg -ksp_monitor_short -pc_mg_galerkin pmat -da_refine 5 -snes_vi_monitor -pc_mg_type full -snes_max_it 100 -snes_converged_reason 586c4762a1bSJed Brown requires: !single 587c4762a1bSJed Brown 588c4762a1bSJed Brown TEST*/ 589