Actual source code: ex29.c
petsc-3.12.2 2019-11-22
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*);
44: int main(int argc,char **argv)
45: {
47: PetscInt its,n,Nx=PETSC_DECIDE,Ny=PETSC_DECIDE,nlocal,i;
48: PetscMPIInt size;
49: PC pc;
50: PetscInt mx,my;
51: Mat A;
52: GridCtx fine_ctx;
53: KSP ksp;
54: PetscBool flg;
56: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
57: /* set up discretization matrix for fine grid */
58: /* ML requires input of fine-grid matrix. It determines nlevels. */
59: fine_ctx.mx = 9; fine_ctx.my = 9;
60: PetscOptionsGetInt(NULL,NULL,"-mx",&mx,&flg);
61: if (flg) fine_ctx.mx = mx;
62: PetscOptionsGetInt(NULL,NULL,"-my",&my,&flg);
63: if (flg) fine_ctx.my = my;
64: PetscPrintf(PETSC_COMM_WORLD,"Fine grid size %D by %D\n",fine_ctx.mx,fine_ctx.my);
65: n = fine_ctx.mx*fine_ctx.my;
67: MPI_Comm_size(PETSC_COMM_WORLD,&size);
68: PetscOptionsGetInt(NULL,NULL,"-Nx",&Nx,NULL);
69: PetscOptionsGetInt(NULL,NULL,"-Ny",&Ny,NULL);
71: 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);
72: DMSetFromOptions(fine_ctx.da);
73: DMSetUp(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: MatCreateAIJ(PETSC_COMM_WORLD,nlocal,nlocal,n,n,5,NULL,3,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) {
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);
98: } else { /* test SAME_NONZERO_PATTERN */
99: KSPSetOperators(ksp,A,A);
100: }
101: KSPSolve(ksp,fine_ctx.b,fine_ctx.x);
102: KSPGetIterationNumber(ksp,&its);
103: if (its > 6) {
104: PetscPrintf(PETSC_COMM_WORLD,"Warning: Number of iterations = %D greater than expected\n",its);
105: }
106: }
108: /* free data structures */
109: VecDestroy(&fine_ctx.x);
110: VecDestroy(&fine_ctx.b);
111: DMDestroy(&fine_ctx.da);
112: VecDestroy(&fine_ctx.localX);
113: VecDestroy(&fine_ctx.localF);
114: MatDestroy(&A);
115: KSPDestroy(&ksp);
116: PetscFinalize();
117: return ierr;
118: }
120: int FormJacobian_Grid(GridCtx *grid,Mat *J)
121: {
122: Mat jac = *J;
123: PetscErrorCode ierr;
124: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
125: PetscInt grow;
126: const PetscInt *ltog;
127: PetscScalar two = 2.0,one = 1.0,v[5],hx,hy,hxdhy,hydhx,value;
128: ISLocalToGlobalMapping ltogm;
130: mx = grid->mx; my = grid->my;
131: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
132: hxdhy = hx/hy; hydhx = hy/hx;
134: /* Get ghost points */
135: DMDAGetCorners(grid->da,&xs,&ys,0,&xm,&ym,0);
136: DMDAGetGhostCorners(grid->da,&Xs,&Ys,0,&Xm,&Ym,0);
137: DMGetLocalToGlobalMapping(grid->da,<ogm);
138: ISLocalToGlobalMappingGetIndices(ltogm,<og);
140: /* Evaluate Jacobian of function */
141: for (j=ys; j<ys+ym; j++) {
142: row = (j - Ys)*Xm + xs - Xs - 1;
143: for (i=xs; i<xs+xm; i++) {
144: row++;
145: grow = ltog[row];
146: if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
147: v[0] = -hxdhy; col[0] = ltog[row - Xm];
148: v[1] = -hydhx; col[1] = ltog[row - 1];
149: v[2] = two*(hydhx + hxdhy); col[2] = grow;
150: v[3] = -hydhx; col[3] = ltog[row + 1];
151: v[4] = -hxdhy; col[4] = ltog[row + Xm];
152: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
153: } else if ((i > 0 && i < mx-1) || (j > 0 && j < my-1)) {
154: value = .5*two*(hydhx + hxdhy);
155: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
156: } else {
157: value = .25*two*(hydhx + hxdhy);
158: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
159: }
160: }
161: }
162: ISLocalToGlobalMappingRestoreIndices(ltogm,<og);
163: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
164: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
165: return 0;
166: }
168: /*TEST
170: test:
171: output_file: output/ex29.out
172: args: -mat_no_inode
173: requires: ml
175: test:
176: suffix: 2
177: nsize: 3
178: requires: ml
180: TEST*/