Actual source code: ex29.c

  1: /*
  2: Added at the request of Marc Garbey.

  4: Inhomogeneous Laplacian in 2D. Modeled by the partial differential equation

  6:    -div \rho grad u = f,  0 < x,y < 1,

  8: with forcing function

 10:    f = e^{-x^2/\nu} e^{-y^2/\nu}

 12: with Dirichlet boundary conditions

 14:    u = f(x,y) for x = 0, x = 1, y = 0, y = 1

 16: or pure Neumman boundary conditions

 18: This uses multigrid to solve the linear system
 19: */

 21: static char help[] = "Solves 2D inhomogeneous Laplacian using multigrid.\n\n";

 23: #include <petscdm.h>
 24: #include <petscdmda.h>
 25: #include <petscksp.h>

 27: extern PetscErrorCode ComputeMatrix(KSP, Mat, Mat, void *);
 28: extern PetscErrorCode ComputeRHS(KSP, Vec, void *);

 30: typedef enum {
 31:   DIRICHLET,
 32:   NEUMANN
 33: } BCType;

 35: typedef struct {
 36:   PetscReal rho;
 37:   PetscReal nu;
 38:   BCType    bcType;
 39: } UserContext;

 41: int main(int argc, char **argv)
 42: {
 43:   KSP         ksp;
 44:   DM          da;
 45:   UserContext user;
 46:   const char *bcTypes[2] = {"dirichlet", "neumann"};
 47:   PetscInt    bc;
 48:   Vec         b, x;
 49:   PetscBool   testsolver = PETSC_FALSE;

 51:   PetscFunctionBeginUser;
 52:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
 53:   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
 54:   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 3, 3, PETSC_DECIDE, PETSC_DECIDE, 1, 1, 0, 0, &da));
 55:   PetscCall(DMSetFromOptions(da));
 56:   PetscCall(DMSetUp(da));
 57:   PetscCall(DMDASetUniformCoordinates(da, 0, 1, 0, 1, 0, 0));
 58:   PetscCall(DMDASetFieldName(da, 0, "Pressure"));

 60:   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for the inhomogeneous Poisson equation", "DMqq");
 61:   user.rho = 1.0;
 62:   PetscCall(PetscOptionsReal("-rho", "The conductivity", "ex29.c", user.rho, &user.rho, NULL));
 63:   user.nu = 0.1;
 64:   PetscCall(PetscOptionsReal("-nu", "The width of the Gaussian source", "ex29.c", user.nu, &user.nu, NULL));
 65:   bc = (PetscInt)DIRICHLET;
 66:   PetscCall(PetscOptionsEList("-bc_type", "Type of boundary condition", "ex29.c", bcTypes, 2, bcTypes[0], &bc, NULL));
 67:   user.bcType = (BCType)bc;
 68:   PetscCall(PetscOptionsBool("-testsolver", "Run solver multiple times, useful for performance studies of solver", "ex29.c", testsolver, &testsolver, NULL));
 69:   PetscOptionsEnd();

 71:   PetscCall(KSPSetComputeRHS(ksp, ComputeRHS, &user));
 72:   PetscCall(KSPSetComputeOperators(ksp, ComputeMatrix, &user));
 73:   PetscCall(KSPSetDM(ksp, da));
 74:   PetscCall(KSPSetFromOptions(ksp));
 75:   PetscCall(KSPSetUp(ksp));
 76:   PetscCall(KSPSolve(ksp, NULL, NULL));

 78:   if (testsolver) {
 79:     PetscCall(KSPGetSolution(ksp, &x));
 80:     PetscCall(KSPGetRhs(ksp, &b));
 81:     PetscCall(KSPSetDMActive(ksp, PETSC_FALSE));
 82:     PetscCall(KSPSolve(ksp, b, x));
 83:     {
 84:       PetscLogStage stage;
 85:       PetscInt      i, n = 20;

 87:       PetscCall(PetscLogStageRegister("Solve only", &stage));
 88:       PetscCall(PetscLogStagePush(stage));
 89:       for (i = 0; i < n; i++) PetscCall(KSPSolve(ksp, b, x));
 90:       PetscCall(PetscLogStagePop());
 91:     }
 92:   }

 94:   PetscCall(DMDestroy(&da));
 95:   PetscCall(KSPDestroy(&ksp));
 96:   PetscCall(PetscFinalize());
 97:   return 0;
 98: }

