Actual source code: ex1.c

  1: /*$Id: ex1.c,v 1.88 2001/03/23 23:23:55 balay 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:   double  norm;         /* norm of solution error */
 32:   int     ierr,i,n = 10,col[3],its,size;
 33:   Scalar  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,PETSC_DECIDE,n,&x);
 50:   VecSetFromOptions(x);
 51:   VecDuplicate(x,&b);
 52:   VecDuplicate(x,&u);

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

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

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

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

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

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

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

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

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

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

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