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 KSP.  The matrix uses simple\n\
  4: bilinear elements on the unit square. Input arguments are:\n\
  5:   -m <size> : problem size\n\n";

 7:  #include petscksp.h

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

 29: int main(int argc,char **args)
 30: {
 31:   Mat         C;
 32:   int         i,m = 2,N,M,its,ierr,idx[4],count,*rows;
 33:   PetscScalar val,zero = 0.0,one = 1.0,none = -1.0,Ke[16],r[4];
 34:   PetscReal   x,y,h,norm;
 35:   Vec         u,ustar,b;
 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 */

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

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

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

 63:   VecCreateSeq(PETSC_COMM_SELF,N,&u);
 64:   VecDuplicate(u,&b);
 65:   VecDuplicate(b,&ustar);
 66:   VecSet(&zero,u);
 67:   VecSet(&zero,b);

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

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

108:   MatZeroRows(C,is,&one);
109:   ISDestroy(is);

111:   /* solve linear system */
112:   KSPCreate(PETSC_COMM_WORLD,&ksp);
113:   KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
114: 
115:   KSPSetFromOptions(ksp);
116:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
117:   KSPSetRhs(ksp,b);
118:   KSPSetSolution(ksp,u);
119:   KSPSolve(ksp);

121:   /* check error */
122:   for (i=0; i<N; i++) {
123:      x = h*(i % (m+1)); y = h*(i/(m+1));
124:      val = y;
125:      VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
126:   }
127:   VecAssemblyBegin(ustar);
128:   VecAssemblyEnd(ustar);

130:   VecAXPY(&none,ustar,u);
131:   VecNorm(u,NORM_2,&norm);
132:   KSPGetIterationNumber(ksp,&its);
133:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %d\n",norm*h,its);

135:   KSPDestroy(ksp);
136:   VecDestroy(ustar);
137:   VecDestroy(u);
138:   VecDestroy(b);
139:   MatDestroy(C);
140:   PetscFinalize();
141:   return 0;
142: }