Actual source code: ex22.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: */
15: static char help[] = "Solves 3D Laplacian using multigrid.\n\n";
17: #include petscda.h
18: #include petscksp.h
19: #include petscdmmg.h
26: int main(int argc,char **argv)
27: {
29: DMMG *dmmg;
30: PetscReal norm;
31: DA da;
33: PetscInitialize(&argc,&argv,(char *)0,help);
35: DMMGCreate(PETSC_COMM_WORLD,3,PETSC_NULL,&dmmg);
36: DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-3,-3,-3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);
37: DMMGSetDM(dmmg,(DM)da);
38: DADestroy(da);
40: DMMGSetKSP(dmmg,ComputeRHS,ComputeJacobian);
42: DMMGSolve(dmmg);
44: MatMult(DMMGGetJ(dmmg),DMMGGetx(dmmg),DMMGGetr(dmmg));
45: VecAXPY(DMMGGetr(dmmg),-1.0,DMMGGetRHS(dmmg));
46: VecNorm(DMMGGetr(dmmg),NORM_2,&norm);
47: /* PetscPrintf(PETSC_COMM_WORLD,"Residual norm %G\n",norm); */
49: DMMGDestroy(dmmg);
50: PetscFinalize();
52: return 0;
53: }
57: PetscErrorCode ComputeRHS(DMMG dmmg,Vec b)
58: {
60: PetscInt mx,my,mz;
61: PetscScalar h;
64: DAGetInfo((DA)dmmg->dm,0,&mx,&my,&mz,0,0,0,0,0,0,0);
65: h = 1.0/((mx-1)*(my-1)*(mz-1));
66: VecSet(b,h);
67: return(0);
68: }
69:
72: PetscErrorCode ComputeJacobian(DMMG dmmg,Mat jac,Mat B)
73: {
74: DA da = (DA)dmmg->dm;
76: PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
77: PetscScalar v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
78: MatStencil row,col[7];
80: DAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0);
81: Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
82: HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
83: DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
84:
85: for (k=zs; k<zs+zm; k++){
86: for (j=ys; j<ys+ym; j++){
87: for(i=xs; i<xs+xm; i++){
88: row.i = i; row.j = j; row.k = k;
89: if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1){
90: v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
91: MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);
92: } else {
93: v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
94: v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
95: v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
96: v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
97: v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
98: v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
99: v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
100: MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);
101: }
102: }
103: }
104: }
105: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
106: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
107: return 0;
108: }