Actual source code: ex18.c

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

  3: #if !defined(PETSC_USE_COMPLEX)

  5: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.n
  6: Input arguments are:n
  7:   -f <input_file> : file to load.  For a 5X5 example of the 5-pt. stencil,n
  8:                     use the file petsc/src/mat/examples/matbinary.exnn";

 10:  #include petscmat.h
 11:  #include petscsles.h

 13: int main(int argc,char **args)
 14: {
 15:   int            ierr,its,m,n,mvec;
 16:   PetscLogDouble time1,time2,time;
 17:   double         norm;
 18:   Scalar         zero = 0.0,none = -1.0;
 19:   Vec            x,b,u;
 20:   Mat            A;
 21:   SLES           sles;
 22:   char           file[128];
 23:   PetscViewer    fd;

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

 27:   /* Read matrix and RHS */
 28:   PetscOptionsGetString(PETSC_NULL,"-f",file,127,PETSC_NULL);
 29:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_BINARY_RDONLY,&fd);
 30:   MatLoad(fd,MATSEQAIJ,&A);
 31:   VecLoad(fd,&b);
 32:   PetscViewerDestroy(fd);

 34:   /* 
 35:      If the load matrix is larger then the vector, due to being padded 
 36:      to match the blocksize then create a new padded vector
 37:   */
 38:   MatGetSize(A,&m,&n);
 39:   VecGetSize(b,&mvec);
 40:   if (m > mvec) {
 41:     Vec    tmp;
 42:     Scalar *bold,*bnew;
 43:     /* create a new vector b by padding the old one */
 44:     VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,m,&tmp);
 45:     VecSetFromOptions(tmp);
 46:     VecGetArray(tmp,&bnew);
 47:     VecGetArray(b,&bold);
 48:     PetscMemcpy(bnew,bold,mvec*sizeof(Scalar));
 49:     VecDestroy(b);
 50:     b = tmp;
 51:   }

 53:   /* Set up solution */
 54:   VecDuplicate(b,&x);
 55:   VecDuplicate(b,&u);
 56:   VecSet(&zero,x);

 58:   /* Solve system */
 59:   PetscLogStagePush(1);
 60:   SLESCreate(PETSC_COMM_WORLD,&sles);
 61:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
 62:   SLESSetFromOptions(sles);
 63:   PetscGetTime(&time1);
 64:   SLESSolve(sles,b,x,&its);
 65:   PetscGetTime(&time2);
 66:   time = time2 - time1;
 67:   PetscLogStagePop();

 69:   /* Show result */
 70:   MatMult(A,x,u);
 71:   VecAXPY(&none,b,u);
 72:   VecNorm(u,NORM_2,&norm);
 73:   PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3dn",its);
 74:   PetscPrintf(PETSC_COMM_WORLD,"Residual norm %An",norm);
 75:   PetscPrintf(PETSC_COMM_WORLD,"Time for solve = %5.2f secondsn",time);

 77:   /* Cleanup */
 78:   SLESDestroy(sles);
 79:   VecDestroy(x);
 80:   VecDestroy(b);
 81:   VecDestroy(u);
 82:   MatDestroy(A);

 84:   PetscFinalize();
 85:   return 0;
 86: }

 88: #else
 89: #include <stdio.h>
 90: int main(int argc,char **args)
 91: {
 92:   fprintf(stdout,"This example does not work for complex numbers.n");
 93:   return 0;
 94: }
 95: #endif