Actual source code: petscsysmod.F

petsc-3.9.4 2018-09-11
Report Typos and Errors

  2:         module petscsysdefdummy
  3: #include <petscconf.h>
  4: #if defined(PETSC_HAVE_MPI_F90MODULE) || defined(PETSC_HAVE_MPIUNI)
  5:         use mpi
  6: #define PETSC_AVOID_MPIF_H
  7: #endif
  8: #include <../src/sys/f90-mod/petscsys.h>
  9: #include <../src/sys/f90-mod/petscdraw.h>
 10: #include <../src/sys/f90-mod/petscviewer.h>
 11:         end module

 13:         module petscsysdef
 14:         use petscsysdefdummy
 15:         interface operator(.ne.)
 16:           function petscviewernotequal(A,B)
 17:             use petscsysdefdummy
 18:             logical petscviewernotequal
 19:             type(tPetscViewer), intent(in) :: A,B
 20:           end function
 21:         end interface operator (.ne.)
 22:         interface operator(.eq.)
 23:           function petscviewerequals(A,B)
 24:             use petscsysdefdummy
 25:             logical petscviewerequals
 26:             type(tPetscViewer), intent(in) :: A,B
 27:           end function
 28:         end interface operator (.eq.)

 30:         interface operator(.ne.)
 31:         function petscrandomnotequal(A,B)
 32:           use petscsysdefdummy
 33:           logical petscrandomnotequal
 34:           type(tPetscRandom), intent(in) :: A,B
 35:         end function
 36:         end interface operator (.ne.)
 37:         interface operator(.eq.)
 38:         function petscrandomequals(A,B)
 39:           use petscsysdefdummy
 40:           logical petscrandomequals
 41:           type(tPetscRandom), intent(in) :: A,B
 42:         end function
 43:         end interface operator (.eq.)
 44:         end module

 46:         function petscviewernotequal(A,B)
 47:           use petscsysdefdummy
 48:           logical petscviewernotequal
 49:           type(tPetscViewer), intent(in) :: A,B
 50:           petscviewernotequal = (A%v .ne. B%v)
 51:         end function
 52:         function petscviewerequals(A,B)
 53:           use petscsysdefdummy
 54:           logical petscviewerequals
 55:           type(tPetscViewer), intent(in) :: A,B
 56:           petscviewerequals = (A%v .eq. B%v)
 57:         end function

 59:         function petscrandomnotequal(A,B)
 60:           use petscsysdefdummy
 61:           logical petscrandomnotequal
 62:           type(tPetscRandom), intent(in) :: A,B
 63:           petscrandomnotequal = (A%v .ne. B%v)
 64:         end function
 65:         function petscrandomequals(A,B)
 66:           use petscsysdefdummy
 67:           logical petscrandomequals
 68:           type(tPetscRandom), intent(in) :: A,B
 69:           petscrandomequals = (A%v .eq. B%v)
 70:         end function

 72:         module petscsys
 73:         use iso_c_binding
 74:         use petscsysdef
 75: #include <../src/sys/f90-mod/petscsys.h90>
 76:         interface
 77: #include <../src/sys/f90-mod/ftn-auto-interfaces/petscsys.h90>
 78:         end interface
 79:         end module