Actual source code: ex16.c

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

  3: /* Usage:  mpirun ex16 [-help] [all PETSc options] */

  5: static char help[] = "Solves a sequence of linear systems with different right-hand-side vectors.n
  6: Input parameters include:n
  7:   -ntimes <ntimes>  : number of linear systems to solven
  8:   -view_exact_sol   : write exact solution vector to stdoutn
  9:   -m <mesh_x>       : number of mesh points in x-directionn
 10:   -n <mesh_n>       : number of mesh points in y-directionnn";

 12: /*T
 13:    Concepts: SLES^repeatedly solving linear systems;
 14:    Concepts: SLES^Laplacian, 2d
 15:    Concepts: Laplacian, 2d
 16:    Processors: n
 17: T*/

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

 29: int main(int argc,char **args)
 30: {
 31:   Vec         x,b,u;  /* approx solution, RHS, exact solution */
 32:   Mat         A;        /* linear system matrix */
 33:   SLES        sles;     /* linear solver context */
 34:   PetscReal   norm;     /* norm of solution error */
 35:   int         ntimes,i,j,k,I,J,Istart,Iend,ierr;
 36:   int         m = 8,n = 7,its;
 37:   PetscTruth  flg;
 38:   PetscScalar v,one = 1.0,neg_one = -1.0,rhs;

 40:   PetscInitialize(&argc,&args,(char *)0,help);
 41:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 42:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 44:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 45:          Compute the matrix for use in solving a series of
 46:          linear systems of the form, A x_i = b_i, for i=1,2,...
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 48:   /* 
 49:      Create parallel matrix, specifying only its global dimensions.
 50:      When using MatCreate(), the matrix format can be specified at
 51:      runtime. Also, the parallel partitioning of the matrix is
 52:      determined by PETSc at runtime.
 53:   */
 54:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&A);
 55:   MatSetFromOptions(A);

 57:   /* 
 58:      Currently, all PETSc parallel matrix formats are partitioned by
 59:      contiguous chunks of rows across the processors.  Determine which
 60:      rows of the matrix are locally owned. 
 61:   */
 62:   MatGetOwnershipRange(A,&Istart,&Iend);

 64:   /* 
 65:      Set matrix elements for the 2-D, five-point stencil in parallel.
 66:       - Each processor needs to insert only elements that it owns
 67:         locally (but any non-local elements will be sent to the
 68:         appropriate processor during matrix assembly). 
 69:       - Always specify global rows and columns of matrix entries.
 70:    */
 71:   for (I=Istart; I<Iend; I++) {
 72:     v = -1.0; i = I/n; j = I - i*n;
 73:     if (i>0)   {J = I - n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 74:     if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 75:     if (j>0)   {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 76:     if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 77:     v = 4.0; MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
 78:   }

 80:   /* 
 81:      Assemble matrix, using the 2-step process:
 82:        MatAssemblyBegin(), MatAssemblyEnd()
 83:      Computations can be done while messages are in transition
 84:      by placing code between these two statements.
 85:   */
 86:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 87:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 89:   /* 
 90:      Create parallel vectors.
 91:       - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
 92:         we specify only the vector's global
 93:         dimension; the parallel partitioning is determined at runtime. 
 94:       - When solving a linear system, the vectors and matrices MUST
 95:         be partitioned accordingly.  PETSc automatically generates
 96:         appropriately partitioned matrices and vectors when MatCreate()
 97:         and VecCreate() are used with the same communicator. 
 98:       - Note: We form 1 vector from scratch and then duplicate as needed.
 99:   */
100:   VecCreate(PETSC_COMM_WORLD,&u);
101:   VecSetSizes(u,PETSC_DECIDE,m*n);
102:   VecSetFromOptions(u);
103:   VecDuplicate(u,&b);
104:   VecDuplicate(b,&x);

106:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
107:                 Create the linear solver and set various options
108:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

110:   /* 
111:      Create linear solver context
112:   */
113:   SLESCreate(PETSC_COMM_WORLD,&sles);

115:   /* 
116:      Set operators. Here the matrix that defines the linear system
117:      also serves as the preconditioning matrix.
118:   */
119:   SLESSetOperators(sles,A,A,SAME_PRECONDITIONER);

121:   /* 
122:     Set runtime options, e.g.,
123:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
124:     These options will override those specified above as long as
125:     SLESSetFromOptions() is called _after_ any other customization
126:     routines.
127:   */
128:   SLESSetFromOptions(sles);

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
131:        Solve several linear systems of the form  A x_i = b_i
132:        I.e., we retain the same matrix (A) for all systems, but
133:        change the right-hand-side vector (b_i) at each step.

135:        In this case, we simply call SLESSolve() multiple times.  The
136:        preconditioner setup operations (e.g., factorization for ILU)
137:        be done during the first call to SLESSolve() only; such operations
138:        will NOT be repeated for successive solves.
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

141:   ntimes = 2;
142:   PetscOptionsGetInt(PETSC_NULL,"-ntimes",&ntimes,PETSC_NULL);
143:   for (k=1; k<ntimes+1; k++) {

145:     /* 
146:        Set exact solution; then compute right-hand-side vector.  We use
147:        an exact solution of a vector with all elements equal to 1.0*k.
148:     */
149:     rhs = one * (PetscReal)k;
150:     VecSet(&rhs,u);
151:     MatMult(A,u,b);

153:     /*
154:        View the exact solution vector if desired
155:     */
156:     PetscOptionsHasName(PETSC_NULL,"-view_exact_sol",&flg);
157:     if (flg) {VecView(u,PETSC_VIEWER_STDOUT_WORLD);}

159:     SLESSolve(sles,b,x,&its);

161:     /* 
162:        Check the error
163:     */
164:     VecAXPY(&neg_one,u,x);
165:     VecNorm(x,NORM_2,&norm);

167:     /*
168:        Print convergence information.  PetscPrintf() produces a single 
169:        print statement from all processes that share a communicator.
170:     */
171:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A System %d: iterations %dn",norm,k,its);
172:   }

174:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
175:                       Clean up
176:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177:   /* 
178:      Free work space.  All PETSc objects should be destroyed when they
179:      are no longer needed.
180:   */
181:   SLESDestroy(sles);
182:   VecDestroy(u);  VecDestroy(x);
183:   VecDestroy(b);  MatDestroy(A);

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