Actual source code: ex22.c
2: /*$Id: ex22.c,v 1.19 2001/08/07 21:30:54 bsmith Exp $*/
3: /*
4: Laplacian in 3D. Modeled by the partial differential equation
6: Laplacian u = 1,0 < x,y,z < 1,
8: with boundary conditions
10: u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
12: This uses multigrid to solve the linear system
14: */
16: static char help[] = "Solves 3D Laplacian using multigrid.nn";
18: #include petscda.h
19: #include petscsles.h
21: extern int ComputeJacobian(DMMG,Mat);
22: extern int ComputeRHS(DMMG,Vec);
24: #undef __FUNCT__
26: int main(int argc,char **argv)
27: {
28: int ierr;
29: DMMG *dmmg;
30: PetscScalar mone = -1.0;
31: PetscReal norm;
32: DA da;
34: PetscInitialize(&argc,&argv,(char *)0,help);
36: DMMGCreate(PETSC_COMM_WORLD,3,PETSC_NULL,&dmmg);
37: DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,3,3,3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);
38: DMMGSetDM(dmmg,(DM)da);
39: DADestroy(da);
41: DMMGSetSLES(dmmg,ComputeRHS,ComputeJacobian);
43: DMMGSolve(dmmg);
45: MatMult(DMMGGetJ(dmmg),DMMGGetx(dmmg),DMMGGetr(dmmg));
46: VecAXPY(&mone,DMMGGetb(dmmg),DMMGGetr(dmmg));
47: VecNorm(DMMGGetr(dmmg),NORM_2,&norm);
48: /* PetscPrintf(PETSC_COMM_WORLD,"Residual norm %gn",norm); */
50: DMMGDestroy(dmmg);
51: PetscFinalize();
53: return 0;
54: }
56: #undef __FUNCT__
58: int ComputeRHS(DMMG dmmg,Vec b)
59: {
60: int ierr,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(&h,b);
67: return(0);
68: }
69:
70: #undef __FUNCT__
72: int ComputeJacobian(DMMG dmmg,Mat jac)
73: {
74: DA da = (DA)dmmg->dm;
75: int ierr,i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
76: PetscScalar v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
77: MatStencil row,col[7];
79: DAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0);
80: Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
81: HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
82: DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
83:
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: row.i = i; row.j = j; row.k = k;
88: if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1){
89: v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
90: MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
91: } else {
92: v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
93: v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
94: v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
95: v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
96: v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
97: v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
98: v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
99: MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);
100: }
101: }
102: }
103: }
104: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
105: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
106: return 0;
107: }