1 /* DMDA/KSP solving a system of linear equations.
2 Poisson equation in 2D:
3
4 div(grad p) = f, 0 < x,y < 1
5 with
6 forcing function f = -cos(m*pi*x)*cos(n*pi*y),
7 Neuman boundary conditions
8 dp/dx = 0 for x = 0, x = 1.
9 dp/dy = 0 for y = 0, y = 1.
10
11 Contributed by Michael Boghosian <boghmic@iit.edu>, 2008,
12 based on petsc/src/ksp/ksp/tutorials/ex29.c and ex32.c
13
14 Compare to ex66.c
15
16 Example of Usage:
17 ./ex50 -da_grid_x 3 -da_grid_y 3 -pc_type mg -da_refine 3 -ksp_monitor -ksp_view -dm_view draw -draw_pause -1
18 ./ex50 -da_grid_x 100 -da_grid_y 100 -pc_type mg -pc_mg_levels 1 -mg_levels_0_pc_type ilu -mg_levels_0_pc_factor_levels 1 -ksp_monitor -ksp_view
19 ./ex50 -da_grid_x 100 -da_grid_y 100 -pc_type mg -pc_mg_levels 1 -mg_levels_0_pc_type lu -mg_levels_0_pc_factor_shift_type NONZERO -ksp_monitor
20 mpiexec -n 4 ./ex50 -da_grid_x 3 -da_grid_y 3 -pc_type mg -da_refine 10 -ksp_monitor -ksp_view -log_view
21 */
22
23 static char help[] = "Solves 2D Poisson equation using multigrid.\n\n";
24
25 #include <petscdm.h>
26 #include <petscdmda.h>
27 #include <petscksp.h>
28 #include <petscsys.h>
29 #include <petscvec.h>
30
31 extern PetscErrorCode ComputeJacobian(KSP,Mat,Mat,void*);
32 extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
33
34 typedef struct {
35 PetscScalar uu, tt;
36 } UserContext;
37
main(int argc,char ** argv)38 int main(int argc,char **argv)
39 {
40 KSP ksp;
41 DM da;
42 UserContext user;
43 PetscErrorCode ierr;
44
45 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
46 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
47 ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
48 ierr = DMSetFromOptions(da);CHKERRQ(ierr);
49 ierr = DMSetUp(da);CHKERRQ(ierr);
50 ierr = KSPSetDM(ksp,(DM)da);CHKERRQ(ierr);
51 ierr = DMSetApplicationContext(da,&user);CHKERRQ(ierr);
52
53 user.uu = 1.0;
54 user.tt = 1.0;
55
56 ierr = KSPSetComputeRHS(ksp,ComputeRHS,&user);CHKERRQ(ierr);
57 ierr = KSPSetComputeOperators(ksp,ComputeJacobian,&user);CHKERRQ(ierr);
58 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
59 ierr = KSPSolve(ksp,NULL,NULL);CHKERRQ(ierr);
60
61 ierr = DMDestroy(&da);CHKERRQ(ierr);
62 ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
63 ierr = PetscFinalize();
64 return ierr;
65 }
66
ComputeRHS(KSP ksp,Vec b,void * ctx)67 PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
68 {
69 UserContext *user = (UserContext*)ctx;
70 PetscErrorCode ierr;
71 PetscInt i,j,M,N,xm,ym,xs,ys;
72 PetscScalar Hx,Hy,pi,uu,tt;
73 PetscScalar **array;
74 DM da;
75 MatNullSpace nullspace;
76
77 PetscFunctionBeginUser;
78 ierr = KSPGetDM(ksp,&da);CHKERRQ(ierr);
79 ierr = DMDAGetInfo(da, 0, &M, &N, 0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
80 uu = user->uu; tt = user->tt;
81 pi = 4*atan(1.0);
82 Hx = 1.0/(PetscReal)(M);
83 Hy = 1.0/(PetscReal)(N);
84
85 ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); /* Fine grid */
86 ierr = DMDAVecGetArray(da, b, &array);CHKERRQ(ierr);
87 for (j=ys; j<ys+ym; j++) {
88 for (i=xs; i<xs+xm; i++) {
89 array[j][i] = -PetscCosScalar(uu*pi*((PetscReal)i+0.5)*Hx)*PetscCosScalar(tt*pi*((PetscReal)j+0.5)*Hy)*Hx*Hy;
90 }
91 }
92 ierr = DMDAVecRestoreArray(da, b, &array);CHKERRQ(ierr);
93 ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
94 ierr = VecAssemblyEnd(b);CHKERRQ(ierr);
95
96 /* force right hand side to be consistent for singular matrix */
97 /* note this is really a hack, normally the model would provide you with a consistent right handside */
98 ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr);
99 ierr = MatNullSpaceRemove(nullspace,b);CHKERRQ(ierr);
100 ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
101 PetscFunctionReturn(0);
102 }
103
ComputeJacobian(KSP ksp,Mat J,Mat jac,void * ctx)104 PetscErrorCode ComputeJacobian(KSP ksp,Mat J, Mat jac,void *ctx)
105 {
106 PetscErrorCode ierr;
107 PetscInt i, j, M, N, xm, ym, xs, ys, num, numi, numj;
108 PetscScalar v[5], Hx, Hy, HydHx, HxdHy;
109 MatStencil row, col[5];
110 DM da;
111 MatNullSpace nullspace;
112
113 PetscFunctionBeginUser;
114 ierr = KSPGetDM(ksp,&da);CHKERRQ(ierr);
115 ierr = DMDAGetInfo(da,0,&M,&N,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
116 Hx = 1.0 / (PetscReal)(M);
117 Hy = 1.0 / (PetscReal)(N);
118 HxdHy = Hx/Hy;
119 HydHx = Hy/Hx;
120 ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
121 for (j=ys; j<ys+ym; j++) {
122 for (i=xs; i<xs+xm; i++) {
123 row.i = i; row.j = j;
124
125 if (i==0 || j==0 || i==M-1 || j==N-1) {
126 num=0; numi=0; numj=0;
127 if (j!=0) {
128 v[num] = -HxdHy; col[num].i = i; col[num].j = j-1;
129 num++; numj++;
130 }
131 if (i!=0) {
132 v[num] = -HydHx; col[num].i = i-1; col[num].j = j;
133 num++; numi++;
134 }
135 if (i!=M-1) {
136 v[num] = -HydHx; col[num].i = i+1; col[num].j = j;
137 num++; numi++;
138 }
139 if (j!=N-1) {
140 v[num] = -HxdHy; col[num].i = i; col[num].j = j+1;
141 num++; numj++;
142 }
143 v[num] = ((PetscReal)(numj)*HxdHy + (PetscReal)(numi)*HydHx); col[num].i = i; col[num].j = j;
144 num++;
145 ierr = MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);CHKERRQ(ierr);
146 } else {
147 v[0] = -HxdHy; col[0].i = i; col[0].j = j-1;
148 v[1] = -HydHx; col[1].i = i-1; col[1].j = j;
149 v[2] = 2.0*(HxdHy + HydHx); col[2].i = i; col[2].j = j;
150 v[3] = -HydHx; col[3].i = i+1; col[3].j = j;
151 v[4] = -HxdHy; col[4].i = i; col[4].j = j+1;
152 ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
153 }
154 }
155 }
156 ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157 ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158
159 ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr);
160 ierr = MatSetNullSpace(J,nullspace);CHKERRQ(ierr);
161 ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
162 PetscFunctionReturn(0);
163 }
164
165
166
167 /*TEST
168
169 build:
170 requires: !complex !single
171
172 test:
173 args: -pc_type mg -pc_mg_type full -ksp_type cg -ksp_monitor_short -da_refine 3 -mg_coarse_pc_type svd -ksp_view
174
175 test:
176 suffix: 2
177 nsize: 4
178 args: -pc_type mg -pc_mg_type full -ksp_type cg -ksp_monitor_short -da_refine 3 -mg_coarse_pc_type redundant -mg_coarse_redundant_pc_type svd -ksp_view
179
180 test:
181 suffix: 3
182 nsize: 2
183 args: -pc_type mg -pc_mg_type full -ksp_monitor_short -da_refine 5 -mg_coarse_ksp_type cg -mg_coarse_ksp_converged_reason -mg_coarse_ksp_rtol 1e-2 -mg_coarse_ksp_max_it 5 -mg_coarse_pc_type none -pc_mg_levels 2 -ksp_type pipefgmres -ksp_pipefgmres_shift 1.5
184
185 test:
186 suffix: tut_1
187 nsize: 1
188 args: -da_grid_x 4 -da_grid_y 4 -mat_view
189
190 test:
191 suffix: tut_2
192 requires: superlu_dist parmetis
193 nsize: 4
194 args: -da_grid_x 120 -da_grid_y 120 -pc_type lu -pc_factor_mat_solver_type superlu_dist -ksp_monitor -ksp_view
195
196 test:
197 suffix: tut_3
198 nsize: 4
199 args: -da_grid_x 1025 -da_grid_y 1025 -pc_type mg -pc_mg_levels 9 -ksp_monitor
200
201 TEST*/
202