Actual source code: ex45.c
petsc-3.8.4 2018-03-24
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/examples/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: {
32: KSP ksp;
33: PetscReal norm;
34: DM da;
35: Vec x,b,r;
36: Mat A;
38: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
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 ierr;
66: }
68: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
69: {
71: PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
72: DM dm;
73: PetscScalar Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
74: PetscScalar ***barray;
77: KSPGetDM(ksp,&dm);
78: DMDAGetInfo(dm,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);
79: Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
80: HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
81: DMDAGetCorners(dm,&xs,&ys,&zs,&xm,&ym,&zm);
82: DMDAVecGetArray(dm,b,&barray);
84: for (k=zs; k<zs+zm; k++) {
85: for (j=ys; j<ys+ym; j++) {
86: for (i=xs; i<xs+xm; i++) {
87: if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
88: barray[k][j][i] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
89: } else {
90: barray[k][j][i] = Hx*Hy*Hz;
91: }
92: }
93: }
94: }
95: DMDAVecRestoreArray(dm,b,&barray);
96: return(0);
97: }
99: PetscErrorCode ComputeInitialGuess(KSP ksp,Vec b,void *ctx)
100: {
104: VecSet(b,0);
105: return(0);
106: }
108: PetscErrorCode ComputeMatrix(KSP ksp,Mat jac,Mat B,void *ctx)
109: {
110: 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); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
120: HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
121: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
123: for (k=zs; k<zs+zm; k++) {
124: for (j=ys; j<ys+ym; j++) {
125: for (i=xs; i<xs+xm; i++) {
126: row.i = i; row.j = j; row.k = k;
127: if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
128: v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
129: MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);
130: } else {
131: v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
132: v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
133: v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
134: v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
135: v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
136: v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
137: v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
138: MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);
139: }
140: }
141: }
142: }
143: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
144: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
145: return(0);
146: }