Actual source code: ex1.c

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

  4: static char help[] = "Solves a tridiagonal linear system with KSP.\n\n";

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

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

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

 25: int main(int argc,char **args)
 26: {
 27:   Vec            x, b, u;      /* approx solution, RHS, exact solution */
 28:   Mat            A;            /* linear system matrix */
 29:   KSP            ksp;         /* linear solver context */
 30:   PC             pc;           /* preconditioner context */
 31:   PetscReal      norm;         /* norm of solution error */
 33:   PetscInt       i,n = 10,col[3],its;
 34:   PetscMPIInt    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,&A);
 73:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 74:   MatSetFromOptions(A);

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

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

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

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

111:   /* 
112:      Set linear solver defaults for this problem (optional).
113:      - By extracting the KSP and PC contexts from the KSP context,
114:        we can then directly call any KSP and PC routines to set
115:        various options.
116:      - The following four statements are optional; all of these
117:        parameters could alternatively be specified at runtime via
118:        KSPSetFromOptions();
119:   */
120:   KSPGetPC(ksp,&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:     KSPSetFromOptions() is called _after_ any other customization
129:     routines.
130:   */
131:   KSPSetFromOptions(ksp);
132: 
133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
134:                       Solve the linear system
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136:   /* 
137:      Solve linear system
138:   */
139:   KSPSolve(ksp,b,x);

141:   /* 
142:      View solver info; we could instead use the option -ksp_view to
143:      print this info to the screen at the conclusion of KSPSolve().
144:   */
145:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

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

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