Actual source code: ex34.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: /*T
  2:    Concepts: KSP^solving a system of linear equations
  3:    Concepts: KSP^Laplacian, 3d
  4:    Processors: n
  5: T*/

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

 10:    div  grad u = f,  0 < x,y,z < 1,

 12: with pure Neumann boundary conditions

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

 16: The functions are cell-centered

 18: This uses multigrid to solve the linear system

 20:        Contributed by Jianming Yang <jianming-yang@uiowa.edu>
 21: */

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

 25: #include <petscdm.h>
 26: #include <petscdmda.h>
 27: #include <petscksp.h>

 29: extern PetscErrorCode ComputeMatrix(KSP,Mat,Mat,void*);
 30: extern PetscErrorCode ComputeRHS(KSP,Vec,void*);

 34: int main(int argc,char **argv)
 35: {
 36:   KSP            ksp;
 37:   DM             da;
 38:   PetscReal      norm;

 41:   PetscInt    i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
 42:   PetscScalar Hx,Hy,Hz;
 43:   PetscScalar ***array;
 44:   Vec         x,b,r;
 45:   Mat         J;

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

 49:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 50:   DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-12,-12,-12,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);
 51:   DMDASetInterpolationType(da, DMDA_Q0);

 53:   KSPSetDM(ksp,da);

 55:   KSPSetComputeRHS(ksp,ComputeRHS,NULL);
 56:   KSPSetComputeOperators(ksp,ComputeMatrix,NULL);
 57:   KSPSetFromOptions(ksp);
 58:   KSPSolve(ksp,NULL,NULL);
 59:   KSPGetSolution(ksp,&x);
 60:   KSPGetRhs(ksp,&b);
 61:   KSPGetOperators(ksp,NULL,&J);
 62:   VecDuplicate(b,&r);

 64:   MatMult(J,x,r);
 65:   VecAXPY(r,-1.0,b);
 66:   VecNorm(r,NORM_2,&norm);
 67:   PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)norm);

 69:   DMDAGetInfo(da, 0, &mx, &my, &mz, 0,0,0,0,0,0,0,0,0);
 70:   Hx   = 1.0 / (PetscReal)(mx);
 71:   Hy   = 1.0 / (PetscReal)(my);
 72:   Hz   = 1.0 / (PetscReal)(mz);
 73:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
 74:   DMDAVecGetArray(da, x, &array);

 76:   for (k=zs; k<zs+zm; k++) {
 77:     for (j=ys; j<ys+ym; j++) {
 78:       for (i=xs; i<xs+xm; i++) {
 79:         array[k][j][i] -=
 80:           PetscCosScalar(2*PETSC_PI*(((PetscReal)i+0.5)*Hx))*
 81:           PetscCosScalar(2*PETSC_PI*(((PetscReal)j+0.5)*Hy))*
 82:           PetscCosScalar(2*PETSC_PI*(((PetscReal)k+0.5)*Hz));
 83:       }
 84:     }
 85:   }
 86:   DMDAVecRestoreArray(da, x, &array);
 87:   VecAssemblyBegin(x);
 88:   VecAssemblyEnd(x);

 90:   VecNorm(x,NORM_INFINITY,&norm);
 91:   PetscPrintf(PETSC_COMM_WORLD,"Error norm %g\n",(double)norm);
 92:   VecNorm(x,NORM_1,&norm);
 93:   PetscPrintf(PETSC_COMM_WORLD,"Error norm %g\n",(double)(norm/((PetscReal)(mx)*(PetscReal)(my)*(PetscReal)(mz))));
 94:   VecNorm(x,NORM_2,&norm);
 95:   PetscPrintf(PETSC_COMM_WORLD,"Error norm %g\n",(double)(norm/((PetscReal)(mx)*(PetscReal)(my)*(PetscReal)(mz))));

 97:   VecDestroy(&r);
 98:   KSPDestroy(&ksp);
 99:   DMDestroy(&da);
100:   PetscFinalize();
101:   return 0;
102: }

106: PetscErrorCode ComputeRHS(KSP ksp,Vec b,void *ctx)
107: {
109:   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
110:   PetscScalar    Hx,Hy,Hz;
111:   PetscScalar    ***array;
112:   DM             da;
113:   MatNullSpace   nullspace;

116:   KSPGetDM(ksp,&da);
117:   DMDAGetInfo(da, 0, &mx, &my, &mz, 0,0,0,0,0,0,0,0,0);
118:   Hx   = 1.0 / (PetscReal)(mx);
119:   Hy   = 1.0 / (PetscReal)(my);
120:   Hz   = 1.0 / (PetscReal)(mz);
121:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
122:   DMDAVecGetArray(da, b, &array);
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:         array[k][j][i] = 12 * PETSC_PI * PETSC_PI
127:                          * PetscCosScalar(2*PETSC_PI*(((PetscReal)i+0.5)*Hx))
128:                          * PetscCosScalar(2*PETSC_PI*(((PetscReal)j+0.5)*Hy))
129:                          * PetscCosScalar(2*PETSC_PI*(((PetscReal)k+0.5)*Hz))
130:                          * Hx * Hy * Hz;
131:       }
132:     }
133:   }
134:   DMDAVecRestoreArray(da, b, &array);
135:   VecAssemblyBegin(b);
136:   VecAssemblyEnd(b);

