Actual source code: ex18.c

  1: /*$Id: ex18.c,v 1.27 2001/08/07 21:30:50 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.ex\n\n";

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

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

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

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

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

 56:   /* Set up solution */
 57:   VecDuplicate(b,&x);
 58:   VecDuplicate(b,&u);
 59:   VecSet(&zero,x);

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

 72:   /* Show result */
 73:   MatMult(A,x,u);
 74:   VecAXPY(&none,b,u);
 75:   VecNorm(u,NORM_2,&norm);
 76:   PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3d\n",its);
 77:   PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A\n",norm);
 78:   PetscPrintf(PETSC_COMM_WORLD,"Time for solve = %5.2f seconds\n",time);

 80:   /* Cleanup */
 81:   SLESDestroy(sles);
 82:   VecDestroy(x);
 83:   VecDestroy(b);
 84:   VecDestroy(u);
 85:   MatDestroy(A);

 87:   PetscFinalize();
 88:   return 0;
 89: }

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