Actual source code: ex1.c

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

  3: /* Program usage:  mpirun ex1 [-help] [all PETSc options] */

  5: static char help[] = "Solves a tridiagonal linear system with SLES.nn";

  7: /*T
  8:    Concepts: SLES^solving a system of linear equations
  9:    Processors: 1
 10: T*/

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

 20:   Note:  The corresponding parallel example is ex23.c
 21: */
 22:  #include petscsles.h

 24: int main(int argc,char **args)
 25: {
 26:   Vec         x, b, u;      /* approx solution, RHS, exact solution */
 27:   Mat         A;            /* linear system matrix */
 28:   SLES        sles;         /* linear solver context */
 29:   PC          pc;           /* preconditioner context */
 30:   KSP         ksp;          /* Krylov subspace method context */
 31:   PetscReal   norm;         /* norm of solution error */
 32:   int         ierr,i,n = 10,col[3],its,size;
 33:   PetscScalar neg_one = -1.0,one = 1.0,value[3];

 35:   PetscInitialize(&argc,&args,(char *)0,help);
 36:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 37:   if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
 38:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 40:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 41:          Compute the matrix and right-hand-side vector that define
 42:          the linear system, Ax = b.
 43:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 45:   /* 
 46:      Create vectors.  Note that we form 1 vector from scratch and
 47:      then duplicate as needed.
 48:   */
 49:   VecCreate(PETSC_COMM_WORLD,&x);
 50:   PetscObjectSetName((PetscObject) x, "Solution");
 51:   VecSetSizes(x,PETSC_DECIDE,n);
 52:   VecSetFromOptions(x);
 53:   VecDuplicate(x,&b);
 54:   VecDuplicate(x,&u);

 56:   /* 
 57:      Create matrix.  When using MatCreate(), the matrix format can
 58:      be specified at runtime.

 60:      Performance tuning note:  For problems of substantial size,
 61:      preallocation of matrix memory is crucial for attaining good 
 62:      performance.  Since preallocation is not possible via the generic
 63:      matrix creation routine MatCreate(), we recommend for practical 
 64:      problems instead to use the creation routine for a particular matrix
 65:      format, e.g.,
 66:          MatCreateSeqAIJ() - sequential AIJ (compressed sparse row)
 67:          MatCreateSeqBAIJ() - block AIJ
 68:      See the matrix chapter of the users manual for details.
 69:   */
 70:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A);
 71:   MatSetFromOptions(A);

 73:   /* 
 74:      Assemble matrix
 75:   */
 76:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 77:   for (i=1; i<n-1; i++) {
 78:     col[0] = i-1; col[1] = i; col[2] = i+1;
 79:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 80:   }
 81:   i = n - 1; col[0] = n - 2; col[1] = n - 1;
 82:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 83:   i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 84:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 85:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 86:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 88:   /* 
 89:      Set exact solution; then compute right-hand-side vector.
 90:   */
 91:   VecSet(&one,u);
 92:   MatMult(A,u,b);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 95:                 Create the linear solver and set various options
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   /* 
 98:      Create linear solver context
 99:   */
100:   SLESCreate(PETSC_COMM_WORLD,&sles);

102:   /* 
103:      Set operators. Here the matrix that defines the linear system
104:      also serves as the preconditioning matrix.
105:   */
106:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);

108:   /* 
109:      Set linear solver defaults for this problem (optional).
110:      - By extracting the KSP and PC contexts from the SLES context,
111:        we can then directly call any KSP and PC routines to set
112:        various options.
113:      - The following four statements are optional; all of these
114:        parameters could alternatively be specified at runtime via
115:        SLESSetFromOptions();
116:   */
117:   SLESGetKSP(sles,&ksp);
118:   SLESGetPC(sles,&pc);
119:   PCSetType(pc,PCJACOBI);
120:   KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

122:   /* 
123:     Set runtime options, e.g.,
124:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
125:     These options will override those specified above as long as
126:     SLESSetFromOptions() is called _after_ any other customization
127:     routines.
128:   */
129:   SLESSetFromOptions(sles);
130: 
131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
132:                       Solve the linear system
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134:   /* 
135:      Solve linear system
136:   */
137:   SLESSolve(sles,b,x,&its);

139:   /* 
140:      View solver info; we could instead use the option -sles_view to
141:      print this info to the screen at the conclusion of SLESSolve().
142:   */
143:   SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
146:                       Check solution and clean up
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:   /* 
149:      Check the error
150:   */
151:   VecAXPY(&neg_one,u,x);
152:   ierr  = VecNorm(x,NORM_2,&norm);
153:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %dn",norm,its);
154:   /* 
155:      Free work space.  All PETSc objects should be destroyed when they
156:      are no longer needed.
157:   */
158:   VecDestroy(x); VecDestroy(u);
159:   VecDestroy(b); MatDestroy(A);
160:   SLESDestroy(sles);

162:   /*
163:      Always call PetscFinalize() before exiting a program.  This routine
164:        - finalizes the PETSc libraries as well as MPI
165:        - provides summary and diagnostic information if certain runtime
166:          options are chosen (e.g., -log_summary).
167:   */
168:   PetscFinalize();
169:   return 0;
170: }