2.4. Simple PETSc Examples

Up: Contents Next: Include Files Previous: Writing PETSc Programs

To help the user start using PETSc immediately, we begin with a simple uniprocessor example in Figure 3 that solves the one-dimensional Laplacian problem with finite differences. This sequential code, which can be found in ${}PETSC_DIR/src/sles/examples/tutorials/ex1.c, illustrates the solution of a linear system with SLES, the interface to the preconditioners, Krylov subspace methods, and direct linear solvers of PETSc. Following the code we highlight a few of the most important parts of this example.


#ifdef PETSC_RCS_HEADER 
static char vcid[] = "$Id: ex1.c,v 1.73 1999/03/19 21:22:11 bsmith Exp $"; 
#endif 
 
/* Program usage:  mpirun ex1 [-help] [all PETSc options] */ 
 
static char help[] = "Solves a tridiagonal linear system with SLES.\n\n"; 
 
/*T 
   Concepts: SLES^Solving a system of linear equations (basic uniprocessor example); 
   Routines: SLESCreate(); SLESSetOperators(); SLESSetFromOptions(); 
   Routines: SLESSolve(); SLESView(); SLESGetKSP(); SLESGetPC(); 
   Routines: KSPSetTolerances(); PCSetType(); 
   Processors: 1 
T*/ 
 
/*  
  Include "sles.h" so that we can use SLES solvers.  Note that this file 
  automatically includes: 
     petsc.h  - base PETSc routines   vec.h - vectors 
     sys.h    - system routines       mat.h - matrices 
     is.h     - index sets            ksp.h - Krylov subspace methods 
     viewer.h - viewers               pc.h  - preconditioners 
*/ 
#include "sles.h" 
 
