Actual source code: ex32.c

petsc-3.4.2 2013-07-02
  1: /*
  2:   Laplacian in 3D. Use for testing BAIJ matrix.
  3:   Modeled by the partial differential equation

  5:    - Laplacian u = 1,0 < x,y,z < 1,

  7:    with boundary conditions
  8:    u = 1 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1.
  9: */

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

 13: #include <petscdmda.h>
 14: #include <petscksp.h>

 16: extern PetscErrorCode ComputeMatrix(DM,Mat);
 17: extern PetscErrorCode ComputeRHS(DM,Vec);

 21: int main(int argc,char **argv)
 22: {
 24:   KSP            ksp;
 25:   PC             pc;
 26:   Vec            x,b;
 27:   DM             da;
 28:   Mat            A,Atrans;
 29:   PetscInt       dof=1,M=-8;
 30:   PetscBool      flg,trans=PETSC_FALSE;

 32:   PetscInitialize(&argc,&argv,(char*)0,help);
 33:   PetscOptionsGetInt(NULL,"-dof",&dof,NULL);
 34:   PetscOptionsGetInt(NULL,"-M",&M,NULL);
 35:   PetscOptionsGetBool(NULL,"-trans",&trans,NULL);

 37:   DMDACreate(PETSC_COMM_WORLD,&da);
 38:   DMDASetDim(da,3);
 39:   DMDASetBoundaryType(da,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE);
 40:   DMDASetStencilType(da,DMDA_STENCIL_STAR);
 41:   DMDASetSizes(da,M,M,M);
 42:   DMDASetNumProcs(da,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
 43:   DMDASetDof(da,dof);
 44:   DMDASetStencilWidth(da,1);
 45:   DMDASetOwnershipRanges(da,NULL,NULL,NULL);
 46:   DMSetFromOptions(da);
 47:   DMSetUp(da);

 49:   DMCreateGlobalVector(da,&x);
 50:   DMCreateGlobalVector(da,&b);
 51:   ComputeRHS(da,b);
 52:   DMCreateMatrix(da,MATBAIJ,&A);
 53:   ComputeMatrix(da,A);


 56:   /* A is non-symmetric. Make A = 0.5*(A + Atrans) symmetric for testing icc and cholesky */
 57:   MatTranspose(A,MAT_INITIAL_MATRIX,&Atrans);
 58:   MatAXPY(A,1.0,Atrans,DIFFERENT_NONZERO_PATTERN);
 59:   MatScale(A,0.5);
 60:   MatDestroy(&Atrans);

 62:   /* Test sbaij matrix */
 63:   flg  = PETSC_FALSE;
 64:   PetscOptionsGetBool(NULL, "-test_sbaij1", &flg,NULL);
 65:   if (flg) {
 66:     Mat       sA;
 67:     PetscBool issymm;
 68:     MatIsTranspose(A,A,0.0,&issymm);
 69:     if (issymm) {
 70:       MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
 71:     } else printf("Warning: A is non-symmetric\n");
 72:     MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA);
 73:     MatDestroy(&A);
 74:     A    = sA;
 75:   }

 77:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 78:   KSPSetFromOptions(ksp);
 79:   KSPSetOperators(ksp,A,A,SAME_NONZERO_PATTERN);
 80:   KSPGetPC(ksp,&pc);
 81:   PCSetDM(pc,(DM)da);

 83:   if (trans) {
 84:     KSPSolveTranspose(ksp,b,x);
 85:   } else {
 86:     KSPSolve(ksp,b,x);
 87:   }

 89:   /* check final residual */
 90:   flg  = PETSC_FALSE;
 91:   PetscOptionsGetBool(NULL, "-check_final_residual", &flg,NULL);
 92:   if (flg) {
 93:     Vec       b1;
 94:     PetscReal norm;
 95:     KSPGetSolution(ksp,&x);
 96:     VecDuplicate(b,&b1);
 97:     MatMult(A,x,b1);
 98:     VecAXPY(b1,-1.0,b);
 99:     VecNorm(b1,NORM_2,&norm);
100:     PetscPrintf(PETSC_COMM_WORLD,"Final residual %g\n",norm);
101:     VecDestroy(&b1);
102:   }

104:   KSPDestroy(&ksp);
105:   VecDestroy(&x);
106:   VecDestroy(&b);
107:   MatDestroy(&A);
108:   DMDestroy(&da);
109:   PetscFinalize();
110:   return 0;
111: }

