Actual source code: petscviewer.h90

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1: interface
  2:   Subroutine PetscViewerDestroy(v,ierr)
  3:     use petscdef
  4:     PetscViewer v
  5:     PetscErrorCode, intent(out):: ierr
  6:   End Subroutine
  7:   Subroutine PetscViewerBinaryOpen(c,n,t,v,ierr)
  8:     use petscdef
  9:     MPI_Comm :: c
 10:     character(len=*), intent(in) :: n
 11:     PetscFileMode, intent(in) :: t
 12:     PetscViewer, intent(out) :: v
 13:     PetscErrorCode, intent(out):: ierr
 14:   End Subroutine
 15:   Subroutine PetscViewerGetSubViewer(v,c,vn,ierr)
 16:     use petscdef
 17:     PetscViewer, intent(in) :: v
 18:     MPI_Comm :: c
 19:     PetscFileMode, intent(out) :: vn
 20:     PetscErrorCode, intent(out):: ierr
 21:   End Subroutine
 22: end interface

 24: Interface PetscViewerBinaryWrite
 25:   Subroutine PetscViewerBinaryWriteInt(v,a,cnt,tmp,ierr)
 26:     use petscdef
 27:     PetscViewer v
 28:     PetscInt a(*)
 29:     PetscInt cnt
 30:     PetscBool  tmp
 31:     PetscErrorCode, intent(out):: ierr
 32:   End Subroutine

 34:   Subroutine PetscViewerBinaryWriteScalar(v,a,cnt,tmp,ierr)
 35:     use petscdef
 36:     PetscViewer v
 37:     PetscScalar a(*)
 38:     PetscInt cnt
 39:     PetscBool  tmp
 40:     PetscErrorCode, intent(out):: ierr
 41:   End Subroutine

 43: #if defined(PETSC_USE_COMPLEX)
 44:   Subroutine PetscViewerBinaryWriteReal(v,a,cnt,tmp,ierr)
 45:     use petscdef
 46:     PetscViewer v
 47:     PetscReal a(*)
 48:     PetscInt cnt
 49:     PetscBool  tmp
 50:     PetscErrorCode, intent(out):: ierr
 51:   End Subroutine
 52: #endif

 54:   Subroutine PetscViewerBinaryReadInt(v,a,cnt,ierr)
 55:     use petscdef
 56:     PetscViewer v
 57:     PetscInt a(*)
 58:     PetscInt cnt
 59:     PetscErrorCode, intent(out):: ierr
 60:   End Subroutine

 62:   Subroutine PetscViewerBinaryReadScalar(v,a,cnt,ierr)
 63:     use petscdef
 64:     PetscViewer v
 65:     PetscScalar a(*)
 66:     PetscInt cnt
 67:     PetscErrorCode, intent(out):: ierr
 68:   End Subroutine

 70: #if defined(PETSC_USE_COMPLEX)
 71:   Subroutine PetscViewerBinaryReadReal(v,a,cnt,ierr)
 72:     use petscdef
 73:     PetscViewer v
 74:     PetscReal a(*)
 75:     PetscInt cnt
 76:     PetscErrorCode, intent(out):: ierr
 77:   End Subroutine
 78: #endif

 80: End Interface