Actual source code: ex29.c

  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
 11:   
 12:             -Laplacian u  = g,  0 < x,y < 1,
 13:   
 14:     with boundary conditions
 15:    
 16:              u = 0  for  x = 0, x = 1, y = 0, y = 1.
 17:   
 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 cascade 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 <petscdmda.h>

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

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

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

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

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

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

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

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

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

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

128:   mx = grid->mx;            my = grid->my;
129:   hx = one/(PetscReal)(mx-1);  hy = one/(PetscReal)(my-1);
130:   hxdhy = hx/hy;            hydhx = hy/hx;

132:   /* Get ghost points */
133:   DMDAGetCorners(grid->da,&xs,&ys,0,&xm,&ym,0);
134:   DMDAGetGhostCorners(grid->da,&Xs,&Ys,0,&Xm,&Ym,0);
135:   DMDAGetGlobalIndices(grid->da,&nloc,&ltog);

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