Actual source code: petscmatmod.F90

  1: module petscmatdef
  2:   use, intrinsic :: ISO_C_binding
  3:   use petscvecdef
  4: #include "petsc/finclude/petscmat.h"
  5: #include "petsc/finclude/petscmatcoarsen.h"
  6: #include "petsc/finclude/petscpartitioner.h"
  7: #include "petsc/finclude/petscmathypre.h"
  8: #include "petsc/finclude/petscmathtool.h"
  9: #include "petsc/finclude/petscmatelemental.h"
 10: #include <../ftn/mat/petscmat.h>
 11: #include <../ftn/mat/petscmatcoarsen.h>
 12: #include <../ftn/mat/petscpartitioner.h>

 14: end module

 16: module petscmat
 17:   use petscmatdef
 18:   use petscvec

 20: #include <../src/mat/ftn-mod/petscmat.h90>
 21: #include <../ftn/mat/petscmat.h90>
 22: #include <../ftn/mat/petscmatcoarsen.h90>
 23: #include <../ftn/mat/petscpartitioner.h90>

 25: !     deprecated functions

 27:   interface MatDenseGetArrayF90
 28:     module procedure MatDenseGetArrayF901d, MatDenseGetArrayF902d
 29:   end interface

 31:   interface MatDenseRestoreArrayF90
 32:     module procedure MatDenseRestoreArrayF901d, MatDenseRestoreArrayF902d
 33:   end interface

 35:   interface MatDenseGetArrayReadF90
 36:     module procedure MatDenseGetArrayReadF901d, MatDenseGetArrayReadF902d
 37:   end interface

 39:   interface MatDenseRestoreArrayReadF90
 40:     module procedure MatDenseRestoreArrayWriteF901d, MatDenseRestoreArrayWriteF902d
 41:   end interface

 43:   interface MatDenseGetArrayWriteF90
 44:     module procedure MatDenseGetArrayWriteF901d, MatDenseGetArrayWriteF902d
 45:   end interface

 47:   interface MatDenseRestoreArrayWriteF90
 48:     module procedure MatDenseRestoreArrayWriteF901d, MatDenseRestoreArrayWriteF902d
 49:   end interface

 51: contains

 53: #include <../ftn/mat/petscmat.hf90>
 54: #include <../ftn/mat/petscmatcoarsen.hf90>
 55: #include <../ftn/mat/petscpartitioner.hf90>

 57: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
 58: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseGetArrayF901d
 59: #endif
 60:   subroutine MatDenseGetArrayF901d(v, array, ierr)
 61:     PetscScalar, pointer :: array(:)
 62:     PetscErrorCode ierr
 63:     Mat v
 64:     call MatDenseGetArray(v, array, ierr)
 65:   end subroutine

 67: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
 68: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseRestoreArrayF901d
 69: #endif
 70:   subroutine MatDenseRestoreArrayF901d(v, array, ierr)
 71:     PetscScalar, pointer :: array(:)
 72:     PetscErrorCode ierr
 73:     Mat v
 74:     call MatDenseRestoreArray(v, array, ierr)
 75:   end subroutine

 77: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
 78: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseGetArrayReadF901d
 79: #endif
 80:   subroutine MatDenseGetArrayReadF901d(v, array, ierr)
 81:     PetscScalar, pointer :: array(:)
 82:     PetscErrorCode ierr
 83:     Mat v
 84:     call MatDenseGetArrayRead(v, array, ierr)
 85:   end subroutine

 87: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
 88: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseRestoreArrayReadF901d
 89: #endif
 90:   subroutine MatDenseRestoreArrayReadF901d(v, array, ierr)
 91:     PetscScalar, pointer :: array(:)
 92:     PetscErrorCode ierr
 93:     Mat v
 94:     call MatDenseRestoreArrayRead(v, array, ierr)
 95:   end subroutine

 97: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
 98: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseGetArrayWriteF901d
 99: #endif
100:   subroutine MatDenseGetArrayWriteF901d(v, array, ierr)
101:     PetscScalar, pointer :: array(:)
102:     PetscErrorCode ierr
103:     Mat v
104:     call MatDenseGetArrayWrite(v, array, ierr)
105:   end subroutine

107: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
108: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseRestoreArrayWriteF901d
109: #endif
110:   subroutine MatDenseRestoreArrayWriteF901d(v, array, ierr)
111:     PetscScalar, pointer :: array(:)
112:     PetscErrorCode ierr
113:     Mat v
114:     call MatDenseRestoreArrayWrite(v, array, ierr)
115:   end subroutine

117: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
118: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseGetArrayF902d
119: #endif
120:   subroutine MatDenseGetArrayF902d(v, array, ierr)
121:     PetscScalar, pointer :: array(:, :)
122:     PetscErrorCode ierr
123:     Mat v
124:     call MatDenseGetArray(v, array, ierr)
125:   end subroutine

127: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
128: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseRestoreArrayF902d
129: #endif
130:   subroutine MatDenseRestoreArrayF902d(v, array, ierr)
131:     PetscScalar, pointer :: array(:, :)
132:     PetscErrorCode ierr
133:     Mat v
134:     call MatDenseRestoreArray(v, array, ierr)
135:   end subroutine

137: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
138: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseGetArrayReadF902d
139: #endif
140:   subroutine MatDenseGetArrayReadF902d(v, array, ierr)
141:     PetscScalar, pointer :: array(:, :)
142:     PetscErrorCode ierr
143:     Mat v
144:     call MatDenseGetArrayRead(v, array, ierr)
145:   end subroutine

147: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
148: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseRestoreArrayReadF90
149: #endif
150:   subroutine MatDenseRestoreArrayReadF902d(v, array, ierr)
151:     PetscScalar, pointer :: array(:, :)
152:     PetscErrorCode ierr
153:     Mat v
154:     call MatDenseRestoreArrayRead(v, array, ierr)
155:   end subroutine

157: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
158: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseGetArrayWriteF90
159: #endif
160:   subroutine MatDenseGetArrayWriteF902d(v, array, ierr)
161:     PetscScalar, pointer :: array(:, :)
162:     PetscErrorCode ierr
163:     Mat v
164:     call MatDenseGetArrayWrite(v, array, ierr)
165:   end subroutine

167: #if defined(_WIN32) && defined(PETSC_USE_SHARED_LIBRARIES)
168: !DEC$ ATTRIBUTES DLLEXPORT:: MatDenseRestoreArrayWriteF90
169: #endif
170:   subroutine MatDenseRestoreArrayWriteF902d(v, array, ierr)
171:     PetscScalar, pointer :: array(:, :)
172:     PetscErrorCode ierr
173:     Mat v
174:     call MatDenseRestoreArrayWrite(v, array, ierr)
175:   end subroutine

177: end module