Actual source code: ex19.c
petsc-3.6.4 2016-04-12
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 <petscdm.h>
25: #include <petscdmda.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: DM 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
48: extern int FormJacobian_Grid(AppCtx*,GridCtx*,Mat*);
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,NULL,help);
67: user.ratio = 2;
68: user.coarse.mx = 5; user.coarse.my = 5;
70: PetscOptionsGetInt(NULL,"-Mx",&user.coarse.mx,NULL);
71: PetscOptionsGetInt(NULL,"-My",&user.coarse.my,NULL);
72: PetscOptionsGetInt(NULL,"-ratio",&user.ratio,NULL);
74: user.fine.mx = user.ratio*(user.coarse.mx-1)+1; user.fine.my = user.ratio*(user.coarse.my-1)+1;
76: PetscPrintf(PETSC_COMM_WORLD,"Coarse grid size %D by %D\n",user.coarse.mx,user.coarse.my);
77: PetscPrintf(PETSC_COMM_WORLD,"Fine grid size %D by %D\n",user.fine.mx,user.fine.my);
79: n = user.fine.mx*user.fine.my; N = user.coarse.mx*user.coarse.my;
81: MPI_Comm_size(PETSC_COMM_WORLD,&size);
82: PetscOptionsGetInt(NULL,"-Nx",&Nx,NULL);
83: PetscOptionsGetInt(NULL,"-Ny",&Ny,NULL);
85: /* Set up distributed array for fine grid */
86: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.fine.mx,
87: user.fine.my,Nx,Ny,1,1,NULL,NULL,&user.fine.da);
88: DMCreateGlobalVector(user.fine.da,&user.fine.x);
89: VecDuplicate(user.fine.x,&user.fine.r);
90: VecDuplicate(user.fine.x,&user.fine.b);
91: VecGetLocalSize(user.fine.x,&nlocal);
92: DMCreateLocalVector(user.fine.da,&user.fine.localX);
93: VecDuplicate(user.fine.localX,&user.fine.localF);
94: MatCreateAIJ(PETSC_COMM_WORLD,nlocal,nlocal,n,n,5,NULL,3,NULL,&user.fine.J);
96: /* Set up distributed array for coarse grid */
97: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.coarse.mx,
98: user.coarse.my,Nx,Ny,1,1,NULL,NULL,&user.coarse.da);
99: DMCreateGlobalVector(user.coarse.da,&user.coarse.x);
100: VecDuplicate(user.coarse.x,&user.coarse.b);
101: VecGetLocalSize(user.coarse.x,&Nlocal);
102: DMCreateLocalVector(user.coarse.da,&user.coarse.localX);
103: VecDuplicate(user.coarse.localX,&user.coarse.localF);
104: MatCreateAIJ(PETSC_COMM_WORLD,Nlocal,Nlocal,N,N,5,NULL,3,NULL,&user.coarse.J);
106: /* Create linear solver */
107: KSPCreate(PETSC_COMM_WORLD,&ksp);
109: /* set two level additive Schwarz preconditioner */
110: KSPGetPC(ksp,&pc);
111: PCSetType(pc,PCMG);
112: PCMGSetLevels(pc,2,NULL);
113: PCMGSetType(pc,PC_MG_ADDITIVE);
115: FormJacobian_Grid(&user,&user.coarse,&user.coarse.J);
116: FormJacobian_Grid(&user,&user.fine,&user.fine.J);
118: /* Create coarse level */
119: PCMGGetCoarseSolve(pc,&user.ksp_coarse);
120: KSPSetOptionsPrefix(user.ksp_coarse,"coarse_");
121: KSPSetFromOptions(user.ksp_coarse);
122: KSPSetOperators(user.ksp_coarse,user.coarse.J,user.coarse.J);
123: PCMGSetX(pc,COARSE_LEVEL,user.coarse.x);
124: PCMGSetRhs(pc,COARSE_LEVEL,user.coarse.b);
126: /* Create fine level */
127: PCMGGetSmoother(pc,FINE_LEVEL,&ksp_fine);
128: KSPSetOptionsPrefix(ksp_fine,"fine_");
129: KSPSetFromOptions(ksp_fine);
130: KSPSetOperators(ksp_fine,user.fine.J,user.fine.J);
131: PCMGSetR(pc,FINE_LEVEL,user.fine.r);
133: /* Create interpolation between the levels */
134: DMCreateInterpolation(user.coarse.da,user.fine.da,&user.Ii,NULL);
135: PCMGSetInterpolation(pc,FINE_LEVEL,user.Ii);
136: PCMGSetRestriction(pc,FINE_LEVEL,user.Ii);
138: KSPSetOperators(ksp,user.fine.J,user.fine.J);
140: VecSet(user.fine.b,one);
141: {
142: PetscRandom rdm;
143: PetscRandomCreate(PETSC_COMM_WORLD,&rdm);
144: PetscRandomSetFromOptions(rdm);
145: VecSetRandom(user.fine.b,rdm);
146: PetscRandomDestroy(&rdm);
147: }
149: /* Set options, then solve nonlinear system */
150: KSPSetFromOptions(ksp);
152: KSPSolve(ksp,user.fine.b,user.fine.x);
153: KSPGetIterationNumber(ksp,&its);
154: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %D\n",its);
156: /* Free data structures */
157: MatDestroy(&user.fine.J);
158: VecDestroy(&user.fine.x);
159: VecDestroy(&user.fine.r);
160: VecDestroy(&user.fine.b);
161: DMDestroy(&user.fine.da);
162: VecDestroy(&user.fine.localX);
163: VecDestroy(&user.fine.localF);
165: MatDestroy(&user.coarse.J);
166: VecDestroy(&user.coarse.x);
167: VecDestroy(&user.coarse.b);
168: DMDestroy(&user.coarse.da);
169: VecDestroy(&user.coarse.localX);
170: VecDestroy(&user.coarse.localF);
172: KSPDestroy(&ksp);
173: MatDestroy(&user.Ii);
174: PetscFinalize();
176: return 0;
177: }
181: int FormJacobian_Grid(AppCtx *user,GridCtx *grid,Mat *J)
182: {
183: Mat jac = *J;
184: PetscErrorCode ierr;
185: PetscInt i,j,row,mx,my,xs,ys,xm,ym,Xs,Ys,Xm,Ym,col[5];
186: PetscInt grow;
187: const PetscInt *ltog;
188: PetscScalar two = 2.0,one = 1.0,v[5],hx,hy,hxdhy,hydhx,value;
189: ISLocalToGlobalMapping ltogm;
191: mx = grid->mx; my = grid->my;
192: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
193: hxdhy = hx/hy; hydhx = hy/hx;
195: /* Get ghost points */
196: DMDAGetCorners(grid->da,&xs,&ys,0,&xm,&ym,0);
197: DMDAGetGhostCorners(grid->da,&Xs,&Ys,0,&Xm,&Ym,0);
198: DMGetLocalToGlobalMapping(grid->da,<ogm);
199: ISLocalToGlobalMappingGetIndices(ltogm,<og);
201: /* Evaluate Jacobian of function */
202: for (j=ys; j<ys+ym; j++) {
203: row = (j - Ys)*Xm + xs - Xs - 1;
204: for (i=xs; i<xs+xm; i++) {
205: row++;
206: grow = ltog[row];
207: if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
208: v[0] = -hxdhy; col[0] = ltog[row - Xm];
209: v[1] = -hydhx; col[1] = ltog[row - 1];
210: v[2] = two*(hydhx + hxdhy); col[2] = grow;
211: v[3] = -hydhx; col[3] = ltog[row + 1];
212: v[4] = -hxdhy; col[4] = ltog[row + Xm];
213: MatSetValues(jac,1,&grow,5,col,v,INSERT_VALUES);
214: } else if ((i > 0 && i < mx-1) || (j > 0 && j < my-1)) {
215: value = .5*two*(hydhx + hxdhy);
216: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
217: } else {
218: value = .25*two*(hydhx + hxdhy);
219: MatSetValues(jac,1,&grow,1,&grow,&value,INSERT_VALUES);
220: }
221: }
222: }
223: ISLocalToGlobalMappingRestoreIndices(ltogm,<og);
224: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
225: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
227: return 0;
228: }