100: PetscErrorCode ComputeRHS(KSP ksp, Vec b, void *ctx)
101: {
102:   UserContext  *user = (UserContext *)ctx;
103:   PetscInt      i, j, mx, my, xm, ym, xs, ys;
104:   PetscScalar   Hx, Hy;
105:   PetscScalar **array;
106:   DM            da;

108:   PetscFunctionBeginUser;
109:   PetscCall(KSPGetDM(ksp, &da));
110:   PetscCall(DMDAGetInfo(da, 0, &mx, &my, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
111:   Hx = 1.0 / (PetscReal)(mx - 1);
112:   Hy = 1.0 / (PetscReal)(my - 1);
113:   PetscCall(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0));
114:   PetscCall(DMDAVecGetArray(da, b, &array));
115:   for (j = ys; j < ys + ym; j++) {
116:     for (i = xs; i < xs + xm; i++) array[j][i] = PetscExpScalar(-((PetscReal)i * Hx) * ((PetscReal)i * Hx) / user->nu) * PetscExpScalar(-((PetscReal)j * Hy) * ((PetscReal)j * Hy) / user->nu) * Hx * Hy;
117:   }
118:   PetscCall(DMDAVecRestoreArray(da, b, &array));
119:   PetscCall(VecAssemblyBegin(b));
120:   PetscCall(VecAssemblyEnd(b));

122:   /* force right-hand side to be consistent for singular matrix */
123:   /* note this is really a hack, normally the model would provide you with a consistent right handside */
124:   if (user->bcType == NEUMANN) {
125:     MatNullSpace nullspace;

127:     PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
128:     PetscCall(MatNullSpaceRemove(nullspace, b));
129:     PetscCall(MatNullSpaceDestroy(&nullspace));
130:   }
131:   PetscFunctionReturn(PETSC_SUCCESS);
132: }

134: PetscErrorCode ComputeRho(PetscInt i, PetscInt j, PetscInt mx, PetscInt my, PetscReal centerRho, PetscReal *rho)
135: {
136:   PetscFunctionBeginUser;
137:   if ((i > mx / 3.0) && (i < 2.0 * mx / 3.0) && (j > my / 3.0) && (j < 2.0 * my / 3.0)) {
138:     *rho = centerRho;
139:   } else {
140:     *rho = 1.0;
141:   }
142:   PetscFunctionReturn(PETSC_SUCCESS);
143: }

145: PetscErrorCode ComputeMatrix(KSP ksp, Mat J, Mat jac, void *ctx)
146: {
147:   UserContext *user = (UserContext *)ctx;
148:   PetscReal    centerRho;
149:   PetscInt     i, j, mx, my, xm, ym, xs, ys;
150:   PetscScalar  v[5];
151:   PetscReal    Hx, Hy, HydHx, HxdHy, rho;
152:   MatStencil   row, col[5];
153:   DM           da;
154:   PetscBool    check_matis = PETSC_FALSE;

156:   PetscFunctionBeginUser;
157:   PetscCall(KSPGetDM(ksp, &da));
158:   centerRho = user->rho;
159:   PetscCall(DMDAGetInfo(da, 0, &mx, &my, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
160:   Hx    = 1.0 / (PetscReal)(mx - 1);
161:   Hy    = 1.0 / (PetscReal)(my - 1);
162:   HxdHy = Hx / Hy;
163:   HydHx = Hy / Hx;
164:   PetscCall(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0));
165:   for (j = ys; j < ys + ym; j++) {
166:     for (i = xs; i < xs + xm; i++) {
167:       row.i = i;
168:       row.j = j;
169:       PetscCall(ComputeRho(i, j, mx, my, centerRho, &rho));
170:       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
171:         if (user->bcType == DIRICHLET) {
172:           v[0] = 2.0 * rho * (HxdHy + HydHx);
173:           PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES));
174:         } else if (user->bcType == NEUMANN) {
175:           PetscInt numx = 0, numy = 0, num = 0;
176:           if (j != 0) {
177:             v[num]     = -rho * HxdHy;
178:             col[num].i = i;
179:             col[num].j = j - 1;
180:             numy++;
181:             num++;
182:           }
183:           if (i != 0) {
184:             v[num]     = -rho * HydHx;
185:             col[num].i = i - 1;
186:             col[num].j = j;
187:             numx++;
188:             num++;
189:           }
190:           if (i != mx - 1) {
191:             v[num]     = -rho * HydHx;
192:             col[num].i = i + 1;
193:             col[num].j = j;
194:             numx++;
195:             num++;
196:           }
197:           if (j != my - 1) {
198:             v[num]     = -rho * HxdHy;
199:             col[num].i = i;
200:             col[num].j = j + 1;
201:             numy++;
202:             num++;
203:           }
204:           v[num]     = numx * rho * HydHx + numy * rho * HxdHy;
205:           col[num].i = i;
206:           col[num].j = j;
207:           num++;
208:           PetscCall(MatSetValuesStencil(jac, 1, &row, num, col, v, INSERT_VALUES));
209:         }
210:       } else {
211:         v[0]     = -rho * HxdHy;
212:         col[0].i = i;
213:         col[0].j = j - 1;
214:         v[1]     = -rho * HydHx;
215:         col[1].i = i - 1;
216:         col[1].j = j;
217:         v[2]     = 2.0 * rho * (HxdHy + HydHx);
218:         col[2].i = i;
219:         col[2].j = j;
220:         v[3]     = -rho * HydHx;
221:         col[3].i = i + 1;
222:         col[3].j = j;
223:         v[4]     = -rho * HxdHy;
224:         col[4].i = i;
225:         col[4].j = j + 1;
226:         PetscCall(MatSetValuesStencil(jac, 1, &row, 5, col, v, INSERT_VALUES));
227:       }
228:     }
229:   }
230:   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
231:   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
232:   PetscCall(MatViewFromOptions(jac, NULL, "-view_mat"));
233:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-check_matis", &check_matis, NULL));
234:   if (check_matis) {
235:     void (*f)(void);
236:     Mat       J2;
237:     MatType   jtype;
238:     PetscReal nrm;

240:     PetscCall(MatGetType(jac, &jtype));
241:     PetscCall(MatConvert(jac, MATIS, MAT_INITIAL_MATRIX, &J2));
242:     PetscCall(MatViewFromOptions(J2, NULL, "-view_conv"));
243:     PetscCall(MatConvert(J2, jtype, MAT_INPLACE_MATRIX, &J2));
244:     PetscCall(MatGetOperation(jac, MATOP_VIEW, &f));
245:     PetscCall(MatSetOperation(J2, MATOP_VIEW, f));
246:     PetscCall(MatSetDM(J2, da));
247:     PetscCall(MatViewFromOptions(J2, NULL, "-view_conv_assembled"));
248:     PetscCall(MatAXPY(J2, -1., jac, DIFFERENT_NONZERO_PATTERN));
249:     PetscCall(MatNorm(J2, NORM_FROBENIUS, &nrm));
250:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MATIS %g\n", (double)nrm));
251:     PetscCall(MatViewFromOptions(J2, NULL, "-view_conv_err"));
252:     PetscCall(MatDestroy(&J2));
253:   }
254:   if (user->bcType == NEUMANN) {
255:     MatNullSpace nullspace;

257:     PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
258:     PetscCall(MatSetNullSpace(J, nullspace));
259:     PetscCall(MatNullSpaceDestroy(&nullspace));
260:   }
261:   PetscFunctionReturn(PETSC_SUCCESS);
262: }

