Actual source code: ex7.c

  1: /*$Id: ex7.c,v 1.17 2001/04/10 19:36:37 bsmith Exp $*/

  3: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.n
  4:  Tests inplace factorization for SeqBAIJ. Input parameters includen
  5:   -f0 <input_file> : first file to load (small system)nn";

  7: /*T
  8:    Concepts: SLES^solving a linear system
  9:    Concepts: PetscLog^profiling multiple stages of code;
 10:    Processors: n
 11: T*/

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

 23: int main(int argc,char **args)
 24: {
 25:   SLES        sles;             /* linear solver context */
 26:   Mat         A,B;                /* matrix */
 27:   Vec         x,b,u;          /* approx solution, RHS, exact solution */
 28:   PetscViewer fd;               /* viewer */
 29:   char        file[2][128];     /* input file name */
 30:   int         ierr,its;
 31:   PetscTruth  flg;
 32:   double      norm;
 33:   Scalar      zero = 0.0,none = -1.0;

 35:   PetscInitialize(&argc,&args,(char *)0,help);

 37:   /* 
 38:      Determine files from which we read the two linear systems
 39:      (matrix and right-hand-side vector).
 40:   */
 41:   PetscOptionsGetString(PETSC_NULL,"-f0",file[0],127,&flg);
 42:   if (!flg) SETERRQ(1,"Must indicate binary file with the -f0 option");


 45:   /* 
 46:        Open binary file.  Note that we use PETSC_BINARY_RDONLY to indicate
 47:        reading from this file.
 48:   */
 49:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],PETSC_BINARY_RDONLY,&fd);

 51:   /*
 52:        Load the matrix and vector; then destroy the viewer.
 53:   */
 54:   MatLoad(fd,MATSEQBAIJ,&A);
 55:   MatConvert(A,MATSAME,&B);
 56:   VecLoad(fd,&b);
 57:   PetscViewerDestroy(fd);

 59:   /* 
 60:        If the loaded matrix is larger than the vector (due to being padded 
 61:        to match the block size of the system), then create a new padded vector.
 62:   */
 63:   {
 64:       int    m,n,j,mvec,start,end,index;
 65:       Vec    tmp;
 66:       Scalar *bold;

 68:       /* Create a new vector b by padding the old one */
 69:       MatGetLocalSize(A,&m,&n);
 70:       VecCreateMPI(PETSC_COMM_WORLD,m,PETSC_DECIDE,&tmp);
 71:       VecGetOwnershipRange(b,&start,&end);
 72:       VecGetLocalSize(b,&mvec);
 73:       VecGetArray(b,&bold);
 74:       for (j=0; j<mvec; j++) {
 75:         index = start+j;
 76:         ierr  = VecSetValues(tmp,1,&index,bold+j,INSERT_VALUES);
 77:       }
 78:       VecRestoreArray(b,&bold);
 79:       VecDestroy(b);
 80:       VecAssemblyBegin(tmp);
 81:       VecAssemblyEnd(tmp);
 82:       b = tmp;
 83:     }
 84:   VecDuplicate(b,&x);
 85:   VecDuplicate(b,&u);
 86:   VecSet(&zero,x);

 88:   /*
 89:       Create linear solver; set operators; set runtime options.
 90:   */
 91:   SLESCreate(PETSC_COMM_WORLD,&sles);
 92:   SLESSetOperators(sles,A,B,DIFFERENT_NONZERO_PATTERN);
 93:   SLESSetFromOptions(sles);

 95:   /* 
 96:        Here we explicitly call SLESSetUp() and SLESSetUpOnBlocks() to
 97:        enable more precise profiling of setting up the preconditioner.
 98:        These calls are optional, since both will be called within
 99:        SLESSolve() if they haven't been called already.
100:   */
101:   SLESSetUp(sles,b,x);
102:   SLESSetUpOnBlocks(sles);

104:   SLESSolve(sles,b,x,&its);

106:   /*
107:             Check error, print output, free data structures.
108:             This stage is not profiled separately.
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   /* 
112:      Check error
113:   */
114:   MatMult(A,x,u);
115:   VecAXPY(&none,b,u);
116:   VecNorm(u,NORM_2,&norm);

118:   PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3dn",its);
119:   PetscPrintf(PETSC_COMM_WORLD,"Residual norm = %An",norm);

121:   /* 
122:        Free work space.  All PETSc objects should be destroyed when they
123:        are no longer needed.
124:   */
125:   MatDestroy(A);
126:   MatDestroy(B);
127:   VecDestroy(b);
128:   VecDestroy(u); VecDestroy(x);
129:   SLESDestroy(sles);


132:   PetscFinalize();
133:   return 0;
134: }