Actual source code: ex19.c
2: static char help[] ="Solvers Laplacian with multigrid, bad way.\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.
21: */
23: #include petscksp.h
24: #include petscda.h
25: #include petscmg.h
27: /* User-defined application contexts */
29: typedef struct {
30: PetscInt mx,my; /* number grid points in x and y direction */
31: Vec localX,localF; /* local vectors with ghost region */
32: DA da;
33: Vec x,b,r; /* global vectors */
34: Mat J; /* Jacobian on grid */
35: } GridCtx;
37: typedef struct {
38: GridCtx fine;
39: GridCtx coarse;
40: KSP ksp_coarse;
41: PetscInt ratio;
42: Mat Ii; /* interpolation from coarse to fine */
43: } AppCtx;
45: #define COARSE_LEVEL 0
46: #define FINE_LEVEL 1
50: /*
51: Mm_ratio - ration of grid lines between fine and coarse grids.
52: */
55: int main(int argc,char **argv)
56: {
57: AppCtx user;
59: PetscInt its,N,n,Nx = PETSC_DECIDE,Ny = PETSC_DECIDE,nlocal,Nlocal;
60: PetscMPIInt size;
61: KSP ksp,ksp_fine;
62: PC pc;
63: PetscScalar one = 1.0;
65: PetscInitialize(&argc,&argv,PETSC_NULL,help);
67: user.ratio = 2;
68: user.coarse.mx = 5; user.coarse.my = 5;
69: PetscOptionsGetInt(PETSC_NULL,"-Mx",&user.coarse.mx,PETSC_NULL);
70: PetscOptionsGetInt(PETSC_NULL,"-My",&user.coarse.my,PETSC_NULL);
71: PetscOptionsGetInt(PETSC_NULL,"-ratio",&user.ratio,PETSC_NULL);
72: user.fine.mx = user.ratio*(user.coarse.mx-1)+1; user.fine.my = user.ratio*(user.coarse.my-1)+1;
74: PetscPrintf(PETSC_COMM_WORLD,"Coarse grid size %D by %D\n",user.coarse.mx,user.coarse.my);
75: PetscPrintf(PETSC_COMM_WORLD,"Fine grid size %D by %D\n",user.fine.mx,user.fine.my);
77: n = user.fine.mx*user.fine.my; N = user.coarse.mx*user.coarse.my;
79: MPI_Comm_size(PETSC_COMM_WORLD,&size);
80: PetscOptionsGetInt(PETSC_NULL,"-Nx",&Nx,PETSC_NULL);
81: PetscOptionsGetInt(PETSC_NULL,"-Ny",&Ny,PETSC_NULL);
83: /* Set up distributed array for fine grid */
84: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.fine.mx,
85: user.fine.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.fine.da);
86: DACreateGlobalVector(user.fine.da,&user.fine.x);
87: VecDuplicate(user.fine.x,&user.fine.r);
88: VecDuplicate(user.fine.x,&user.fine.b);
89: VecGetLocalSize(user.fine.x,&nlocal);
90: DACreateLocalVector(user.fine.da,&user.fine.localX);
91: VecDuplicate(user.fine.localX,&user.fine.localF);
92: MatCreateMPIAIJ(PETSC_COMM_WORLD,nlocal,nlocal,n,n,5,PETSC_NULL,3,PETSC_NULL,&user.fine.J);
94: /* Set up distributed array for coarse grid */
95: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,user.coarse.mx,
96: user.coarse.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&user.coarse.da);
97: DACreateGlobalVector(user.coarse.da,&user.coarse.x);
98: VecDuplicate(user.coarse.x,&user.coarse.b);
99: VecGetLocalSize(user.coarse.x,&Nlocal);
100: DACreateLocalVector(user.coarse.da,&user.coarse.localX);
101: VecDuplicate(user.coarse.localX,&user.coarse.localF);
102: MatCreateMPIAIJ(PETSC_COMM_WORLD,Nlocal,Nlocal,N,N,5,PETSC_NULL,3,PETSC_NULL,&user.coarse.J);
104: /* Create linear solver */
105: KSPCreate(PETSC_COMM_WORLD,&ksp);
107: /* set two level additive Schwarz preconditioner */
108: KSPGetPC(ksp,&pc);
109: PCSetType(pc,PCMG);
110: PCMGSetLevels(pc,2,PETSC_NULL);
111: PCMGSetType(pc,PC_MG_ADDITIVE);
113: FormJacobian_Grid(&user,&user.coarse,&user.coarse.J);
114: FormJacobian_Grid(&user,&user.fine,&user.fine.J);
116: /* Create coarse level */
117: PCMGGetCoarseSolve(pc,&user.ksp_coarse);
118: KSPSetOptionsPrefix(user.ksp_coarse,"coarse_");
119: KSPSetFromOptions(user.ksp_coarse);
120: KSPSetOperators(user.ksp_coarse,user.coarse.J,user.coarse.J,DIFFERENT_NONZERO_PATTERN);
121: PCMGSetX(pc,COARSE_LEVEL,user.coarse.x);
122: PCMGSetRhs(pc,COARSE_LEVEL,user.coarse.b);
124: /* Create fine level */
125: PCMGGetSmoother(pc,FINE_LEVEL,&ksp_fine);
126: KSPSetOptionsPrefix(ksp_fine,"fine_");
127: KSPSetFromOptions(ksp_fine);
128: KSPSetOperators(ksp_fine,user.fine.J,user.fine.J,DIFFERENT_NONZERO_PATTERN);
129: PCMGSetR(pc,FINE_LEVEL,user.fine.r);
130: PCMGSetResidual(pc,FINE_LEVEL,PCMGDefaultResidual,user.fine.J);
132: /* Create interpolation between the levels */
133: DAGetInterpolation(user.coarse.da,user.fine.da,&user.Ii,PETSC_NULL);
134: PCMGSetInterpolation(pc,FINE_LEVEL,user.Ii);
135: PCMGSetRestriction(pc,FINE_LEVEL,user.Ii);
137: KSPSetOperators(ksp,user.fine.J,user.fine.J,DIFFERENT_NONZERO_PATTERN);
139: VecSet(user.fine.b,one);
140: {
141: PetscRandom rdm;
142: PetscRandomCreate(PETSC_COMM_WORLD,&rdm);
143: PetscRandomSetFromOptions(rdm);
144: VecSetRandom(user.fine.b,rdm);
145: PetscRandomDestroy(rdm);
146: }
148: /* Set options, then solve nonlinear system */
149: KSPSetFromOptions(ksp);
151: KSPSolve(ksp,user.fine.b,user.fine.x);
152: KSPGetIterationNumber(ksp,&its);
153: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %D\n",its);
155: /* Free data structures */
156: MatDestroy(user.fine.J);
157: VecDestroy(user.fine.x);
158: VecDestroy(user.fine.r);
159: VecDestroy(user.fine.b);
160: DADestroy(user.fine.da);
161: VecDestroy(user.fine.localX);
162: VecDestroy(user.fine.localF);
164: MatDestroy(user.coarse.J);
165: VecDestroy(user.coarse.x);
166: VecDestroy(user.coarse.b);
167: DADestroy(user.coarse.da);
168: VecDestroy(user.coarse.localX);
169: VecDestroy(user.coarse.localF);
171: KSPDestroy(ksp);
172: MatDestroy(user.Ii);
173: PetscFinalize();
175: return 0;
176: }
180: int FormJacobian_Grid(AppCtx *user,GridCtx *grid,Mat *J)
181: {
182: Mat jac = *J;
184: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
185: PetscInt nloc,*ltog,grow;
186: PetscScalar two = 2.0,one = 1.0,v[5],hx,hy,hxdhy,hydhx,value;
188: mx = grid->mx; my = grid->my;
189: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
190: hxdhy = hx/hy; hydhx = hy/hx;
192: /* Get ghost points */
193: DAGetCorners(grid->da,&xs,&ys,0,&xm,&ym,0);
194: DAGetGhostCorners(grid->da,&Xs,&Ys,0,&Xm,&Ym,0);
195: DAGetGlobalIndices(grid->da,&nloc,<og);
197: /* Evaluate Jacobian of function */
198: for (j=ys; j<ys+ym; j++) {
199: row = (j - Ys)*Xm + xs - Xs - 1;
200: for (i=xs; i<xs+xm; i++) {
201: row++;
202: grow = ltog[row];
203: if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
204: v[0] = -hxdhy; col[0] = ltog[row - Xm];
205: v[1] = -hydhx; col[1] = ltog[row - 1];
206: v[2] = two*(hydhx + hxdhy); col[2] = grow;
207: v[3] = -hydhx; col[3] = ltog[row + 1];
208: v[4] = -hxdhy; col[4] = ltog[row + Xm];
209: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
210: } else if ((i > 0 && i < mx-1) || (j > 0 && j < my-1)){
211: value = .5*two*(hydhx + hxdhy);
212: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
213: } else {
214: value = .25*two*(hydhx + hxdhy);
215: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
216: }
217: }
218: }
219: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
220: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
222: return 0;
223: }