1 /*T
2    Concepts: KSP^solving a system of linear equations
3    Concepts: KSP^Laplacian, 2d
4    Processors: n
5 T*/
6 
7 /*
8 Added at the request of Marc Garbey.
9 
10 Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation
11 
12    -div \rho grad u = f,  0 < x,y < 1,
13 
14 with forcing function
15 
16    f = e^{-x^2/\nu} e^{-y^2/\nu}
17 
18 with Dirichlet boundary conditions
19 
20    u = f(x,y) for x = 0, x = 1, y = 0, y = 1
21 
22 or pure Neumman boundary conditions
23 
24 This uses multigrid to solve the linear system
25 */
26 
27 static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";
28 
29 #include <petscdm.h>
30 #include <petscdmda.h>
31 #include <petscksp.h>
32 
33 extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
34 extern PetscErrorCode ComputeRHS(KSP,Vec,void*);
35 
36 typedef enum {DIRICHLET, NEUMANN} BCType;
37 
38 typedef struct {
39   PetscReal rho;
40   PetscReal nu;
41   BCType    bcType;
42 } UserContext;
43 
main(int argc,char ** argv)44 int main(int argc,char **argv)
45 {
46   KSP            ksp;
47   DM             da;
48   UserContext    user;
49   const char     *bcTypes[2] = {"dirichlet","neumann"};
50   PetscErrorCode ierr;
51   PetscInt       bc;
52   Vec            b,x;
53   PetscBool      testsolver = PETSC_FALSE;
54 
55   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
56   ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
57   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,3,3,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da);CHKERRQ(ierr);
58   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
59   ierr = DMSetUp(da);CHKERRQ(ierr);
60   ierr = DMDASetUniformCoordinates(da,0,1,0,1,0,0);CHKERRQ(ierr);
61   ierr = DMDASetFieldName(da,0,"Pressure");CHKERRQ(ierr);
62 
63   ierr        = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");CHKERRQ(ierr);
64   user.rho    = 1.0;
65   ierr        = PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL);CHKERRQ(ierr);
66   user.nu     = 0.1;
67   ierr        = PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL);CHKERRQ(ierr);
68   bc          = (PetscInt)DIRICHLET;
69   ierr        = PetscOptionsEList("-bc_type","Type of boundary condition","ex29.c",bcTypes,2,bcTypes[0],&bc,NULL);CHKERRQ(ierr);
70   user.bcType = (BCType)bc;
71   ierr        = PetscOptionsBool("-testsolver", "Run solver multiple times, useful for performance studies of solver", "ex29.c", testsolver, &testsolver, NULL);CHKERRQ(ierr);
72   ierr        = PetscOptionsEnd();CHKERRQ(ierr);
73 
74   ierr = KSPSetComputeRHS(ksp,ComputeRHS,&user);CHKERRQ(ierr);
75   ierr = KSPSetComputeOperators(ksp,ComputeMatrix,&user);CHKERRQ(ierr);
76   ierr = KSPSetDM(ksp,da);CHKERRQ(ierr);
77   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
78   ierr = KSPSetUp(ksp);CHKERRQ(ierr);
79   ierr = KSPSolve(ksp,NULL,NULL);CHKERRQ(ierr);
80 
81   if (testsolver) {
82     ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
83     ierr = KSPGetRhs(ksp,&b);CHKERRQ(ierr);
84     KSPSetDMActive(ksp,PETSC_FALSE);
85     ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
86     {
87 #if defined(PETSC_USE_LOG)
88       PetscLogStage stage;
89 #endif
90       PetscInt      i,n = 20;
91 
92       ierr = PetscLogStageRegister("Solve only",&stage);CHKERRQ(ierr);
93       ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
94       for (i=0; i<n; i++) {
95         ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
96       }
97       ierr = PetscLogStagePop();CHKERRQ(ierr);
98     }
99   }
100 
101   ierr = DMDestroy(&da);CHKERRQ(ierr);
102   ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
103   ierr = PetscFinalize();
104   return ierr;
105 }
106 
ComputeRHS(KSP ksp,Vec b,void * ctx)107 PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
108 {
109   UserContext    *user = (UserContext*)ctx;
110   PetscErrorCode ierr;
111   PetscInt       i,j,mx,my,xm,ym,xs,ys;
112   PetscScalar    Hx,Hy;
113   PetscScalar    **array;
114   DM             da;
115 
116   PetscFunctionBeginUser;
117   ierr = KSPGetDM(ksp,&da);CHKERRQ(ierr);
118   ierr = DMDAGetInfo(da, 0, &mx, &my, 0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
119   Hx   = 1.0 / (PetscReal)(mx-1);
120   Hy   = 1.0 / (PetscReal)(my-1);
121   ierr = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
122   ierr = DMDAVecGetArray(da, b, &array);CHKERRQ(ierr);
123   for (j=ys; j<ys+ym; j++) {
124     for (i=xs; i<xs+xm; i++) {
125       array[j][i] = PetscExpScalar(-((PetscReal)i*Hx)*((PetscReal)i*Hx)/user->nu)*PetscExpScalar(-((PetscReal)j*Hy)*((PetscReal)j*Hy)/user->nu)*Hx*Hy;
126     }
127   }
128   ierr = DMDAVecRestoreArray(da, b, &array);CHKERRQ(ierr);
129   ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
130   ierr = VecAssemblyEnd(b);CHKERRQ(ierr);
131 
132   /* force right hand side to be consistent for singular matrix */
133   /* note this is really a hack, normally the model would provide you with a consistent right handside */
134   if (user->bcType == NEUMANN) {
135     MatNullSpace nullspace;
136 
137     ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr);
138     ierr = MatNullSpaceRemove(nullspace,b);CHKERRQ(ierr);
139     ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
140   }
141   PetscFunctionReturn(0);
142 }
143 
144 
ComputeRho(PetscInt i,PetscInt j,PetscInt mx,PetscInt my,PetscReal centerRho,PetscReal * rho)145 PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
146 {
147   PetscFunctionBeginUser;
148   if ((i > mx/3.0) && (i < 2.0*mx/3.0) && (j > my/3.0) && (j < 2.0*my/3.0)) {
149     *rho = centerRho;
150   } else {
151     *rho = 1.0;
152   }
153   PetscFunctionReturn(0);
154 }
155 
ComputeMatrix(KSP ksp,Mat J,Mat jac,void * ctx)156 PetscErrorCode ComputeMatrix(KSP ksp,Mat J,Mat jac,void *ctx)
157 {
158   UserContext    *user = (UserContext*)ctx;
159   PetscReal      centerRho;
160   PetscErrorCode ierr;
161   PetscInt       i,j,mx,my,xm,ym,xs,ys;
162   PetscScalar    v[5];
163   PetscReal      Hx,Hy,HydHx,HxdHy,rho;
164   MatStencil     row, col[5];
165   DM             da;
166   PetscBool      check_matis = PETSC_FALSE;
167 
168   PetscFunctionBeginUser;
169   ierr      = KSPGetDM(ksp,&da);CHKERRQ(ierr);
170   centerRho = user->rho;
171   ierr      = DMDAGetInfo(da,0,&mx,&my,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
172   Hx        = 1.0 / (PetscReal)(mx-1);
173   Hy        = 1.0 / (PetscReal)(my-1);
174   HxdHy     = Hx/Hy;
175   HydHx     = Hy/Hx;
176   ierr      = DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
177   for (j=ys; j<ys+ym; j++) {
178     for (i=xs; i<xs+xm; i++) {
179       row.i = i; row.j = j;
180       ierr  = ComputeRho(i, j, mx, my, centerRho, &rho);CHKERRQ(ierr);
181       if (i==0 || j==0 || i==mx-1 || j==my-1) {
182         if (user->bcType == DIRICHLET) {
183           v[0] = 2.0*rho*(HxdHy + HydHx);
184           ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
185         } else if (user->bcType == NEUMANN) {
186           PetscInt numx = 0, numy = 0, num = 0;
187           if (j!=0) {
188             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j-1;
189             numy++; num++;
190           }
191           if (i!=0) {
192             v[num] = -rho*HydHx;              col[num].i = i-1; col[num].j = j;
193             numx++; num++;
194           }
195           if (i!=mx-1) {
196             v[num] = -rho*HydHx;              col[num].i = i+1; col[num].j = j;
197             numx++; num++;
198           }
199           if (j!=my-1) {
200             v[num] = -rho*HxdHy;              col[num].i = i;   col[num].j = j+1;
201             numy++; num++;
202           }
203           v[num] = numx*rho*HydHx + numy*rho*HxdHy; col[num].i = i;   col[num].j = j;
204           num++;
205           ierr = MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);CHKERRQ(ierr);
206         }
207       } else {
208         v[0] = -rho*HxdHy;              col[0].i = i;   col[0].j = j-1;
209         v[1] = -rho*HydHx;              col[1].i = i-1; col[1].j = j;
210         v[2] = 2.0*rho*(HxdHy + HydHx); col[2].i = i;   col[2].j = j;
211         v[3] = -rho*HydHx;              col[3].i = i+1; col[3].j = j;
212         v[4] = -rho*HxdHy;              col[4].i = i;   col[4].j = j+1;
213         ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
214       }
215     }
216   }
217   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
218   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
219   ierr = MatViewFromOptions(jac,NULL,"-view_mat");CHKERRQ(ierr);
220   ierr = PetscOptionsGetBool(NULL,NULL,"-check_matis",&check_matis,NULL);CHKERRQ(ierr);
221   if (check_matis) {
222     void      (*f)(void);
223     Mat       J2;
224     MatType   jtype;
225     PetscReal nrm;
226 
227     ierr = MatGetType(jac,&jtype);CHKERRQ(ierr);
228     ierr = MatConvert(jac,MATIS,MAT_INITIAL_MATRIX,&J2);CHKERRQ(ierr);
229     ierr = MatViewFromOptions(J2,NULL,"-view_conv");CHKERRQ(ierr);
230     ierr = MatConvert(J2,jtype,MAT_INPLACE_MATRIX,&J2);CHKERRQ(ierr);
231     ierr = MatGetOperation(jac,MATOP_VIEW,&f);CHKERRQ(ierr);
232     ierr = MatSetOperation(J2,MATOP_VIEW,f);CHKERRQ(ierr);
233     ierr = MatSetDM(J2,da);CHKERRQ(ierr);
234     ierr = MatViewFromOptions(J2,NULL,"-view_conv_assembled");CHKERRQ(ierr);
235     ierr = MatAXPY(J2,-1.,jac,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
236     ierr = MatNorm(J2,NORM_FROBENIUS,&nrm);CHKERRQ(ierr);
237     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MATIS %g\n",(double)nrm);CHKERRQ(ierr);
238     ierr = MatViewFromOptions(J2,NULL,"-view_conv_err");CHKERRQ(ierr);
239     ierr = MatDestroy(&J2);CHKERRQ(ierr);
240   }
241   if (user->bcType == NEUMANN) {
242     MatNullSpace nullspace;
243 
244     ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);CHKERRQ(ierr);
245     ierr = MatSetNullSpace(J,nullspace);CHKERRQ(ierr);
246     ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
247   }
248   PetscFunctionReturn(0);
249 }
250 
251 
252 /*TEST
253 
254    test:
255       args: -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -ksp_rtol 1.e-3
256 
257    test:
258       suffix: 2
259       args: -bc_type neumann -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -mg_coarse_pc_factor_shift_type nonzero
260       requires: !single
261 
262    test:
263       suffix: telescope
264       nsize: 4
265       args: -ksp_monitor_short -da_grid_x 257 -da_grid_y 257 -pc_type mg -pc_mg_galerkin pmat -pc_mg_levels 4 -ksp_type richardson -mg_levels_ksp_type chebyshev -mg_levels_pc_type jacobi -mg_coarse_pc_type telescope -mg_coarse_pc_telescope_ignore_kspcomputeoperators -mg_coarse_telescope_pc_type mg -mg_coarse_telescope_pc_mg_galerkin pmat -mg_coarse_telescope_pc_mg_levels 3 -mg_coarse_telescope_mg_levels_ksp_type chebyshev -mg_coarse_telescope_mg_levels_pc_type jacobi -mg_coarse_pc_telescope_reduction_factor 4
266 
267    test:
268       suffix: 3
269       args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_pc_type jacobi
270 
271    test:
272       suffix: 4
273       args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_ksp_max_it 3 -mg_levels_ksp_max_it 4
274 
275    test:
276       suffix: 5
277       nsize: 2
278       requires: hypre !complex
279       args: -pc_type mg  -da_refine 2 -ksp_monitor  -matptap_via hypre -pc_mg_galerkin both
280 
281 TEST*/
282