Actual source code: ex3.c

  2: static char help[] = "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: /*T
  9:    Concepts: KSP^basic parallel example
 10:    Concepts: Matrices^inserting elements by blocks
 11:    Processors: n
 12: T*/

 14: /* 
 15:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 16:   automatically includes:
 17:      petsc.h       - base PETSc routines   petscvec.h - vectors
 18:      petscsys.h    - system routines       petscmat.h - matrices
 19:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 20:      petscviewer.h - viewers               petscpc.h  - preconditioners
 21: */
 22:  #include petscksp.h

 24: /* Declare user-defined routines */

 30: int main(int argc,char **args)
 31: {
 32:   Vec            u,b,ustar; /* approx solution, RHS, exact solution */
 33:   Mat            A;           /* linear system matrix */
 34:   KSP            ksp;         /* Krylov subspace method context */
 35:   IS             is;          /* index set - used for boundary conditions */
 36:   PetscInt       N;           /* dimension of system (global) */
 37:   PetscInt       M;           /* number of elements (global) */
 38:   PetscMPIInt    rank;        /* processor rank */
 39:   PetscMPIInt    size;        /* size of communicator */
 40:   PetscScalar    Ke[16];      /* element matrix */
 41:   PetscScalar    r[4];        /* element vector */
 42:   PetscReal      h;           /* mesh width */
 43:   PetscReal      norm;        /* norm of solution error */
 44:   PetscReal      x,y;
 45:   PetscScalar    val,zero = 0.0,one = 1.0,none = -1.0;
 47:   PetscInt       idx[4],count,*rows,i,m = 5,start,end,its;

 49:   PetscInitialize(&argc,&args,(char *)0,help);
 50:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 51:   N = (m+1)*(m+1);
 52:   M = m*m;
 53:   h = 1.0/m;
 54:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 55:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 58:          Compute the matrix and right-hand-side vector that define
 59:          the linear system, Au = b.
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   /* 
 63:      Create stiffness matrix
 64:   */
 65:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
 66:   MatSetFromOptions(A);
 67:   start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
 68:   end   = start + M/size + ((M%size) > rank);

 70:   /*
 71:      Assemble matrix
 72:   */
 73:   FormElementStiffness(h*h,Ke);
 74:   for (i=start; i<end; i++) {
 75:      /* location of lower left corner of element */
 76:      x = h*(i % m); y = h*(i/m);
 77:      /* node numbers for the four corners of element */
 78:      idx[0] = (m+1)*(i/m) + (i % m);
 79:      idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
 80:      MatSetValues(A,4,idx,4,idx,Ke,ADD_VALUES);
 81:   }
 82:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 85:   /*
 86:      Create right-hand-side and solution vectors
 87:   */
 88:   VecCreate(PETSC_COMM_WORLD,&u);
 89:   VecSetSizes(u,PETSC_DECIDE,N);
 90:   VecSetFromOptions(u);
 91:   PetscObjectSetName((PetscObject)u,"Approx. Solution");
 92:   VecDuplicate(u,&b);
 93:   PetscObjectSetName((PetscObject)b,"Right hand side");
 94:   VecDuplicate(b,&ustar);
 95:   VecSet(&zero,u);
 96:   VecSet(&zero,b);

 98:   /* 
 99:      Assemble right-hand-side vector
100:   */
101:   for (i=start; i<end; i++) {
102:      /* location of lower left corner of element */
103:      x = h*(i % m); y = h*(i/m);
104:      /* node numbers for the four corners of element */
105:      idx[0] = (m+1)*(i/m) + (i % m);
106:      idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
107:      FormElementRhs(x,y,h*h,r);
108:      VecSetValues(b,4,idx,r,ADD_VALUES);
109:   }
110:   VecAssemblyBegin(b);
111:   VecAssemblyEnd(b);

113:   /* 
114:      Modify matrix and right-hand-side for Dirichlet boundary conditions
115:   */
116:   PetscMalloc(4*m*sizeof(PetscInt),&rows);
117:   for (i=0; i<m+1; i++) {
118:     rows[i] = i; /* bottom */
119:     rows[3*m - 1 +i] = m*(m+1) + i; /* top */
120:   }
121:   count = m+1; /* left side */
122:   for (i=m+1; i<m*(m+1); i+= m+1) {
123:     rows[count++] = i;
124:   }
125:   count = 2*m; /* left side */
126:   for (i=2*m+1; i<m*(m+1); i+= m+1) {
127:     rows[count++] = i;
128:   }
129:   ISCreateGeneral(PETSC_COMM_SELF,4*m,rows,&is);
130:   for (i=0; i<4*m; i++) {
131:      x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
132:      val = y;
133:      VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
134:      VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
135:   }
136:   PetscFree(rows);
137:   VecAssemblyBegin(u);
138:   VecAssemblyEnd(u);
139:   VecAssemblyBegin(b);
140:   VecAssemblyEnd(b);

142:   MatZeroRows(A,is,&one);
143:   ISDestroy(is);

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
146:                 Create the linear solver and set various options
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

149:   KSPCreate(PETSC_COMM_WORLD,&ksp);
150:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
151:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
152:   KSPSetFromOptions(ksp);

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
155:                       Solve the linear system
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

158:   KSPSolve(ksp,b,u);

160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
161:                       Check solution and clean up
162:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

164:   /* Check error */
165:   VecGetOwnershipRange(ustar,&start,&end);
166:   for (i=start; i<end; i++) {
167:      x = h*(i % (m+1)); y = h*(i/(m+1));
168:      val = y;
169:      VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
170:   }
171:   VecAssemblyBegin(ustar);
172:   VecAssemblyEnd(ustar);
173:   VecAXPY(&none,ustar,u);
174:   VecNorm(u,NORM_2,&norm);
175:   KSPGetIterationNumber(ksp,&its);
176:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %D\n",norm*h,its);

178:   /* 
179:      Free work space.  All PETSc objects should be destroyed when they
180:      are no longer needed.
181:   */
182:   KSPDestroy(ksp); VecDestroy(u);
183:   VecDestroy(ustar); VecDestroy(b);
184:   MatDestroy(A);

186:   /*
187:      Always call PetscFinalize() before exiting a program.  This routine
188:        - finalizes the PETSc libraries as well as MPI
189:        - provides summary and diagnostic information if certain runtime
190:          options are chosen (e.g., -log_summary).
191:   */
192:   PetscFinalize();
193:   return 0;
194: }

196: /* --------------------------------------------------------------------- */
199:    /* element stiffness for Laplacian */
200: PetscErrorCode FormElementStiffness(PetscReal H,PetscScalar *Ke)
201: {
203:   Ke[0]  = H/6.0;    Ke[1]  = -.125*H; Ke[2]  = H/12.0;   Ke[3]  = -.125*H;
204:   Ke[4]  = -.125*H;  Ke[5]  = H/6.0;   Ke[6]  = -.125*H;  Ke[7]  = H/12.0;
205:   Ke[8]  = H/12.0;   Ke[9]  = -.125*H; Ke[10] = H/6.0;    Ke[11] = -.125*H;
206:   Ke[12] = -.125*H;  Ke[13] = H/12.0;  Ke[14] = -.125*H;  Ke[15] = H/6.0;
207:   return(0);
208: }
209: /* --------------------------------------------------------------------- */
212: PetscErrorCode FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
213: {
215:   r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
216:   return(0);
217: }