Actual source code: ex29.c

petsc-3.7.1 2016-05-15
Report Typos and Errors
  2: static char help[] ="Tests ML interface. Modified from ~src/ksp/ksp/examples/tests/ex19.c \n\
  3:   -mx <xg>, where <xg> = number of grid points in the x-direction\n\
  4:   -my <yg>, where <yg> = number of grid points in the y-direction\n\
  5:   -Nx <npx>, where <npx> = number of processors in the x-direction\n\
  6:   -Ny <npy>, where <npy> = number of processors in the y-direction\n\n";

  8: /*
  9:     This problem is modeled by
 10:     the partial differential equation

 12:             -Laplacian u  = g,  0 < x,y < 1,

 14:     with boundary conditions

 16:              u = 0  for  x = 0, x = 1, y = 0, y = 1.

 18:     A finite difference approximation with the usual 5-point stencil
 19:     is used to discretize the boundary value problem to obtain a nonlinear
 20:     system of equations.

 22:     Usage: ./ex29 -ksp_type gmres -ksp_monitor
 23:            -pc_mg_type <multiplicative> (one of) additive multiplicative full kascade
 24:            -mg_coarse_ksp_max_it 10 -mg_levels_3_ksp_max_it 10 -mg_levels_2_ksp_max_it 10
 25:            -mg_levels_1_ksp_max_it 10 -mg_fine_ksp_max_it 10
 26: */

 28: #include <petscksp.h>
 29: #include <petscdm.h>
 30: #include <petscdmda.h>

 32: /* User-defined application contexts */
 33: typedef struct {
 34:   PetscInt mx,my;              /* number grid points in x and y direction */
 35:   Vec      localX,localF;      /* local vectors with ghost region */
 36:   DM       da;
 37:   Vec      x,b,r;              /* global vectors */
 38:   Mat      J;                  /* Jacobian on grid */
 39:   Mat      A,P,R;
 40:   KSP      ksp;
 41: } GridCtx;
 42: extern int FormJacobian_Grid(GridCtx*,Mat*);

 46: int main(int argc,char **argv)
 47: {
 49:   PetscInt       its,n,Nx=PETSC_DECIDE,Ny=PETSC_DECIDE,nlocal,i;
 50:   PetscMPIInt    size;
 51:   PC             pc;
 52:   PetscInt       mx,my;
 53:   Mat            A;
 54:   GridCtx        fine_ctx;
 55:   KSP            ksp;
 56:   PetscBool      flg;

 58:   PetscInitialize(&argc,&argv,NULL,help);
 59:   /* set up discretization matrix for fine grid */
 60:   /* ML requires input of fine-grid matrix. It determines nlevels. */
 61:   fine_ctx.mx = 9; fine_ctx.my = 9;
 62:   PetscOptionsGetInt(NULL,NULL,"-mx",&mx,&flg);
 63:   if (flg) fine_ctx.mx = mx;
 64:   PetscOptionsGetInt(NULL,NULL,"-my",&my,&flg);
 65:   if (flg) fine_ctx.my = my;
 66:   PetscPrintf(PETSC_COMM_WORLD,"Fine grid size %D by %D\n",fine_ctx.mx,fine_ctx.my);
 67:   n    = fine_ctx.mx*fine_ctx.my;

 69:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 70:   PetscOptionsGetInt(NULL,NULL,"-Nx",&Nx,NULL);
 71:   PetscOptionsGetInt(NULL,NULL,"-Ny",&Ny,NULL);

 73:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,fine_ctx.mx,
 74:                       fine_ctx.my,Nx,Ny,1,1,NULL,NULL,&fine_ctx.da);
 75:   DMCreateGlobalVector(fine_ctx.da,&fine_ctx.x);
 76:   VecDuplicate(fine_ctx.x,&fine_ctx.b);
 77:   VecGetLocalSize(fine_ctx.x,&nlocal);
 78:   DMCreateLocalVector(fine_ctx.da,&fine_ctx.localX);
 79:   VecDuplicate(fine_ctx.localX,&fine_ctx.localF);
 80:   MatCreateAIJ(PETSC_COMM_WORLD,nlocal,nlocal,n,n,5,NULL,3,NULL,&A);
 81:   FormJacobian_Grid(&fine_ctx,&A);

 83:   /* create linear solver */
 84:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 85:   KSPGetPC(ksp,&pc);
 86:   PCSetType(pc,PCML);

 88:   /* set options, then solve system */
 89:   KSPSetFromOptions(ksp); /* calls PCSetFromOptions_MG/ML */

 91:   for (i=0; i<3; i++) {
 92:     if (i<2) {
 93:       /* set values for rhs vector */
 94:       VecSet(fine_ctx.b,i+1.0);
 95:       /* modify A */
 96:       MatShift(A,1.0);
 97:       MatScale(A,2.0);
 98:       KSPSetOperators(ksp,A,A);
 99:     } else {  /* test SAME_NONZERO_PATTERN */
100:       KSPSetOperators(ksp,A,A);
101:     }
102:     KSPSolve(ksp,fine_ctx.b,fine_ctx.x);
103:     KSPGetIterationNumber(ksp,&its);
104:     PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %D\n",its);
105:   }