138:   /* force right hand side to be consistent for singular matrix */
139:   /* note this is really a hack, normally the model would provide you with a consistent right handside */

141:   MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
142:   MatNullSpaceRemove(nullspace,b);
143:   MatNullSpaceDestroy(&nullspace);
144:   return(0);
145: }


150: PetscErrorCode ComputeMatrix(KSP ksp, Mat J,Mat jac, void *ctx)
151: {
153:   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs,num, numi, numj, numk;
154:   PetscScalar    v[7],Hx,Hy,Hz,HyHzdHx,HxHzdHy,HxHydHz;
155:   MatStencil     row, col[7];
156:   DM             da;
157:   MatNullSpace   nullspace;

160:   KSPGetDM(ksp,&da);
161:   DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);
162:   Hx      = 1.0 / (PetscReal)(mx);
163:   Hy      = 1.0 / (PetscReal)(my);
164:   Hz      = 1.0 / (PetscReal)(mz);
165:   HyHzdHx = Hy*Hz/Hx;
166:   HxHzdHy = Hx*Hz/Hy;
167:   HxHydHz = Hx*Hy/Hz;
168:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
169:   for (k=zs; k<zs+zm; k++) {
170:     for (j=ys; j<ys+ym; j++) {
171:       for (i=xs; i<xs+xm; i++) {
172:         row.i = i; row.j = j; row.k = k;
173:         if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) {
174:           num = 0; numi=0; numj=0; numk=0;
175:           if (k!=0) {
176:             v[num]     = -HxHydHz;
177:             col[num].i = i;
178:             col[num].j = j;
179:             col[num].k = k-1;
180:             num++; numk++;
181:           }
182:           if (j!=0) {
183:             v[num]     = -HxHzdHy;
184:             col[num].i = i;
185:             col[num].j = j-1;
186:             col[num].k = k;
187:             num++; numj++;
188:             }
189:           if (i!=0) {
190:             v[num]     = -HyHzdHx;
191:             col[num].i = i-1;
192:             col[num].j = j;
193:             col[num].k = k;
194:             num++; numi++;
195:           }
196:           if (i!=mx-1) {
197:             v[num]     = -HyHzdHx;
198:             col[num].i = i+1;
199:             col[num].j = j;
200:             col[num].k = k;
201:             num++; numi++;
202:           }
203:           if (j!=my-1) {
204:             v[num]     = -HxHzdHy;
205:             col[num].i = i;
206:             col[num].j = j+1;
207:             col[num].k = k;
208:             num++; numj++;
209:           }
210:           if (k!=mz-1) {
211:             v[num]     = -HxHydHz;
212:             col[num].i = i;
213:             col[num].j = j;
214:             col[num].k = k+1;
215:             num++; numk++;
216:           }
217:           v[num]     = (PetscReal)(numk)*HxHydHz + (PetscReal)(numj)*HxHzdHy + (PetscReal)(numi)*HyHzdHx;
218:           col[num].i = i;   col[num].j = j;   col[num].k = k;
219:           num++;
220:           MatSetValuesStencil(jac,1,&row,num,col,v,INSERT_VALUES);
221:         } else {
222:           v[0] = -HxHydHz;                          col[0].i = i;   col[0].j = j;   col[0].k = k-1;
223:           v[1] = -HxHzdHy;                          col[1].i = i;   col[1].j = j-1; col[1].k = k;
224:           v[2] = -HyHzdHx;                          col[2].i = i-1; col[2].j = j;   col[2].k = k;
225:           v[3] = 2.0*(HyHzdHx + HxHzdHy + HxHydHz); col[3].i = i;   col[3].j = j;   col[3].k = k;
226:           v[4] = -HyHzdHx;                          col[4].i = i+1; col[4].j = j;   col[4].k = k;
227:           v[5] = -HxHzdHy;                          col[5].i = i;   col[5].j = j+1; col[5].k = k;
228:           v[6] = -HxHydHz;                          col[6].i = i;   col[6].j = j;   col[6].k = k+1;
229:           MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);
230:         }
231:       }
232:     }
233:   }
234:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
235:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
236:   MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,0,&nullspace);
237:   MatSetNullSpace(J,nullspace);
238:   MatNullSpaceDestroy(&nullspace);
239:   return(0);
240: }