Actual source code: ex19.c
petsc-3.4.2 2013-07-02
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
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.
21: */
23: #include <petscksp.h>
24: #include <petscdmda.h>
26: /* User-defined application contexts */
28: typedef struct {
29: PetscInt mx,my; /* number grid points in x and y direction */
30: Vec localX,localF; /* local vectors with ghost region */
31: DM da;
32: Vec x,b,r; /* global vectors */
33: Mat J; /* Jacobian on grid */
34: } GridCtx;
36: typedef struct {
37: GridCtx fine;
38: GridCtx coarse;
39: KSP ksp_coarse;
40: PetscInt ratio;
41: Mat Ii; /* interpolation from coarse to fine */
42: } AppCtx;
44: #define COARSE_LEVEL 0
45: #define FINE_LEVEL 1
47: extern int FormJacobian_Grid(AppCtx*,GridCtx*,Mat*);
49: /*
50: Mm_ratio - ration of grid lines between fine and coarse grids.
51: */
54: int main(int argc,char **argv)
55: {
56: AppCtx user;
58: PetscInt its,N,n,Nx = PETSC_DECIDE,Ny = PETSC_DECIDE,nlocal,Nlocal;
59: PetscMPIInt size;
60: KSP ksp,ksp_fine;
61: PC pc;
62: PetscScalar one = 1.0;
64: PetscInitialize(&argc,&argv,NULL,help);
66: user.ratio = 2;
67: user.coarse.mx = 5; user.coarse.my = 5;
69: PetscOptionsGetInt(NULL,"-Mx",&user.coarse.mx,NULL);
70: PetscOptionsGetInt(NULL,"-My",&user.coarse.my,NULL);
71: PetscOptionsGetInt(NULL,"-ratio",&user.ratio,NULL);
73: user.fine.mx = user.ratio*(user.coarse.mx-1)+1; user.fine.my = user.ratio*(user.coarse.my-1)+1;
75: PetscPrintf(PETSC_COMM_WORLD,"Coarse grid size %D by %D\n",user.coarse.mx,user.coarse.my);
76: PetscPrintf(PETSC_COMM_WORLD,"Fine grid size %D by %D\n",user.fine.mx,user.fine.my);
78: n = user.fine.mx*user.fine.my; N = user.coarse.mx*user.coarse.my;
80: MPI_Comm_size(PETSC_COMM_WORLD,&size);
81: PetscOptionsGetInt(NULL,"-Nx",&Nx,NULL);
82: PetscOptionsGetInt(NULL,"-Ny",&Ny,NULL);
84: /* Set up distributed array for fine grid */
85: DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.fine.mx,
86: user.fine.my,Nx,Ny,1,1,NULL,NULL,&user.fine.da);
87: DMCreateGlobalVector(user.fine.da,&user.fine.x);
88: VecDuplicate(user.fine.x,&user.fine.r);
89: VecDuplicate(user.fine.x,&user.fine.b);
90: VecGetLocalSize(user.fine.x,&nlocal);
91: DMCreateLocalVector(user.fine.da,&user.fine.localX);
92: VecDuplicate(user.fine.localX,&user.fine.localF);
93: MatCreateAIJ(PETSC_COMM_WORLD,nlocal,nlocal,n,n,5,NULL,3,NULL,&user.fine.J);
95: /* Set up distributed array for coarse grid */
96: DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.coarse.mx,
97: user.coarse.my,Nx,Ny,1,1,NULL,NULL,&user.coarse.da);
98: DMCreateGlobalVector(user.coarse.da,&user.coarse.x);
99: VecDuplicate(user.coarse.x,&user.coarse.b);
100: VecGetLocalSize(user.coarse.x,&Nlocal);
101: DMCreateLocalVector(user.coarse.da,&user.coarse.localX);
102: VecDuplicate(user.coarse.localX,&user.coarse.localF);
103: MatCreateAIJ(PETSC_COMM_WORLD,Nlocal,Nlocal,N,N,5,NULL,3,NULL,&user.coarse.J);
105: /* Create linear solver */
106: KSPCreate(PETSC_COMM_WORLD,&ksp);
108: /* set two level additive Schwarz preconditioner */
109: KSPGetPC(ksp,&pc);
110: PCSetType(pc,PCMG);
111: PCMGSetLevels(pc,2,NULL);
112: PCMGSetType(pc,PC_MG_ADDITIVE);
114: FormJacobian_Grid(&user,&user.coarse,&user.coarse.J);
115: FormJacobian_Grid(&user,&user.fine,&user.fine.J);
117: /* Create coarse level */
118: PCMGGetCoarseSolve(pc,&user.ksp_coarse);
119: KSPSetOptionsPrefix(user.ksp_coarse,"coarse_");
120: KSPSetFromOptions(user.ksp_coarse);
121: KSPSetOperators(user.ksp_coarse,user.coarse.J,user.coarse.J,DIFFERENT_NONZERO_PATTERN);
122: PCMGSetX(pc,COARSE_LEVEL,user.coarse.x);
123: PCMGSetRhs(pc,COARSE_LEVEL,user.coarse.b);
125: /* Create fine level */
126: PCMGGetSmoother(pc,FINE_LEVEL,&ksp_fine);
127: KSPSetOptionsPrefix(ksp_fine,"fine_");
128: KSPSetFromOptions(ksp_fine);
129: KSPSetOperators(ksp_fine,user.fine.J,user.fine.J,DIFFERENT_NONZERO_PATTERN);
130: PCMGSetR(pc,FINE_LEVEL,user.fine.r);
132: /* Create interpolation between the levels */
133: DMCreateInterpolation(user.coarse.da,user.fine.da,&user.Ii,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: DMDestroy(&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: DMDestroy(&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: DMDAGetCorners(grid->da,&xs,&ys,0,&xm,&ym,0);
194: DMDAGetGhostCorners(grid->da,&Xs,&Ys,0,&Xm,&Ym,0);
195: DMDAGetGlobalIndices(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: }