Actual source code: ex45.c

  1: /*
  2: Laplacian in 3D. Modeled by the partial differential equation

  4:    - Laplacian u = 1,0 < x,y,z < 1,

  6: with boundary conditions

  8:    u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.

 10:    This uses multigrid to solve the linear system

 12:    See src/snes/tutorials/ex50.c

 14:    Can also be run with -pc_type exotic -ksp_type fgmres

 16: */

 18: static char help[] = "Solves 3D Laplacian using multigrid.\n\n";

 20: #include <petscksp.h>
 21: #include <petscdm.h>
 22: #include <petscdmda.h>

 24: extern PetscErrorCode ComputeMatrix(KSP, Mat, Mat, void *);
 25: extern PetscErrorCode ComputeRHS(KSP, Vec, void *);
 26: extern PetscErrorCode ComputeInitialGuess(KSP, Vec, void *);

 28: int main(int argc, char **argv)
 29: {
 30:   KSP       ksp;
 31:   PetscReal norm;
 32:   DM        da;
 33:   Vec       x, b, r;
 34:   Mat       A;

 36:   PetscFunctionBeginUser;
 37:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));

 39:   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
 40:   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 7, 7, 7, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, 0, 0, 0, &da));
 41:   PetscCall(DMSetFromOptions(da));
 42:   PetscCall(DMSetUp(da));
 43:   PetscCall(KSPSetDM(ksp, da));
 44:   PetscCall(KSPSetComputeInitialGuess(ksp, ComputeInitialGuess, NULL));
 45:   PetscCall(KSPSetComputeRHS(ksp, ComputeRHS, NULL));
 46:   PetscCall(KSPSetComputeOperators(ksp, ComputeMatrix, NULL));
 47:   PetscCall(DMDestroy(&da));

 49:   PetscCall(KSPSetFromOptions(ksp));
 50:   PetscCall(KSPSolve(ksp, NULL, NULL));
 51:   PetscCall(KSPGetSolution(ksp, &x));
 52:   PetscCall(KSPGetRhs(ksp, &b));
 53:   PetscCall(VecDuplicate(b, &r));
 54:   PetscCall(KSPGetOperators(ksp, &A, NULL));

 56:   PetscCall(MatMult(A, x, r));
 57:   PetscCall(VecAXPY(r, -1.0, b));
 58:   PetscCall(VecNorm(r, NORM_2, &norm));
 59:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Residual norm %g\n", (double)norm));

 61:   PetscCall(VecDestroy(&r));
 62:   PetscCall(KSPDestroy(&ksp));
 63:   PetscCall(PetscFinalize());
 64:   return 0;
 65: }

 67: PetscErrorCode ComputeRHS(KSP ksp, Vec b, void *ctx)
 68: {
 69:   PetscInt       i, j, k, mx, my, mz, xm, ym, zm, xs, ys, zs;
 70:   DM             dm;
 71:   PetscScalar    Hx, Hy, Hz, HxHydHz, HyHzdHx, HxHzdHy;
 72:   PetscScalar ***barray;

 74:   PetscFunctionBeginUser;
 75:   PetscCall(KSPGetDM(ksp, &dm));
 76:   PetscCall(DMDAGetInfo(dm, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
 77:   Hx      = 1.0 / (PetscReal)(mx - 1);
 78:   Hy      = 1.0 / (PetscReal)(my - 1);
 79:   Hz      = 1.0 / (PetscReal)(mz - 1);
 80:   HxHydHz = Hx * Hy / Hz;
 81:   HxHzdHy = Hx * Hz / Hy;
 82:   HyHzdHx = Hy * Hz / Hx;
 83:   PetscCall(DMDAGetCorners(dm, &xs, &ys, &zs, &xm, &ym, &zm));
 84:   PetscCall(DMDAVecGetArray(dm, b, &barray));

 86:   for (k = zs; k < zs + zm; k++) {
 87:     for (j = ys; j < ys + ym; j++) {
 88:       for (i = xs; i < xs + xm; i++) {
 89:         if (i == 0 || j == 0 || k == 0 || i == mx - 1 || j == my - 1 || k == mz - 1) {
 90:           barray[k][j][i] = 2.0 * (HxHydHz + HxHzdHy + HyHzdHx);
 91:         } else {
 92:           barray[k][j][i] = Hx * Hy * Hz;
 93:         }
 94:       }
 95:     }
 96:   }
 97:   PetscCall(DMDAVecRestoreArray(dm, b, &barray));
 98:   PetscFunctionReturn(PETSC_SUCCESS);
 99: }

101: PetscErrorCode ComputeInitialGuess(KSP ksp, Vec b, void *ctx)
102: {
103:   PetscFunctionBeginUser;
104:   PetscCall(VecSet(b, 0));
105:   PetscFunctionReturn(PETSC_SUCCESS);
106: }

108: PetscErrorCode ComputeMatrix(KSP ksp, Mat jac, Mat B, void *ctx)
109: {
110:   DM          da;
111:   PetscInt    i, j, k, mx, my, mz, xm, ym, zm, xs, ys, zs;
112:   PetscScalar v[7], Hx, Hy, Hz, HxHydHz, HyHzdHx, HxHzdHy;
113:   MatStencil  row, col[7];

115:   PetscFunctionBeginUser;
116:   PetscCall(KSPGetDM(ksp, &da));
117:   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
118:   Hx      = 1.0 / (PetscReal)(mx - 1);
119:   Hy      = 1.0 / (PetscReal)(my - 1);
120:   Hz      = 1.0 / (PetscReal)(mz - 1);
121:   HxHydHz = Hx * Hy / Hz;
122:   HxHzdHy = Hx * Hz / Hy;
123:   HyHzdHx = Hy * Hz / Hx;
124:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));

126:   for (k = zs; k < zs + zm; k++) {
127:     for (j = ys; j < ys + ym; j++) {
128:       for (i = xs; i < xs + xm; i++) {
129:         row.i = i;
130:         row.j = j;
131:         row.k = k;
132:         if (i == 0 || j == 0 || k == 0 || i == mx - 1 || j == my - 1 || k == mz - 1) {
133:           v[0] = 2.0 * (HxHydHz + HxHzdHy + HyHzdHx);
134:           PetscCall(MatSetValuesStencil(B, 1, &row, 1, &row, v, INSERT_VALUES));
135:         } else {
136:           v[0]     = -HxHydHz;
137:           col[0].i = i;
138:           col[0].j = j;
139:           col[0].k = k - 1;
140:           v[1]     = -HxHzdHy;
141:           col[1].i = i;
142:           col[1].j = j - 1;
143:           col[1].k = k;
144:           v[2]     = -HyHzdHx;
145:           col[2].i = i - 1;
146:           col[2].j = j;
147:           col[2].k = k;
148:           v[3]     = 2.0 * (HxHydHz + HxHzdHy + HyHzdHx);
149:           col[3].i = row.i;
150:           col[3].j = row.j;
151:           col[3].k = row.k;
152:           v[4]     = -HyHzdHx;
153:           col[4].i = i + 1;
154:           col[4].j = j;
155:           col[4].k = k;
156:           v[5]     = -HxHzdHy;
157:           col[5].i = i;
158:           col[5].j = j + 1;
159:           col[5].k = k;
160:           v[6]     = -HxHydHz;
161:           col[6].i = i;
162:           col[6].j = j;
163:           col[6].k = k + 1;
164:           PetscCall(MatSetValuesStencil(B, 1, &row, 7, col, v, INSERT_VALUES));
165:         }
166:       }
167:     }
168:   }
169:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
170:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
171:   PetscFunctionReturn(PETSC_SUCCESS);
172: }

