Actual source code: ex4.c

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

  3: static char help[] = "Solves a linear system with SLES.  The matrix uses simplen
  4: bilinear elements on the unit square. Input arguments are:n
  5:   -m <size> : problem sizenn";

 7:  #include petscsles.h

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

 23: int main(int argc,char **args)
 24: {
 25:   Mat         C;
 26:   int         i,m = 2,N,M,its,ierr,idx[4],count,*rows;
 27:   PetscScalar val,zero = 0.0,one = 1.0,none = -1.0,Ke[16],r[4];
 28:   PetscReal   x,y,h,norm;
 29:   Vec         u,ustar,b;
 30:   SLES        sles;
 31:   KSP         ksp;
 32:   IS          is;

 34:   PetscInitialize(&argc,&args,(char *)0,help);
 35:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 36:   N = (m+1)*(m+1); /* dimension of matrix */
 37:   M = m*m; /* number of elements */
 38:   h = 1.0/m;       /* mesh width */

 40:   /* create stiffness matrix */
 41:   MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,9,PETSC_NULL,&C);

 43:   /* forms the element stiffness for the Laplacian */
 44:   FormElementStiffness(h*h,Ke);
 45:   for (i=0; i<M; i++) {
 46:      /* location of lower left corner of element */
 47:      x = h*(i % m); y = h*(i/m);
 48:      /* node numbers for the four corners of element */
 49:      idx[0] = (m+1)*(i/m) + (i % m);
 50:      idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 51:      MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
 52:   }
 53:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 54:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 56:   /* create right hand side and solution */

 58:   VecCreateSeq(PETSC_COMM_SELF,N,&u);
 59:   VecDuplicate(u,&b);
 60:   VecDuplicate(b,&ustar);
 61:   VecSet(&zero,u);
 62:   VecSet(&zero,b);

 64:   for (i=0; i<M; i++) {
 65:      /* location of lower left corner of element */
 66:      x = h*(i % m); y = h*(i/m);
 67:      /* node numbers for the four corners of element */
 68:      idx[0] = (m+1)*(i/m) + (i % m);
 69:      idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 70:      FormElementRhs(x,y,h*h,r);
 71:      VecSetValues(b,4,idx,r,ADD_VALUES);
 72:   }
 73:   VecAssemblyBegin(b);
 74:   VecAssemblyEnd(b);

 76:   /* modify matrix and rhs for Dirichlet boundary conditions */
 77:   PetscMalloc((4*m+1)*sizeof(int),&rows);
 78:   for (i=0; i<m+1; i++) {
 79:     rows[i] = i; /* bottom */
 80:     rows[3*m - 1 +i] = m*(m+1) + i; /* top */
 81:   }
 82:   count = m+1; /* left side */
 83:   for (i=m+1; i<m*(m+1); i+= m+1) {
 84:     rows[count++] = i;
 85:   }
 86:   count = 2*m; /* left side */
 87:   for (i=2*m+1; i<m*(m+1); i+= m+1) {
 88:     rows[count++] = i;
 89:   }
 90:   ISCreateGeneral(PETSC_COMM_SELF,4*m,rows,&is);
 91:   for (i=0; i<4*m; i++) {
 92:      x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
 93:      val = y;
 94:      VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
 95:      VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
 96:   }
 97:   PetscFree(rows);
 98:   VecAssemblyBegin(u);
 99:   VecAssemblyEnd(u);
100:   VecAssemblyBegin(b);
101:   VecAssemblyEnd(b);

103:   MatZeroRows(C,is,&one);
104:   ISDestroy(is);

106:   /* solve linear system */
107:   SLESCreate(PETSC_COMM_WORLD,&sles);
108:   SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
109: 
110:   SLESSetFromOptions(sles);
111:   SLESGetKSP(sles,&ksp);
112:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
113:   SLESSolve(sles,b,u,&its);

115:   /* check error */
116:   for (i=0; i<N; i++) {
117:      x = h*(i % (m+1)); y = h*(i/(m+1));
118:      val = y;
119:      VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
120:   }
121:   VecAssemblyBegin(ustar);
122:   VecAssemblyEnd(ustar);

124:   VecAXPY(&none,ustar,u);
125:   VecNorm(u,NORM_2,&norm);
126:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %dn",norm*h,its);

128:   SLESDestroy(sles);
129:   VecDestroy(ustar);
130:   VecDestroy(u);
131:   VecDestroy(b);
132:   MatDestroy(C);
133:   PetscFinalize();
134:   return 0;
135: }