Actual source code: ex42.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 2 level multigrid with the wirebasket based coarse problem to solve the linear system
13: */
15: static char help[] = "Solves 3D Laplacian using wirebasket based multigrid.\n\n";
17: #include petscda.h
18: #include petscksp.h
19: #include petscmg.h
27: int main(int argc,char **argv)
28: {
30: KSP ksp;
31: PC pc;
32: Vec x,b;
33: DA da;
34: Mat A;
36: PetscInitialize(&argc,&argv,(char *)0,help);
38: DACreate3d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,-8,-8,-8,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da);
39: DACreateGlobalVector(da,&x);
40: DACreateGlobalVector(da,&b);
41: ComputeRHS(da,b);
42: DAGetMatrix(da,MATAIJ,&A);
43: ComputeMatrix(da,A);
45: KSPCreate(PETSC_COMM_WORLD,&ksp);
46: KSPSetFromOptions(ksp);
47: KSPSetOperators(ksp,A,A,SAME_NONZERO_PATTERN);
48: KSPGetPC(ksp,&pc);
49: PCExoticSetDA(pc,da);
50:
51: KSPSolve(ksp,b,x);
52:
53: KSPDestroy(ksp);
54: VecDestroy(x);
55: VecDestroy(b);
56: MatDestroy(A);
57: DADestroy(da);
58: PetscFinalize();
60: return 0;
61: }
65: PetscErrorCode ComputeRHS(DA da,Vec b)
66: {
68: PetscInt mx,my,mz;
69: PetscScalar h;
72: DAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0);
73: h = 1.0/((mx-1)*(my-1)*(mz-1));
74: VecSet(b,h);
75: return(0);
76: }
77:
80: PetscErrorCode ComputeMatrix(DA da,Mat B)
81: {
83: PetscInt i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs;
84: PetscScalar v[7],Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
85: MatStencil row,col[7];
87: DAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0);
88: Hx = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
89: HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;
90: DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
91:
92: for (k=zs; k<zs+zm; k++){
93: for (j=ys; j<ys+ym; j++){
94: for(i=xs; i<xs+xm; i++){
95: row.i = i; row.j = j; row.k = k;
96: if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1){
97: v[0] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
98: MatSetValuesStencil(B,1,&row,1,&row,v,INSERT_VALUES);
99: } else {
100: v[0] = -HxHydHz;col[0].i = i; col[0].j = j; col[0].k = k-1;
101: v[1] = -HxHzdHy;col[1].i = i; col[1].j = j-1; col[1].k = k;
102: v[2] = -HyHzdHx;col[2].i = i-1; col[2].j = j; col[2].k = k;
103: v[3] = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);col[3].i = row.i; col[3].j = row.j; col[3].k = row.k;
104: v[4] = -HyHzdHx;col[4].i = i+1; col[4].j = j; col[4].k = k;
105: v[5] = -HxHzdHy;col[5].i = i; col[5].j = j+1; col[5].k = k;
106: v[6] = -HxHydHz;col[6].i = i; col[6].j = j; col[6].k = k+1;
107: MatSetValuesStencil(B,1,&row,7,col,v,INSERT_VALUES);
108: }
109: }
110: }
111: }
112: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
113: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
114: return 0;
115: }
121: PetscErrorCode ComputeInterpolation(PC pc,void *ida)
122: {
124: DA da = (DA)ida;
125: Mat A,P;
128: PCGetOperators(pc,PETSC_NULL,&A,PETSC_NULL);
129: DAGetWireBasketInterpolation(da,A,MAT_INITIAL_MATRIX,&P);
130: PCMGSetInterpolation(pc,1,P);
131: MatDestroy(P);
133: return(0);
134: }