174: /*TEST

176:    test:
177:       nsize: 4
178:       args: -pc_type exotic -ksp_monitor_short -ksp_type fgmres -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 1 -mg_levels_pc_type bjacobi
179:       output_file: output/ex45_1.out

181:    test:
182:       suffix: 2
183:       nsize: 4
184:       args: -ksp_monitor_short -da_grid_x 21 -da_grid_y 21 -da_grid_z 21 -pc_type mg -pc_mg_levels 3 -mg_levels_ksp_type richardson -mg_levels_ksp_max_it 1 -mg_levels_pc_type bjacobi

186:    test:
187:       suffix: telescope
188:       nsize: 4
189:       args: -ksp_type fgmres -ksp_monitor_short -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type jacobi -pc_mg_levels 2 -da_grid_x 65 -da_grid_y 65 -da_grid_z 65 -mg_coarse_pc_type telescope -mg_coarse_pc_telescope_ignore_kspcomputeoperators -mg_coarse_pc_telescope_reduction_factor 4 -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 richardson -mg_coarse_telescope_mg_levels_pc_type jacobi -mg_levels_ksp_type richardson -mg_coarse_telescope_mg_levels_ksp_type richardson -ksp_rtol 1.0e-4

191:    test:
192:       suffix: telescope_2
193:       nsize: 4
194:       args: -ksp_type fgmres -ksp_monitor_short -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type jacobi -pc_mg_levels 2 -da_grid_x 65 -da_grid_y 65 -da_grid_z 65 -mg_coarse_pc_type telescope -mg_coarse_pc_telescope_reduction_factor 2 -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 richardson -mg_coarse_telescope_mg_levels_pc_type jacobi -mg_levels_ksp_type richardson -mg_coarse_telescope_mg_levels_ksp_type richardson -ksp_rtol 1.0e-4

196: TEST*/