Actual source code: ex45.c


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

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

  7: with boundary conditions

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

 11:    This uses multigrid to solve the linear system

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

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

 17: */

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

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

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

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

 38:   PetscInitialize(&argc, &argv, (char *)0, help);

 40:   KSPCreate(PETSC_COMM_WORLD, &ksp);
 41:   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);
 42:   DMSetFromOptions(da);
 43:   DMSetUp(da);
 44:   KSPSetDM(ksp, da);
 45:   KSPSetComputeInitialGuess(ksp, ComputeInitialGuess, NULL);
 46:   KSPSetComputeRHS(ksp, ComputeRHS, NULL);
 47:   KSPSetComputeOperators(ksp, ComputeMatrix, NULL);
 48:   DMDestroy(&da);

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

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

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

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

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

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

102: PetscErrorCode ComputeInitialGuess(KSP ksp, Vec b, void *ctx)
103: {
105:   VecSet(b, 0);
106:   return 0;
107: }

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

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

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

175: /*TEST

177:    test:
178:       nsize: 4
179:       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
180:       output_file: output/ex45_1.out

182:    test:
183:       suffix: 2
184:       nsize: 4
185:       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

187:    test:
188:       suffix: telescope
189:       nsize: 4
190:       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

192:    test:
193:       suffix: telescope_2
194:       nsize: 4
195:       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

197: TEST*/