Actual source code: ex45.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  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*);

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

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

 42:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 43:   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);
 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();

 66:   return 0;
 67: }

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

 80:   KSPGetDM(ksp,&dm);
 81:   DMDAGetInfo(dm,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);
 82:   Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
 83:   HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; 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: }

104: PetscErrorCode ComputeInitialGuess(KSP ksp,Vec b,void *ctx)
105: {

109:   VecSet(b,0);
110:   return(0);
111: }

115: PetscErrorCode ComputeMatrix(KSP ksp,Mat jac,Mat B,void *ctx)
116: {
117:   DM             da;
119:   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
120:   PetscScalar    v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
121:   MatStencil     row,col[7];

124:   KSPGetDM(ksp,&da);
125:   DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);
126:   Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
127:   HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
128:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);

130:   for (k=zs; k<zs+zm; k++) {
131:     for (j=ys; j<ys+ym; j++) {
132:       for (i=xs; i<xs+xm; i++) {
133:         row.i = i; row.j = j; row.k = k;
134:         if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
135:           v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
136:           MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);
137:         } else {
138:           v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
139:           v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
140:           v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
141:           v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
142:           v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
143:           v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
144:           v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
145:           MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);
146:         }
147:       }
148:     }
149:   }
150:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
151:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
152:   return(0);
153: }