#undef __FUNC__ 
#define __FUNC__ "main" 
int main(int argc,char **args) 
{ 
  Vec     x, b, u;      /* approx solution, RHS, exact solution */ 
  Mat     A;            /* linear system matrix */ 
  SLES    sles;         /* linear solver context */ 
  PC      pc;           /* preconditioner context */ 
  KSP     ksp;          /* Krylov subspace method context */ 
  double  norm;         /* norm of solution error */ 
  int     ierr, i, n = 10, col[3], its, flg, size; 
  Scalar  neg_one = -1.0, one = 1.0, value[3]; 
 
  PetscInitialize(&argc,&args,(char *)0,help); 
  MPI_Comm_size(PETSC_COMM_WORLD,&size); 
  if (size != 1) SETERRA(1,0,"This is a uniprocessor example only!"); 
  ierr = OptionsGetInt(PETSC_NULL,"-n",&n,&flg); CHKERRA(ierr); 
 
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
         Compute the matrix and right-hand-side vector that define 
         the linear system, Ax = b. 
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
 
  /*  
     Create matrix.  When using MatCreate(), the matrix format can 
     be specified at runtime. 
 
     Performance tuning note:  For problems of substantial size, 
     preallocation of matrix memory is crucial for attaining good  
     performance.  Since preallocation is not possible via the generic 
     matrix creation routine MatCreate(), we recommend for practical  
     problems instead to use the creation routine for a particular matrix 
     format, e.g., 
         MatCreateSeqAIJ() - sequential AIJ (compressed sparse row) 
         MatCreateSeqBAIJ() - block AIJ 
     See the matrix chapter of the users manual for details. 
  */ 
  ierr = MatCreate(PETSC_COMM_WORLD,n,n,&A); CHKERRA(ierr); 
 
  /*  
     Assemble matrix 
  */ 
  value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 
  for (i=1; i<n-1; i++ ) { 
    col[0] = i-1; col[1] = i; col[2] = i+1; 
    ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES); CHKERRA(ierr); 
  } 
  i = n - 1; col[0] = n - 2; col[1] = n - 1; 
  ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES); CHKERRA(ierr); 
  i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0; 
  ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES); CHKERRA(ierr); 
  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRA(ierr); 
  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRA(ierr); 
 
  /*  
     Create vectors.  Note that we form 1 vector from scratch and 
     then duplicate as needed. 
  */ 
  ierr = VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,n,&x); CHKERRA(ierr); 
  ierr = VecSetFromOptions(x);CHKERRA(ierr); 
  ierr = VecDuplicate(x,&b); CHKERRA(ierr); 
  ierr = VecDuplicate(x,&u); CHKERRA(ierr); 
 
  /*  
     Set exact solution; then compute right-hand-side vector. 
  */ 
  ierr = VecSet(&one,u); CHKERRA(ierr); 
  ierr = MatMult(A,u,b); CHKERRA(ierr); 
 
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
                Create the linear solver and set various options 
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
  /*  
     Create linear solver context 
  */ 
  ierr = SLESCreate(PETSC_COMM_WORLD,&sles); CHKERRA(ierr); 
 
  /*  
     Set operators. Here the matrix that defines the linear system 
     also serves as the preconditioning matrix. 
  */ 
  ierr = SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN); CHKERRA(ierr); 
 
  /*  
     Set linear solver defaults for this problem (optional). 
     - By extracting the KSP and PC contexts from the SLES context, 
       we can then directly call any KSP and PC routines to set 
       various options. 
     - The following four statements are optional; all of these 
       parameters could alternatively be specified at runtime via 
       SLESSetFromOptions(); 
  */ 
  ierr = SLESGetKSP(sles,&ksp); CHKERRA(ierr); 
  ierr = SLESGetPC(sles,&pc); CHKERRA(ierr); 
  ierr = PCSetType(pc,PCJACOBI); CHKERRA(ierr); 
  ierr = KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT, 
         PETSC_DEFAULT); CHKERRA(ierr); 
 
  /*  
    Set runtime options, e.g., 
        -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol> 
    These options will override those specified above as long as 
    SLESSetFromOptions() is called _after_ any other customization 
    routines. 
  */ 
  ierr = SLESSetFromOptions(sles); CHKERRA(ierr); 
  
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
                      Solve the linear system 
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
  /*  
     Solve linear system 
  */ 
  ierr = SLESSolve(sles,b,x,&its); CHKERRA(ierr);  
 
  /*  
     View solver info; we could instead use the option -sles_view to 
     print this info to the screen at the conclusion of SLESSolve(). 
  */ 
  ierr = SLESView(sles,VIEWER_STDOUT_WORLD); CHKERRA(ierr); 
 
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
                      Check solution and clean up 
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
  /*  
     Check the error 
  */ 
  ierr = VecAXPY(&neg_one,u,x); CHKERRA(ierr); 
  ierr  = VecNorm(x,NORM_2,&norm); CHKERRA(ierr); 
  if (norm > 1.e-12)  
    PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %d\n",norm,its); 
  else  
    PetscPrintf(PETSC_COMM_WORLD,"Norm of error < 1.e-12, Iterations %d\n",its); 
 
  /*  
     Free work space.  All PETSc objects should be destroyed when they 
     are no longer needed. 
  */ 
  ierr = VecDestroy(x); CHKERRA(ierr); ierr = VecDestroy(u); CHKERRA(ierr); 
  ierr = VecDestroy(b); CHKERRA(ierr); ierr = MatDestroy(A); CHKERRA(ierr); 
  ierr = SLESDestroy(sles); CHKERRA(ierr); 
 
  /* 
     Always call PetscFinalize() before exiting a program.  This routine 
       - finalizes the PETSc libraries as well as MPI 
       - provides summary and diagnostic information if certain runtime 
         options are chosen (e.g., -log_summary). 
  */ 
  PetscFinalize(); 
  return 0; 
} 

Figure 3: Example of Uniprocessor PETSc Code


Up: Contents Next: Include Files Previous: Writing PETSc Programs