Actual source code: ex6.c

  1: /*$Id: ex6.c,v 1.71 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: Input arguments are:n
  5:   -f <input_file> : file to load.  For a 5X5 example of the 5-pt. stencil,n
  6:                     use the file petsc/src/mat/examples/matbinary.exnn";

 8:  #include petscsles.h
 9:  #include petsclog.h

 11: int main(int argc,char **args)
 12: {
 13:   int            ierr,its;
 14:   double         norm;
 15:   PetscLogDouble tsetup1,tsetup2,tsetup,tsolve1,tsolve2,tsolve;
 16:   Scalar         zero = 0.0,none = -1.0;
 17:   Vec            x,b,u;
 18:   Mat            A;
 19:   SLES           sles;
 20:   char           file[128];
 21:   PetscViewer    fd;
 22:   PetscTruth     table,flg;

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

 26:   PetscOptionsHasName(PETSC_NULL,"-table",&table);

 28: #if defined(PETSC_USE_COMPLEX)
 29:   SETERRQ(1,"This example does not work with complex numbers");
 30: #else

 32:   /* Read matrix and RHS */
 33:   PetscOptionsGetString(PETSC_NULL,"-f",file,127,&flg);
 34:   if (!flg) SETERRQ(1,"Must indicate binary file with the -f option");
 35:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_BINARY_RDONLY,&fd);

 37:   MatLoad(fd,MATSEQAIJ,&A);
 38:   VecLoad(fd,&b);
 39:   PetscViewerDestroy(fd);

 41:   /* 
 42:    If the load matrix is larger then the vector, due to being padded 
 43:    to match the blocksize then create a new padded vector
 44:   */
 45:   {
 46:     int    m,n,j,mvec,start,end,index;
 47:     Vec    tmp;
 48:     Scalar *bold;

 50:     MatGetLocalSize(A,&m,&n);
 51:     VecCreateMPI(PETSC_COMM_WORLD,m,PETSC_DECIDE,&tmp);
 52:     VecGetOwnershipRange(b,&start,&end);
 53:     VecGetLocalSize(b,&mvec);
 54:     VecGetArray(b,&bold);
 55:     for (j=0; j<mvec; j++) {
 56:       index = start+j;
 57:       ierr  = VecSetValues(tmp,1,&index,bold+j,INSERT_VALUES);
 58:     }
 59:     VecRestoreArray(b,&bold);
 60:     VecDestroy(b);
 61:     VecAssemblyBegin(tmp);
 62:     VecAssemblyEnd(tmp);
 63:     b = tmp;
 64:   }
 65:   VecDuplicate(b,&x);
 66:   VecDuplicate(b,&u);

 68:   VecSet(&zero,x);
 69:   PetscBarrier((PetscObject)A);

 71:   PetscLogStagePush(1);
 72:   PetscGetTime(&tsetup1);
 73:   SLESCreate(PETSC_COMM_WORLD,&sles);
 74:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
 75:   SLESSetFromOptions(sles);
 76:   SLESSetUp(sles,b,x);
 77:   SLESSetUpOnBlocks(sles);
 78:   PetscGetTime(&tsetup2);
 79:   tsetup = tsetup2 -tsetup1;
 80:   PetscLogStagePop();
 81:   PetscBarrier((PetscObject)A);

 83:   PetscLogStagePush(2);
 84:   PetscGetTime(&tsolve1);
 85:   SLESSolve(sles,b,x,&its);
 86:   PetscGetTime(&tsolve2);
 87:   tsolve = tsolve2 - tsolve1;
 88:   PetscLogStagePop();

 90:   /* Show result */
 91:   MatMult(A,x,u);
 92:   VecAXPY(&none,b,u);
 93:   VecNorm(u,NORM_2,&norm);
 94:   /*  matrix PC   KSP   Options       its    residual setuptime solvetime  */
 95:   if (table) {
 96:     char   *matrixname,slesinfo[120];
 97:     PetscViewer viewer;
 98:     PetscViewerStringOpen(PETSC_COMM_WORLD,slesinfo,120,&viewer);
 99:     SLESView(sles,viewer);
100:     PetscStrrchr(file,'/',&matrixname);
101:     PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3d %2.0e %2.1e %2.1e %2.1e %s n",
102:                        matrixname,its,norm,tsetup+tsolve,tsetup,tsolve,slesinfo);
103:     PetscViewerDestroy(viewer);
104:   } else {
105:     PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3dn",its);
106:     PetscPrintf(PETSC_COMM_WORLD,"Residual norm = %An",norm);
107:   }

109:   /* Cleanup */
110:   SLESDestroy(sles);
111:   VecDestroy(x);
112:   VecDestroy(b);
113:   VecDestroy(u);
114:   MatDestroy(A);

116:   PetscFinalize();
117: #endif
118:   return 0;
119: }