Actual source code: ex29.c
2: static char help[] ="Tests ML interface. Modified from ~src/ksp/ksp/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*);
44: int main(int argc,char **argv)
45: {
46: PetscInt its,n,Nx=PETSC_DECIDE,Ny=PETSC_DECIDE,nlocal,i;
47: PetscMPIInt size;
48: PC pc;
49: PetscInt mx,my;
50: Mat A;
51: GridCtx fine_ctx;
52: KSP ksp;
53: PetscBool flg;
55: PetscInitialize(&argc,&argv,NULL,help);
56: /* set up discretization matrix for fine grid */
57: /* ML requires input of fine-grid matrix. It determines nlevels. */
58: fine_ctx.mx = 9; fine_ctx.my = 9;
59: PetscOptionsGetInt(NULL,NULL,"-mx",&mx,&flg);
60: if (flg) fine_ctx.mx = mx;
61: PetscOptionsGetInt(NULL,NULL,"-my",&my,&flg);
62: if (flg) fine_ctx.my = my;
63: PetscPrintf(PETSC_COMM_WORLD,"Fine grid size %D by %D\n",fine_ctx.mx,fine_ctx.my);
64: n = fine_ctx.mx*fine_ctx.my;
66: MPI_Comm_size(PETSC_COMM_WORLD,&size);
67: PetscOptionsGetInt(NULL,NULL,"-Nx",&Nx,NULL);
68: PetscOptionsGetInt(NULL,NULL,"-Ny",&Ny,NULL);
70: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,fine_ctx.mx,fine_ctx.my,Nx,Ny,1,1,NULL,NULL,&fine_ctx.da);
71: DMSetFromOptions(fine_ctx.da);
72: DMSetUp(fine_ctx.da);
73: DMCreateGlobalVector(fine_ctx.da,&fine_ctx.x);
74: VecDuplicate(fine_ctx.x,&fine_ctx.b);
75: VecGetLocalSize(fine_ctx.x,&nlocal);
76: DMCreateLocalVector(fine_ctx.da,&fine_ctx.localX);
77: VecDuplicate(fine_ctx.localX,&fine_ctx.localF);
78: MatCreateAIJ(PETSC_COMM_WORLD,nlocal,nlocal,n,n,5,NULL,3,NULL,&A);
79: FormJacobian_Grid(&fine_ctx,&A);
81: /* create linear solver */
82: KSPCreate(PETSC_COMM_WORLD,&ksp);
83: KSPGetPC(ksp,&pc);
84: PCSetType(pc,PCML);
86: /* set options, then solve system */
87: KSPSetFromOptions(ksp); /* calls PCSetFromOptions_MG/ML */
89: for (i=0; i<3; i++) {
90: if (i<2) {
91: /* set values for rhs vector */
92: VecSet(fine_ctx.b,i+1.0);
93: /* modify A */
94: MatShift(A,1.0);
95: MatScale(A,2.0);
96: KSPSetOperators(ksp,A,A);
97: } else { /* test SAME_NONZERO_PATTERN */
98: KSPSetOperators(ksp,A,A);
99: }
100: KSPSolve(ksp,fine_ctx.b,fine_ctx.x);
101: KSPGetIterationNumber(ksp,&its);
102: if (its > 6) {
103: PetscPrintf(PETSC_COMM_WORLD,"Warning: Number of iterations = %D greater than expected\n",its);
104: }
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: }
119: int FormJacobian_Grid(GridCtx *grid,Mat *J)
120: {
121: Mat jac = *J;
122: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
123: PetscInt grow;
124: const PetscInt *ltog;
125: PetscScalar two = 2.0,one = 1.0,v[5],hx,hy,hxdhy,hydhx,value;
126: ISLocalToGlobalMapping ltogm;
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: DMGetLocalToGlobalMapping(grid->da,<ogm);
136: ISLocalToGlobalMappingGetIndices(ltogm,<og);
138: /* Evaluate Jacobian of function */
139: for (j=ys; j<ys+ym; j++) {
140: row = (j - Ys)*Xm + xs - Xs - 1;
141: for (i=xs; i<xs+xm; i++) {
142: row++;
143: grow = ltog[row];
144: if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
145: v[0] = -hxdhy; col[0] = ltog[row - Xm];
146: v[1] = -hydhx; col[1] = ltog[row - 1];
147: v[2] = two*(hydhx + hxdhy); col[2] = grow;
148: v[3] = -hydhx; col[3] = ltog[row + 1];
149: v[4] = -hxdhy; col[4] = ltog[row + Xm];
150: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
151: } else if ((i > 0 && i < mx-1) || (j > 0 && j < my-1)) {
152: value = .5*two*(hydhx + hxdhy);
153: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
154: } else {
155: value = .25*two*(hydhx + hxdhy);
156: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
157: }
158: }
159: }
160: ISLocalToGlobalMappingRestoreIndices(ltogm,<og);
161: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
162: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
163: return 0;
164: }
166: /*TEST
168: test:
169: output_file: output/ex29.out
170: args: -mat_no_inode
171: requires: ml
173: test:
174: suffix: 2
175: nsize: 3
176: requires: ml
178: TEST*/