Actual source code: ex3.c

  1: /*$Id: ex3.c,v 1.69 2001/04/10 19:36:37 bsmith Exp $*/

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

 9:  #include petscsles.h

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

 27: int main(int argc,char **args)
 28: {
 29:   Mat         C;
 30:   int         i,m = 5,rank,size,N,start,end,M,its;
 31:   Scalar      val,zero = 0.0,one = 1.0,none = -1.0,Ke[16],r[4];
 32:   double      x,y,h,norm;
 33:   int         ierr,idx[4],count,*rows;
 34:   Vec         u,ustar,b;
 35:   SLES        sles;
 36:   KSP         ksp;
 37:   IS          is;

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

 47:   /* Create stiffness matrix */
 48:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&C);
 49:   MatSetFromOptions(C);
 50:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 51:   end   = start + M/size + ((M%size) > rank);

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

 66:   /* Create right-hand-side and solution vectors */
 67:   VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,N,&u);
 68:   VecSetFromOptions(u);
 69:   PetscObjectSetName((PetscObject)u,"Approx. Solution");
 70:   VecDuplicate(u,&b);
 71:   PetscObjectSetName((PetscObject)b,"Right hand side");
 72:   VecDuplicate(b,&ustar);
 73:   VecSet(&zero,u);
 74:   VecSet(&zero,b);

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

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

116:   MatZeroRows(C,is,&one);
117:   ISDestroy(is);


120:   { Mat A;
121:   MatConvert(C,MATSAME,&A);
122:   MatDestroy(C);
123:   MatConvert(A,MATSAME,&C);
124:   MatDestroy(A);
125:   }

127:   /* Solve linear system */
128:   SLESCreate(PETSC_COMM_WORLD,&sles);
129:   SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
130:   SLESSetFromOptions(sles);
131:   SLESGetKSP(sles,&ksp);
132:   KSPSetInitialGuessNonzero(ksp);
133:   SLESSolve(sles,b,u,&its);

135:   /* Check error */
136:   VecGetOwnershipRange(ustar,&start,&end);
137:   for (i=start; i<end; i++) {
138:      x = h*(i % (m+1)); y = h*(i/(m+1));
139:      val = y;
140:      VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
141:   }
142:   VecAssemblyBegin(ustar);
143:   VecAssemblyEnd(ustar);
144:   VecAXPY(&none,ustar,u);
145:   VecNorm(u,NORM_2,&norm);
146:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %dn",norm*h,its);

148:   /* Free work space */
149:   SLESDestroy(sles);
150:   VecDestroy(ustar);
151:   VecDestroy(u);
152:   VecDestroy(b);
153:   MatDestroy(C);
154:   PetscFinalize();
155:   return 0;
156: }