264: /*TEST

266:    test:
267:       args: -pc_type mg -pc_mg_type full -ksp_type fgmres -ksp_monitor_short -da_refine 8 -ksp_rtol 1.e-3

269:    test:
270:       suffix: 2
271:       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
272:       requires: !single

274:    test:
275:       suffix: telescope
276:       nsize: 4
277:       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

279:    test:
280:       suffix: 3
281:       args: -ksp_view -da_refine 2 -pc_type mg -pc_mg_distinct_smoothup -mg_levels_up_pc_type jacobi

283:    test:
284:       suffix: 4
285:       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

287:    testset:
288:      suffix: aniso
289:      args: -da_grid_x 10 -da_grid_y 2 -da_refine 2 -pc_type mg -ksp_monitor_short -mg_levels_ksp_max_it 6 -mg_levels_pc_type jacobi
290:      test:
291:        suffix: first
292:        args: -mg_levels_ksp_chebyshev_kind first
293:      test:
294:        suffix: fourth
295:        args: -mg_levels_ksp_chebyshev_kind fourth
296:      test:
297:        suffix: opt_fourth
298:        args: -mg_levels_ksp_chebyshev_kind opt_fourth

300:    test:
301:       suffix: 5
302:       nsize: 2
303:       requires: hypre !complex
304:       args: -pc_type mg -da_refine 2 -ksp_monitor -matptap_via hypre -pc_mg_galerkin both

306:    test:
307:       suffix: 6
308:       args: -pc_type svd -pc_svd_monitor ::all

310: TEST*/