Actual source code: ex41.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: static char help[] = "Reads a PETSc matrix and vector from a socket connection,  solves a linear system and sends the result back.\n";

  4: /*T
  5:    Concepts: KSP^solving a linear system
  6:    Processors: n
  7: T*/

  9: /*
 10:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 11:   automatically includes:
 12:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 13:      petscmat.h - matrices
 14:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 15:      petscviewer.h - viewers               petscpc.h  - preconditioners
 16: */
 17:  #include <petscksp.h>

 19: int main(int argc,char **args)
 20: {
 21:   KSP            ksp;             /* linear solver context */
 22:   Mat            A;            /* matrix */
 23:   Vec            x,b;          /* approx solution, RHS, exact solution */
 24:   PetscViewer    fd;               /* viewer */

 27:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 28:   fd = PETSC_VIEWER_SOCKET_WORLD;

 30:   VecCreate(PETSC_COMM_WORLD,&b);
 31:   VecLoad(b,fd);
 32:   MatCreate(PETSC_COMM_WORLD,&A);
 33:   MatLoad(A,fd);
 34:   VecDuplicate(b,&x);

 36:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 37:   KSPSetOperators(ksp,A,A);
 38:   KSPSetFromOptions(ksp);
 39:   KSPSetUp(ksp);
 40:   KSPSolve(ksp,b,x);
 41:   VecView(x,fd);
 42:   MatDestroy(&A);
 43:   VecDestroy(&b);
 44:   VecDestroy(&x);
 45:   KSPDestroy(&ksp);

 47:   PetscFinalize();
 48:   return ierr;
 49: }