107:   /* free data structures */
108:   VecDestroy(&fine_ctx.x);
109:   VecDestroy(&fine_ctx.b);
110:   DMDestroy(&fine_ctx.da);
111:   VecDestroy(&fine_ctx.localX);
112:   VecDestroy(&fine_ctx.localF);
113:   MatDestroy(&A);
114:   KSPDestroy(&ksp);
115:   PetscFinalize();
116:   return 0;
117: }

121: int FormJacobian_Grid(GridCtx *grid,Mat *J)
122: {
123:   Mat                    jac = *J;
124:   PetscErrorCode         ierr;
125:   PetscInt               i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
126:   PetscInt               grow;
127:   const PetscInt         *ltog;
128:   PetscScalar            two = 2.0,one = 1.0,v[5],hx,hy,hxdhy,hydhx,value;
129:   ISLocalToGlobalMapping ltogm;

131:   mx    = grid->mx;               my = grid->my;
132:   hx    = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
133:   hxdhy = hx/hy;               hydhx = hy/hx;

135:   /* Get ghost points */
136:   DMDAGetCorners(grid->da,&xs,&ys,0,&xm,&ym,0);
137:   DMDAGetGhostCorners(grid->da,&Xs,&Ys,0,&Xm,&Ym,0);
138:   DMGetLocalToGlobalMapping(grid->da,&ltogm);
139:   ISLocalToGlobalMappingGetIndices(ltogm,&ltog);

141:   /* Evaluate Jacobian of function */
142:   for (j=ys; j<ys+ym; j++) {
143:     row = (j - Ys)*Xm + xs - Xs - 1;
144:     for (i=xs; i<xs+xm; i++) {
145:       row++;
146:       grow = ltog[row];
147:       if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
148:         v[0] = -hxdhy; col[0] = ltog[row - Xm];
149:         v[1] = -hydhx; col[1] = ltog[row - 1];
150:         v[2] = two*(hydhx + hxdhy); col[2] = grow;
151:         v[3] = -hydhx; col[3] = ltog[row + 1];
152:         v[4] = -hxdhy; col[4] = ltog[row + Xm];
153:         MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
154:       } else if ((i > 0 && i < mx-1) || (j > 0 && j < my-1)) {
155:         value = .5*two*(hydhx + hxdhy);
156:         MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
157:       } else {
158:         value = .25*two*(hydhx + hxdhy);
159:         MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
160:       }
161:     }
162:   }
163:   ISLocalToGlobalMappingRestoreIndices(ltogm,&ltog);
164:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
165:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
166:   return 0;
167: }