Actual source code: ex3.c

  1: /*$Id: ex3.c,v 1.33 2001/08/07 21:30:54 bsmith Exp $*/

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

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

 25: /* Declare user-defined routines */
 26: extern int FormElementStiffness(PetscReal,PetscScalar*);
 27: extern int FormElementRhs(PetscReal,PetscReal,PetscReal,PetscScalar*);

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

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

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

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

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

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

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

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

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

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

150:   SLESCreate(PETSC_COMM_WORLD,&sles);
151:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
152:   SLESGetKSP(sles,&ksp);
153:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
154:   SLESSetFromOptions(sles);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
157:                       Solve the linear system
158:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

160:   SLESSolve(sles,b,u,&its);

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
163:                       Check solution and clean up
164:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

197: /* --------------------------------------------------------------------- */
198: #undef __FUNCT__
200:    /* element stiffness for Laplacian */
201: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
202: {
204:   Ke[0]  = H/6.0;    Ke[1]  = -.125*H; Ke[2]  = H/12.0;   Ke[3]  = -.125*H;
205:   Ke[4]  = -.125*H;  Ke[5]  = H/6.0;   Ke[6]  = -.125*H;  Ke[7]  = H/12.0;
206:   Ke[8]  = H/12.0;   Ke[9]  = -.125*H; Ke[10] = H/6.0;    Ke[11] = -.125*H;
207:   Ke[12] = -.125*H;  Ke[13] = H/12.0;  Ke[14] = -.125*H;  Ke[15] = H/6.0;
208:   return(0);
209: }
210: /* --------------------------------------------------------------------- */
211: #undef __FUNCT__
213: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
214: {
216:   r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
217:   return(0);
218: }