115: PetscErrorCode ComputeRHS(DM da,Vec b)
116: {
118:   PetscInt       mx,my,mz;
119:   PetscScalar    h;

122:   DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);
123:   h    = 1.0/((mx-1)*(my-1)*(mz-1));
124:   VecSet(b,h);
125:   return(0);
126: }

130: PetscErrorCode ComputeMatrix(DM da,Mat B)
131: {
133:   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs,dof,k1,k2,k3;
134:   PetscScalar    *v,*v_neighbor,Hx,Hy,Hz,HxHydHz,HyHzdHx,HxHzdHy;
135:   MatStencil     row,col;

138:   DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);
139:   /* For simplicity, this example only works on mx=my=mz */
140:   if (mx != my || mx != mz) SETERRQ3(PETSC_COMM_SELF,1,"This example only works with mx %d = my %d = mz %d\n",mx,my,mz);

142:   Hx      = 1.0 / (PetscReal)(mx-1); Hy = 1.0 / (PetscReal)(my-1); Hz = 1.0 / (PetscReal)(mz-1);
143:   HxHydHz = Hx*Hy/Hz; HxHzdHy = Hx*Hz/Hy; HyHzdHx = Hy*Hz/Hx;

145:   PetscMalloc((2*dof*dof+1)*sizeof(PetscScalar),&v);
146:   v_neighbor = v + dof*dof;
147:   PetscMemzero(v,(2*dof*dof+1)*sizeof(PetscScalar));
148:   k3         = 0;
149:   for (k1=0; k1<dof; k1++) {
150:     for (k2=0; k2<dof; k2++) {
151:       if (k1 == k2) {
152:         v[k3]          = 2.0*(HxHydHz + HxHzdHy + HyHzdHx);
153:         v_neighbor[k3] = -HxHydHz;
154:       } else {
155:         v[k3]          = k1/(dof*dof);;
156:         v_neighbor[k3] = k2/(dof*dof);
157:       }
158:       k3++;
159:     }
160:   }
161:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);

163:   for (k=zs; k<zs+zm; k++) {
164:     for (j=ys; j<ys+ym; j++) {
165:       for (i=xs; i<xs+xm; i++) {
166:         row.i = i; row.j = j; row.k = k;
167:         if (i==0 || j==0 || k==0 || i==mx-1 || j==my-1 || k==mz-1) { /* boudary points */
168:           MatSetValuesBlockedStencil(B,1,&row,1,&row,v,INSERT_VALUES);
169:         } else { /* interior points */
170:           /* center */
171:           col.i = i; col.j = j; col.k = k;
172:           MatSetValuesBlockedStencil(B,1,&row,1,&col,v,INSERT_VALUES);

174:           /* x neighbors */
175:           col.i = i-1; col.j = j; col.k = k;
176:           MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES);
177:           col.i = i+1; col.j = j; col.k = k;
178:           MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES);

180:           /* y neighbors */
181:           col.i = i; col.j = j-1; col.k = k;
182:           MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES);
183:           col.i = i; col.j = j+1; col.k = k;
184:           MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES);

186:           /* z neighbors */
187:           col.i = i; col.j = j; col.k = k-1;
188:           MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES);
189:           col.i = i; col.j = j; col.k = k+1;
190:           MatSetValuesBlockedStencil(B,1,&row,1,&col,v_neighbor,INSERT_VALUES);
191:         }
192:       }
193:     }
194:   }
195:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
196:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
197:   PetscFree(v);
198:   return(0);
199: }