Actual source code: ex3.c

  1: /*$Id: ex3.c,v 1.73 2001/08/07 21:30:50 bsmith Exp $*/

  3: static char help[] = "This example solves a linear system in parallel with KSP.  The matrix\n\
  4: uses simple bilinear elements on the unit square.  To test the parallel\n\
  5: matrix assembly, the matrix is intentionally laid out across processors\n\
  6: differently from the way it is assembled.  Input arguments are:\n\
  7:   -m <size> : problem size\n\n";

 9:  #include petscksp.h

 13: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
 14: {
 16:   Ke[0]  = H/6.0;    Ke[1]  = -.125*H; Ke[2]  = H/12.0;   Ke[3]  = -.125*H;
 17:   Ke[4]  = -.125*H;  Ke[5]  = H/6.0;   Ke[6]  = -.125*H;  Ke[7]  = H/12.0;
 18:   Ke[8]  = H/12.0;   Ke[9]  = -.125*H; Ke[10] = H/6.0;    Ke[11] = -.125*H;
 19:   Ke[12] = -.125*H;  Ke[13] = H/12.0;  Ke[14] = -.125*H;  Ke[15] = H/6.0;
 20:   return(0);
 21: }
 24: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
 25: {
 27:   r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
 28:   return(0);
 29: }

 33: int main(int argc,char **args)
 34: {
 35:   Mat         C;
 36:   int         i,m = 5,rank,size,N,start,end,M,its;
 37:   PetscScalar val,zero = 0.0,one = 1.0,none = -1.0,Ke[16],r[4];
 38:   PetscReal   x,y,h,norm;
 39:   int         ierr,idx[4],count,*rows;
 40:   Vec         u,ustar,b;
 41:   KSP         ksp;
 42:   IS          is;

 44:   PetscInitialize(&argc,&args,(char *)0,help);
 45:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 46:   N = (m+1)*(m+1); /* dimension of matrix */
 47:   M = m*m; /* number of elements */
 48:   h = 1.0/m;       /* mesh width */
 49:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 50:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 52:   /* Create stiffness matrix */
 53:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&C);
 54:   MatSetFromOptions(C);
 55:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 56:   end   = start + M/size + ((M%size) > rank);

 58:   /* Assemble matrix */
 59:   FormElementStiffness(h*h,Ke);   /* element stiffness for Laplacian */
 60:   for (i=start; i<end; i++) {
 61:      /* location of lower left corner of element */
 62:      x = h*(i % m); y = h*(i/m);
 63:      /* node numbers for the four corners of element */
 64:      idx[0] = (m+1)*(i/m) + (i % m);
 65:      idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 66:      MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
 67:   }
 68:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 69:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 71:   /* Create right-hand-side and solution vectors */
 72:   VecCreate(PETSC_COMM_WORLD,&u);
 73:   VecSetSizes(u,PETSC_DECIDE,N);
 74:   VecSetFromOptions(u);
 75:   PetscObjectSetName((PetscObject)u,"Approx. Solution");
 76:   VecDuplicate(u,&b);
 77:   PetscObjectSetName((PetscObject)b,"Right hand side");
 78:   VecDuplicate(b,&ustar);
 79:   VecSet(&zero,u);
 80:   VecSet(&zero,b);

 82:   /* Assemble right-hand-side vector */
 83:   for (i=start; i<end; i++) {
 84:      /* location of lower left corner of element */
 85:      x = h*(i % m); y = h*(i/m);
 86:      /* node numbers for the four corners of element */
 87:      idx[0] = (m+1)*(i/m) + (i % m);
 88:      idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 89:      FormElementRhs(x,y,h*h,r);
 90:      VecSetValues(b,4,idx,r,ADD_VALUES);
 91:   }
 92:   VecAssemblyBegin(b);
 93:   VecAssemblyEnd(b);

 95:   /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
 96:   PetscMalloc(4*m*sizeof(int),&rows);
 97:   for (i=0; i<m+1; i++) {
 98:     rows[i] = i; /* bottom */
 99:     rows[3*m - 1 +i] = m*(m+1) + i; /* top */
100:   }
101:   count = m+1; /* left side */
102:   for (i=m+1; i<m*(m+1); i+= m+1) {
103:     rows[count++] = i;
104:   }
105:   count = 2*m; /* left side */
106:   for (i=2*m+1; i<m*(m+1); i+= m+1) {
107:     rows[count++] = i;
108:   }
109:   ISCreateGeneral(PETSC_COMM_SELF,4*m,rows,&is);
110:   for (i=0; i<4*m; i++) {
111:      x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
112:      val = y;
113:      VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
114:      VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
115:   }
116:   PetscFree(rows);
117:   VecAssemblyBegin(u);
118:   VecAssemblyEnd(u);
119:   VecAssemblyBegin(b);
120:   VecAssemblyEnd(b);

122:   MatZeroRows(C,is,&one);
123:   ISDestroy(is);


126:   { Mat A;
127:   MatConvert(C,MATSAME,&A);
128:   MatDestroy(C);
129:   MatConvert(A,MATSAME,&C);
130:   MatDestroy(A);
131:   }

133:   /* Solve linear system */
134:   KSPCreate(PETSC_COMM_WORLD,&ksp);
135:   KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
136:   KSPSetFromOptions(ksp);
137:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
138:   KSPSetRhs(ksp,b);
139:   KSPSetSolution(ksp,u);
140:   KSPSolve(ksp);

142:   /* Check error */
143:   VecGetOwnershipRange(ustar,&start,&end);
144:   for (i=start; i<end; i++) {
145:      x = h*(i % (m+1)); y = h*(i/(m+1));
146:      val = y;
147:      VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
148:   }
149:   VecAssemblyBegin(ustar);
150:   VecAssemblyEnd(ustar);
151:   VecAXPY(&none,ustar,u);
152:   VecNorm(u,NORM_2,&norm);
153:   KSPGetIterationNumber(ksp,&its);
154:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %d\n",norm*h,its);

156:   /* Free work space */
157:   KSPDestroy(ksp);
158:   VecDestroy(ustar);
159:   VecDestroy(u);
160:   VecDestroy(b);
161:   MatDestroy(C);
162:   PetscFinalize();
163:   return 0;
164: }