Actual source code: ex3.c

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

 8:  #include petscksp.h

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

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

 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,&C);
 54:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 55:   MatSetFromOptions(C);
 56:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 57:   end   = start + M/size + ((M%size) > rank);

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

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

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

 96:   /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
 97:   PetscMalloc(4*m*sizeof(PetscInt),&rows);
 98:   for (i=0; i<m+1; i++) {
 99:     rows[i] = i; /* bottom */
100:     rows[3*m - 1 +i] = m*(m+1) + i; /* top */
101:   }
102:   count = m+1; /* left side */
103:   for (i=m+1; i<m*(m+1); i+= m+1) {
104:     rows[count++] = i;
105:   }
106:   count = 2*m; /* left side */
107:   for (i=2*m+1; i<m*(m+1); i+= m+1) {
108:     rows[count++] = i;
109:   }
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:   MatZeroRows(C,4*m,rows,1.0);

118:   PetscFree(rows);
119:   VecAssemblyBegin(u);
120:   VecAssemblyEnd(u);
121:   VecAssemblyBegin(b);
122:   VecAssemblyEnd(b);

124:   { Mat A;
125:   MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
126:   MatDestroy(C);
127:   MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);
128:   MatDestroy(A);
129:   }

131:   /* Solve linear system */
132:   KSPCreate(PETSC_COMM_WORLD,&ksp);
133:   KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
134:   KSPSetFromOptions(ksp);
135:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
136:   KSPSolve(ksp,b,u);

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

152:   /* Free work space */
153:   KSPDestroy(ksp);
154:   VecDestroy(ustar);
155:   VecDestroy(u);
156:   VecDestroy(b);
157:   MatDestroy(C);
158:   PetscFinalize();
159:   return 0;
160: }