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

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

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

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

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

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

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

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

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

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

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

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

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

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

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