Actual source code: ex32.c

petsc-dev 2014-02-02
Report Typos and Errors
  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:   DMSetMatType(da,MATBAIJ);
 53:   DMCreateMatrix(da,&A);
 54:   ComputeMatrix(da,A);


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

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

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

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

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

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

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

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

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

139:   DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,&dof,0,0,0,0,0);
140:   /* For simplicity, this example only works on mx=my=mz */
141:   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);

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

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

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

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

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

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