Actual source code: ex26.c
1: static char help[] ="Solvers Laplacian with multigrid, bad way.\n\
2: -mx <xg>, where <xg> = number of grid points in the x-direction\n\
3: -my <yg>, where <yg> = number of grid points in the y-direction\n\
4: -Nx <npx>, where <npx> = number of processors in the x-direction\n\
5: -Ny <npy>, where <npy> = number of processors in the y-direction\n\n";
7: /* Modified from ~src/ksp/examples/tests/ex19.c. Used for testing ML 6.2 interface.
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: ./ex26 -ksp_monitor_short -pc_type ml
23: -mg_coarse_ksp_max_it 10
24: -mg_levels_1_ksp_max_it 10 -mg_levels_2_ksp_max_it 10
25: -mg_fine_ksp_max_it 10
26: */
28: #include petscksp.h
29: #include petscda.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: DA 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;
49: PetscMPIInt size;
50: PetscScalar one = 1.0;
51: PetscInt mx,my;
52: Mat A;
53: GridCtx fine_ctx;
54: KSP ksp;
55: PetscTruth flg;
57: PetscInitialize(&argc,&argv,(char *)0,help);
58: /* set up discretization matrix for fine grid */
59: fine_ctx.mx = 9; fine_ctx.my = 9;
60: PetscOptionsGetInt(PETSC_NULL,"-mx",&mx,&flg);
61: if (flg) fine_ctx.mx = mx;
62: PetscOptionsGetInt(PETSC_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(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
69: PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
71: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,fine_ctx.mx,
72: fine_ctx.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&fine_ctx.da);
73: DACreateGlobalVector(fine_ctx.da,&fine_ctx.x);
74: VecDuplicate(fine_ctx.x,&fine_ctx.b);
75: VecGetLocalSize(fine_ctx.x,&nlocal);
76: DACreateLocalVector(fine_ctx.da,&fine_ctx.localX);
77: VecDuplicate(fine_ctx.localX,&fine_ctx.localF);
78: MatCreateMPIAIJ(PETSC_COMM_WORLD,nlocal,nlocal,n,n,5,PETSC_NULL,3,PETSC_NULL,&A);
79: FormJacobian_Grid(&fine_ctx,&A);
81: /* create linear solver */
82: KSPCreate(PETSC_COMM_WORLD,&ksp);
84: /* set values for rhs vector */
85: VecSet(fine_ctx.b,one);
86: {
87: PetscRandom rdm;
88: PetscRandomCreate(PETSC_COMM_WORLD,&rdm);
89: PetscRandomSetFromOptions(rdm);
90: VecSetRandom(fine_ctx.b,rdm);
91: PetscRandomDestroy(rdm);
92: }
94: /* set options, then solve system */
95: KSPSetFromOptions(ksp); /* calls PCSetFromOptions_ML if 'pc_type=ml' */
96: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
97: KSPSolve(ksp,fine_ctx.b,fine_ctx.x);
98: KSPGetIterationNumber(ksp,&its);
99: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %D\n",its);
101: /* free data structures */
102: VecDestroy(fine_ctx.x);
103: VecDestroy(fine_ctx.b);
104: DADestroy(fine_ctx.da);
105: VecDestroy(fine_ctx.localX);
106: VecDestroy(fine_ctx.localF);
107: MatDestroy(A);
108: KSPDestroy(ksp);
110: PetscFinalize();
111: return 0;
112: }
116: int FormJacobian_Grid(GridCtx *grid,Mat *J)
117: {
118: Mat jac = *J;
120: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
121: PetscInt nloc,*ltog,grow;
122: PetscScalar two = 2.0,one = 1.0,v[5],hx,hy,hxdhy,hydhx,value;
124: mx = grid->mx; my = grid->my;
125: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
126: hxdhy = hx/hy; hydhx = hy/hx;
128: /* Get ghost points */
129: DAGetCorners(grid->da,&xs,&ys,0,&xm,&ym,0);
130: DAGetGhostCorners(grid->da,&Xs,&Ys,0,&Xm,&Ym,0);
131: DAGetGlobalIndices(grid->da,&nloc,<og);
133: /* Evaluate Jacobian of function */
134: for (j=ys; j<ys+ym; j++) {
135: row = (j - Ys)*Xm + xs - Xs - 1;
136: for (i=xs; i<xs+xm; i++) {
137: row++;
138: grow = ltog[row];
139: if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
140: v[0] = -hxdhy; col[0] = ltog[row - Xm];
141: v[1] = -hydhx; col[1] = ltog[row - 1];
142: v[2] = two*(hydhx + hxdhy); col[2] = grow;
143: v[3] = -hydhx; col[3] = ltog[row + 1];
144: v[4] = -hxdhy; col[4] = ltog[row + Xm];
145: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
146: } else if ((i > 0 && i < mx-1) || (j > 0 && j < my-1)){
147: value = .5*two*(hydhx + hxdhy);
148: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
149: } else {
150: value = .25*two*(hydhx + hxdhy);
151: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
152: }
153: }
154: }
155: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
156: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
157: return 0;
158: }