xref: /petsc/src/sys/ftn-src/f90_fwrap.F90 (revision fe66ebcc023cb303106674d426ee542bea707d38)
16dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
26dd63270SBarry Smith!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
36dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
46dd63270SBarry Smith#include <petsc/finclude/petscsys.h>
56dd63270SBarry Smith      subroutine F90Array1dCreateScalar(array,start,len1,ptr)
6*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
76dd63270SBarry Smith      implicit none
86dd63270SBarry Smith      PetscInt start,len1
96dd63270SBarry Smith      PetscScalar, target ::                                                      &
106dd63270SBarry Smith     &     array(start:start+len1-1)
116dd63270SBarry Smith      PetscScalar, pointer :: ptr(:)
126dd63270SBarry Smith
136dd63270SBarry Smith      ptr => array
146dd63270SBarry Smith      end subroutine
156dd63270SBarry Smith
166dd63270SBarry Smith      subroutine F90Array1dCreateReal(array,start,len1,ptr)
17*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
186dd63270SBarry Smith      implicit none
196dd63270SBarry Smith      PetscInt start,len1
206dd63270SBarry Smith      PetscReal, target ::                                                        &
216dd63270SBarry Smith     &     array(start:start+len1-1)
226dd63270SBarry Smith      PetscReal, pointer :: ptr(:)
236dd63270SBarry Smith
246dd63270SBarry Smith      ptr => array
256dd63270SBarry Smith      end subroutine
266dd63270SBarry Smith
276dd63270SBarry Smith      subroutine F90Array1dCreateInt(array,start,len1,ptr)
28*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
296dd63270SBarry Smith      implicit none
306dd63270SBarry Smith      PetscInt start,len1
316dd63270SBarry Smith      PetscInt, target ::                                                         &
326dd63270SBarry Smith     &     array(start:start+len1-1)
336dd63270SBarry Smith      PetscInt, pointer :: ptr(:)
346dd63270SBarry Smith
356dd63270SBarry Smith      ptr => array
366dd63270SBarry Smith      end subroutine
376dd63270SBarry Smith
386dd63270SBarry Smith      subroutine F90Array1dCreateMPIInt(array,start,len1,ptr)
39*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
406dd63270SBarry Smith      implicit none
416dd63270SBarry Smith      PetscInt start,len1
426dd63270SBarry Smith      PetscMPIInt, target ::                                                      &
436dd63270SBarry Smith      &     array(start:start+len1-1)
446dd63270SBarry Smith      PetscMPIInt, pointer :: ptr(:)
456dd63270SBarry Smith
466dd63270SBarry Smith      ptr => array
476dd63270SBarry Smith      end subroutine
486dd63270SBarry Smith
496dd63270SBarry Smith      subroutine F90Array1dCreateFortranAddr(array,start,len1,ptr)
50*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
516dd63270SBarry Smith      implicit none
526dd63270SBarry Smith      PetscInt start,len1
536dd63270SBarry Smith      PetscFortranAddr, target ::                                                 &
546dd63270SBarry Smith     &     array(start:start+len1-1)
556dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:)
566dd63270SBarry Smith
576dd63270SBarry Smith      ptr => array
586dd63270SBarry Smith      end subroutine
596dd63270SBarry Smith
606dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
616dd63270SBarry Smith      subroutine F90Array1dAccessScalar(ptr,address)
62*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
636dd63270SBarry Smith      implicit none
646dd63270SBarry Smith      PetscScalar, pointer :: ptr(:)
656dd63270SBarry Smith      PetscFortranAddr address
666dd63270SBarry Smith      PetscInt start
676dd63270SBarry Smith
686dd63270SBarry Smith      if (associated(ptr) .eqv. .false.) then
696dd63270SBarry Smith        address = 0
706dd63270SBarry Smith      else
716dd63270SBarry Smith        start = lbound(ptr,1)
726dd63270SBarry Smith        call F90Array1dGetAddrScalar(ptr(start),address)
736dd63270SBarry Smith      endif
746dd63270SBarry Smith      end subroutine
756dd63270SBarry Smith
766dd63270SBarry Smith      subroutine F90Array1dAccessReal(ptr,address)
77*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
786dd63270SBarry Smith      implicit none
796dd63270SBarry Smith      PetscReal, pointer :: ptr(:)
806dd63270SBarry Smith      PetscFortranAddr address
816dd63270SBarry Smith      PetscInt start
826dd63270SBarry Smith
836dd63270SBarry Smith      if (associated(ptr) .eqv. .false.) then
846dd63270SBarry Smith        address = 0
856dd63270SBarry Smith      else
866dd63270SBarry Smith        start = lbound(ptr,1)
876dd63270SBarry Smith        call F90Array1dGetAddrReal(ptr(start),address)
886dd63270SBarry Smith      endif
896dd63270SBarry Smith      end subroutine
906dd63270SBarry Smith
916dd63270SBarry Smith      subroutine F90Array1dAccessInt(ptr,address)
92*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
936dd63270SBarry Smith      implicit none
946dd63270SBarry Smith      PetscInt, pointer :: ptr(:)
956dd63270SBarry Smith      PetscFortranAddr address
966dd63270SBarry Smith      PetscInt start
976dd63270SBarry Smith
986dd63270SBarry Smith      if (associated(ptr) .eqv. .false.) then
996dd63270SBarry Smith        address = 0
1006dd63270SBarry Smith      else
1016dd63270SBarry Smith        start = lbound(ptr,1)
1026dd63270SBarry Smith        call F90Array1dGetAddrInt(ptr(start),address)
1036dd63270SBarry Smith      endif
1046dd63270SBarry Smith      end subroutine
1056dd63270SBarry Smith
1066dd63270SBarry Smith      subroutine F90Array1dAccessMPIInt(ptr,address)
107*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1086dd63270SBarry Smith        implicit none
1096dd63270SBarry Smith        PetscMPIInt, pointer :: ptr(:)
1106dd63270SBarry Smith        PetscFortranAddr address
1116dd63270SBarry Smith        PetscInt start
1126dd63270SBarry Smith
1136dd63270SBarry Smith        if (associated(ptr) .eqv. .false.) then
1146dd63270SBarry Smith          address = 0
1156dd63270SBarry Smith        else
1166dd63270SBarry Smith          start = lbound(ptr,1)
1176dd63270SBarry Smith          call F90Array1dGetAddrMPIInt(ptr(start),address)
1186dd63270SBarry Smith        endif
1196dd63270SBarry Smith        end subroutine
1206dd63270SBarry Smith
1216dd63270SBarry Smith      subroutine F90Array1dAccessFortranAddr(ptr,address)
122*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1236dd63270SBarry Smith      implicit none
1246dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:)
1256dd63270SBarry Smith      PetscFortranAddr address
1266dd63270SBarry Smith      PetscInt start
1276dd63270SBarry Smith
1286dd63270SBarry Smith      if (associated(ptr) .eqv. .false.) then
1296dd63270SBarry Smith        address = 0
1306dd63270SBarry Smith      else
1316dd63270SBarry Smith        start = lbound(ptr,1)
1326dd63270SBarry Smith        call F90Array1dGetAddrFortranAddr(ptr(start),address)
1336dd63270SBarry Smith      endif
1346dd63270SBarry Smith      end subroutine
1356dd63270SBarry Smith
1366dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1376dd63270SBarry Smith      subroutine F90Array1dDestroyScalar(ptr)
138*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1396dd63270SBarry Smith      implicit none
1406dd63270SBarry Smith      PetscScalar, pointer :: ptr(:)
1416dd63270SBarry Smith
1426dd63270SBarry Smith      nullify(ptr)
1436dd63270SBarry Smith      end subroutine
1446dd63270SBarry Smith
1456dd63270SBarry Smith      subroutine F90Array1dDestroyReal(ptr)
146*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1476dd63270SBarry Smith      implicit none
1486dd63270SBarry Smith      PetscReal, pointer :: ptr(:)
1496dd63270SBarry Smith
1506dd63270SBarry Smith      nullify(ptr)
1516dd63270SBarry Smith      end subroutine
1526dd63270SBarry Smith
1536dd63270SBarry Smith      subroutine F90Array1dDestroyInt(ptr)
154*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1556dd63270SBarry Smith      implicit none
1566dd63270SBarry Smith      PetscInt, pointer :: ptr(:)
1576dd63270SBarry Smith
1586dd63270SBarry Smith      nullify(ptr)
1596dd63270SBarry Smith      end subroutine
1606dd63270SBarry Smith
1616dd63270SBarry Smith      subroutine F90Array1dDestroyMPIInt(ptr)
162*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1636dd63270SBarry Smith      implicit none
1646dd63270SBarry Smith      PetscMPIInt, pointer :: ptr(:)
1656dd63270SBarry Smith
1666dd63270SBarry Smith      nullify(ptr)
1676dd63270SBarry Smith      end subroutine
1686dd63270SBarry Smith
1696dd63270SBarry Smith      subroutine F90Array1dDestroyFortranAddr(ptr)
170*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1716dd63270SBarry Smith      implicit none
1726dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:)
1736dd63270SBarry Smith
1746dd63270SBarry Smith      nullify(ptr)
1756dd63270SBarry Smith      end subroutine
1766dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1776dd63270SBarry Smith!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
1786dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1796dd63270SBarry Smith      subroutine F90Array2dCreateScalar(array,start1,len1,                        &
1806dd63270SBarry Smith     &     start2,len2,ptr)
181*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1826dd63270SBarry Smith      implicit none
1836dd63270SBarry Smith      PetscInt start1,len1
1846dd63270SBarry Smith      PetscInt start2,len2
1856dd63270SBarry Smith      PetscScalar, target ::                                                      &
1866dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1)
1876dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:)
1886dd63270SBarry Smith
1896dd63270SBarry Smith      ptr => array
1906dd63270SBarry Smith      end subroutine
1916dd63270SBarry Smith
1926dd63270SBarry Smith      subroutine F90Array2dCreateReal(array,start1,len1,                          &
1936dd63270SBarry Smith     &     start2,len2,ptr)
194*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
1956dd63270SBarry Smith      implicit none
1966dd63270SBarry Smith      PetscInt start1,len1
1976dd63270SBarry Smith      PetscInt start2,len2
1986dd63270SBarry Smith      PetscReal, target ::                                                        &
1996dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1)
2006dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:)
2016dd63270SBarry Smith
2026dd63270SBarry Smith      ptr => array
2036dd63270SBarry Smith      end subroutine
2046dd63270SBarry Smith
2056dd63270SBarry Smith      subroutine F90Array2dCreateInt(array,start1,len1,                           &
2066dd63270SBarry Smith     &     start2,len2,ptr)
207*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2086dd63270SBarry Smith      implicit none
2096dd63270SBarry Smith      PetscInt start1,len1
2106dd63270SBarry Smith      PetscInt start2,len2
2116dd63270SBarry Smith      PetscInt, target ::                                                         &
2126dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1)
2136dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:)
2146dd63270SBarry Smith
2156dd63270SBarry Smith      ptr => array
2166dd63270SBarry Smith      end subroutine
2176dd63270SBarry Smith
2186dd63270SBarry Smith      subroutine F90Array2dCreateFortranAddr(array,start1,len1,                   &
2196dd63270SBarry Smith     &     start2,len2,ptr)
220*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2216dd63270SBarry Smith      implicit none
2226dd63270SBarry Smith      PetscInt start1,len1
2236dd63270SBarry Smith      PetscInt start2,len2
2246dd63270SBarry Smith      PetscFortranAddr, target ::                                                 &
2256dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1)
2266dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:)
2276dd63270SBarry Smith
2286dd63270SBarry Smith      ptr => array
2296dd63270SBarry Smith      end subroutine
2306dd63270SBarry Smith
2316dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2326dd63270SBarry Smith      subroutine F90Array2dAccessScalar(ptr,address)
233*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2346dd63270SBarry Smith      implicit none
2356dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:)
2366dd63270SBarry Smith      PetscFortranAddr address
2376dd63270SBarry Smith      PetscInt start1,start2
2386dd63270SBarry Smith
2396dd63270SBarry Smith      start1 = lbound(ptr,1)
2406dd63270SBarry Smith      start2 = lbound(ptr,2)
2416dd63270SBarry Smith      call F90Array2dGetAddrScalar(ptr(start1,start2),address)
2426dd63270SBarry Smith      end subroutine
2436dd63270SBarry Smith
2446dd63270SBarry Smith      subroutine F90Array2dAccessReal(ptr,address)
245*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2466dd63270SBarry Smith      implicit none
2476dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:)
2486dd63270SBarry Smith      PetscFortranAddr address
2496dd63270SBarry Smith      PetscInt start1,start2
2506dd63270SBarry Smith
2516dd63270SBarry Smith      start1 = lbound(ptr,1)
2526dd63270SBarry Smith      start2 = lbound(ptr,2)
2536dd63270SBarry Smith      call F90Array2dGetAddrReal(ptr(start1,start2),address)
2546dd63270SBarry Smith      end subroutine
2556dd63270SBarry Smith
2566dd63270SBarry Smith      subroutine F90Array2dAccessInt(ptr,address)
257*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2586dd63270SBarry Smith      implicit none
2596dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:)
2606dd63270SBarry Smith      PetscFortranAddr address
2616dd63270SBarry Smith      PetscInt start1,start2
2626dd63270SBarry Smith
2636dd63270SBarry Smith      start1 = lbound(ptr,1)
2646dd63270SBarry Smith      start2 = lbound(ptr,2)
2656dd63270SBarry Smith      call F90Array2dGetAddrInt(ptr(start1,start2),address)
2666dd63270SBarry Smith      end subroutine
2676dd63270SBarry Smith
2686dd63270SBarry Smith      subroutine F90Array2dAccessFortranAddr(ptr,address)
269*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2706dd63270SBarry Smith      implicit none
2716dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:)
2726dd63270SBarry Smith      PetscFortranAddr address
2736dd63270SBarry Smith      PetscInt start1,start2
2746dd63270SBarry Smith
2756dd63270SBarry Smith      start1 = lbound(ptr,1)
2766dd63270SBarry Smith      start2 = lbound(ptr,2)
2776dd63270SBarry Smith      call F90Array2dGetAddrFortranAddr(ptr(start1,start2),address)
2786dd63270SBarry Smith      end subroutine
2796dd63270SBarry Smith
2806dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2816dd63270SBarry Smith      subroutine F90Array2dDestroyScalar(ptr)
282*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2836dd63270SBarry Smith      implicit none
2846dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:)
2856dd63270SBarry Smith
2866dd63270SBarry Smith      nullify(ptr)
2876dd63270SBarry Smith      end subroutine
2886dd63270SBarry Smith
2896dd63270SBarry Smith      subroutine F90Array2dDestroyReal(ptr)
290*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2916dd63270SBarry Smith      implicit none
2926dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:)
2936dd63270SBarry Smith
2946dd63270SBarry Smith      nullify(ptr)
2956dd63270SBarry Smith      end subroutine
2966dd63270SBarry Smith
2976dd63270SBarry Smith      subroutine F90Array2dDestroyInt(ptr)
298*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
2996dd63270SBarry Smith      implicit none
3006dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:)
3016dd63270SBarry Smith
3026dd63270SBarry Smith      nullify(ptr)
3036dd63270SBarry Smith      end subroutine
3046dd63270SBarry Smith
3056dd63270SBarry Smith      subroutine F90Array2dDestroyFortranAddr(ptr)
306*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
3076dd63270SBarry Smith      implicit none
3086dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:)
3096dd63270SBarry Smith
3106dd63270SBarry Smith      nullify(ptr)
3116dd63270SBarry Smith      end subroutine
3126dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3136dd63270SBarry Smith!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
3146dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3156dd63270SBarry Smith      subroutine F90Array3dCreateScalar(array,start1,len1,                        &
3166dd63270SBarry Smith     &     start2,len2,start3,len3,ptr)
317*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
3186dd63270SBarry Smith      implicit none
3196dd63270SBarry Smith      PetscInt start1,len1
3206dd63270SBarry Smith      PetscInt start2,len2
3216dd63270SBarry Smith      PetscInt start3,len3
3226dd63270SBarry Smith      PetscScalar, target ::                                                      &
3236dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
3246dd63270SBarry Smith     &           start3:start3+len3-1)
3256dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:,:)
3266dd63270SBarry Smith
3276dd63270SBarry Smith      ptr => array
3286dd63270SBarry Smith      end subroutine
3296dd63270SBarry Smith
3306dd63270SBarry Smith      subroutine F90Array3dCreateReal(array,start1,len1,                          &
3316dd63270SBarry Smith     &     start2,len2,start3,len3,ptr)
332*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
3336dd63270SBarry Smith      implicit none
3346dd63270SBarry Smith      PetscInt start1,len1
3356dd63270SBarry Smith      PetscInt start2,len2
3366dd63270SBarry Smith      PetscInt start3,len3
3376dd63270SBarry Smith      PetscReal, target ::                                                        &
3386dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
3396dd63270SBarry Smith     &           start3:start3+len3-1)
3406dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:,:)
3416dd63270SBarry Smith
3426dd63270SBarry Smith      ptr => array
3436dd63270SBarry Smith      end subroutine
3446dd63270SBarry Smith
3456dd63270SBarry Smith      subroutine F90Array3dCreateInt(array,start1,len1,                           &
3466dd63270SBarry Smith     &     start2,len2,start3,len3,ptr)
347*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
3486dd63270SBarry Smith      implicit none
3496dd63270SBarry Smith      PetscInt start1,len1
3506dd63270SBarry Smith      PetscInt start2,len2
3516dd63270SBarry Smith      PetscInt start3,len3
3526dd63270SBarry Smith      PetscInt, target ::                                                         &
3536dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
3546dd63270SBarry Smith     &           start3:start3+len3-1)
3556dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:,:)
3566dd63270SBarry Smith
3576dd63270SBarry Smith      ptr => array
3586dd63270SBarry Smith      end subroutine
3596dd63270SBarry Smith
3606dd63270SBarry Smith      subroutine F90Array3dCreateFortranAddr(array,start1,len1,                   &
3616dd63270SBarry Smith     &     start2,len2,start3,len3,ptr)
362*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
3636dd63270SBarry Smith      implicit none
3646dd63270SBarry Smith      PetscInt start1,len1
3656dd63270SBarry Smith      PetscInt start2,len2
3666dd63270SBarry Smith      PetscInt start3,len3
3676dd63270SBarry Smith      PetscFortranAddr, target ::                                                 &
3686dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
3696dd63270SBarry Smith     &           start3:start3+len3-1)
3706dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:,:)
3716dd63270SBarry Smith
3726dd63270SBarry Smith      ptr => array
3736dd63270SBarry Smith      end subroutine
3746dd63270SBarry Smith
3756dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3766dd63270SBarry Smith      subroutine F90Array3dAccessScalar(ptr,address)
377*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
3786dd63270SBarry Smith      implicit none
3796dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:,:)
3806dd63270SBarry Smith      PetscFortranAddr address
3816dd63270SBarry Smith      PetscInt start1,start2,start3
3826dd63270SBarry Smith
3836dd63270SBarry Smith      start1 = lbound(ptr,1)
3846dd63270SBarry Smith      start2 = lbound(ptr,2)
3856dd63270SBarry Smith      start3 = lbound(ptr,3)
3866dd63270SBarry Smith      call F90Array3dGetAddrScalar(ptr(start1,start2,start3),address)
3876dd63270SBarry Smith      end subroutine
3886dd63270SBarry Smith
3896dd63270SBarry Smith      subroutine F90Array3dAccessReal(ptr,address)
390*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
3916dd63270SBarry Smith      implicit none
3926dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:,:)
3936dd63270SBarry Smith      PetscFortranAddr address
3946dd63270SBarry Smith      PetscInt start1,start2,start3
3956dd63270SBarry Smith
3966dd63270SBarry Smith      start1 = lbound(ptr,1)
3976dd63270SBarry Smith      start2 = lbound(ptr,2)
3986dd63270SBarry Smith      start3 = lbound(ptr,3)
3996dd63270SBarry Smith      call F90Array3dGetAddrReal(ptr(start1,start2,start3),address)
4006dd63270SBarry Smith      end subroutine
4016dd63270SBarry Smith
4026dd63270SBarry Smith      subroutine F90Array3dAccessInt(ptr,address)
403*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4046dd63270SBarry Smith      implicit none
4056dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:,:)
4066dd63270SBarry Smith      PetscFortranAddr address
4076dd63270SBarry Smith      PetscInt start1,start2,start3
4086dd63270SBarry Smith
4096dd63270SBarry Smith      start1 = lbound(ptr,1)
4106dd63270SBarry Smith      start2 = lbound(ptr,2)
4116dd63270SBarry Smith      start3 = lbound(ptr,3)
4126dd63270SBarry Smith      call F90Array3dGetAddrInt(ptr(start1,start2,start3),address)
4136dd63270SBarry Smith      end subroutine
4146dd63270SBarry Smith
4156dd63270SBarry Smith      subroutine F90Array3dAccessFortranAddr(ptr,address)
416*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4176dd63270SBarry Smith      implicit none
4186dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:,:)
4196dd63270SBarry Smith      PetscFortranAddr address
4206dd63270SBarry Smith      PetscInt start1,start2,start3
4216dd63270SBarry Smith
4226dd63270SBarry Smith      start1 = lbound(ptr,1)
4236dd63270SBarry Smith      start2 = lbound(ptr,2)
4246dd63270SBarry Smith      start3 = lbound(ptr,3)
4256dd63270SBarry Smith      call F90Array3dGetAddrFortranAddr(ptr(start1,start2,start3),        &
4266dd63270SBarry Smith     &                                  address)
4276dd63270SBarry Smith      end subroutine
4286dd63270SBarry Smith
4296dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4306dd63270SBarry Smith      subroutine F90Array3dDestroyScalar(ptr)
431*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4326dd63270SBarry Smith      implicit none
4336dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:,:)
4346dd63270SBarry Smith
4356dd63270SBarry Smith      nullify(ptr)
4366dd63270SBarry Smith      end subroutine
4376dd63270SBarry Smith
4386dd63270SBarry Smith      subroutine F90Array3dDestroyReal(ptr)
439*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4406dd63270SBarry Smith      implicit none
4416dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:,:)
4426dd63270SBarry Smith
4436dd63270SBarry Smith      nullify(ptr)
4446dd63270SBarry Smith      end subroutine
4456dd63270SBarry Smith
4466dd63270SBarry Smith      subroutine F90Array3dDestroyInt(ptr)
447*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4486dd63270SBarry Smith      implicit none
4496dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:,:)
4506dd63270SBarry Smith
4516dd63270SBarry Smith      nullify(ptr)
4526dd63270SBarry Smith      end subroutine
4536dd63270SBarry Smith
4546dd63270SBarry Smith      subroutine F90Array3dDestroyFortranAddr(ptr)
455*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4566dd63270SBarry Smith      implicit none
4576dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:,:)
4586dd63270SBarry Smith
4596dd63270SBarry Smith      nullify(ptr)
4606dd63270SBarry Smith      end subroutine
4616dd63270SBarry Smith
4626dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4636dd63270SBarry Smith      subroutine F90Array4dCreateScalar(array,start1,len1,                        &
4646dd63270SBarry Smith     &     start2,len2,start3,len3,start4,len4,ptr)
465*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4666dd63270SBarry Smith      implicit none
4676dd63270SBarry Smith      PetscInt start1,len1
4686dd63270SBarry Smith      PetscInt start2,len2
4696dd63270SBarry Smith      PetscInt start3,len3
4706dd63270SBarry Smith      PetscInt start4,len4
4716dd63270SBarry Smith      PetscScalar, target ::                                                      &
4726dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
4736dd63270SBarry Smith     &           start3:start3+len3-1,start4:start4+len4-1)
4746dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:,:,:)
4756dd63270SBarry Smith
4766dd63270SBarry Smith      ptr => array
4776dd63270SBarry Smith      end subroutine
4786dd63270SBarry Smith
4796dd63270SBarry Smith      subroutine F90Array4dCreateReal(array,start1,len1,                          &
4806dd63270SBarry Smith     &     start2,len2,start3,len3,start4,len4,ptr)
481*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4826dd63270SBarry Smith      implicit none
4836dd63270SBarry Smith      PetscInt start1,len1
4846dd63270SBarry Smith      PetscInt start2,len2
4856dd63270SBarry Smith      PetscInt start3,len3
4866dd63270SBarry Smith      PetscInt start4,len4
4876dd63270SBarry Smith      PetscReal, target ::                                                        &
4886dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
4896dd63270SBarry Smith     &           start3:start3+len3-1,start4:start4+len4-1)
4906dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:,:,:)
4916dd63270SBarry Smith
4926dd63270SBarry Smith      ptr => array
4936dd63270SBarry Smith      end subroutine
4946dd63270SBarry Smith
4956dd63270SBarry Smith      subroutine F90Array4dCreateInt(array,start1,len1,                           &
4966dd63270SBarry Smith     &     start2,len2,start3,len3,start4,len4,ptr)
497*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
4986dd63270SBarry Smith      implicit none
4996dd63270SBarry Smith      PetscInt start1,len1
5006dd63270SBarry Smith      PetscInt start2,len2
5016dd63270SBarry Smith      PetscInt start3,len3
5026dd63270SBarry Smith      PetscInt start4,len4
5036dd63270SBarry Smith      PetscInt, target ::                                                         &
5046dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
5056dd63270SBarry Smith     &           start3:start3+len3-1,start4:start4+len4-1)
5066dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:,:,:)
5076dd63270SBarry Smith
5086dd63270SBarry Smith      ptr => array
5096dd63270SBarry Smith      end subroutine
5106dd63270SBarry Smith
5116dd63270SBarry Smith      subroutine F90Array4dCreateFortranAddr(array,start1,len1,                   &
5126dd63270SBarry Smith     &     start2,len2,start3,len3,start4,len4,ptr)
513*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
5146dd63270SBarry Smith      implicit none
5156dd63270SBarry Smith      PetscInt start1,len1
5166dd63270SBarry Smith      PetscInt start2,len2
5176dd63270SBarry Smith      PetscInt start3,len3
5186dd63270SBarry Smith      PetscInt start4,len4
5196dd63270SBarry Smith      PetscFortranAddr, target ::                                                 &
5206dd63270SBarry Smith     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
5216dd63270SBarry Smith     &           start3:start3+len3-1,start4:start4+len4-1)
5226dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:,:,:)
5236dd63270SBarry Smith
5246dd63270SBarry Smith      ptr => array
5256dd63270SBarry Smith      end subroutine
5266dd63270SBarry Smith
5276dd63270SBarry Smith      subroutine F90Array4dAccessScalar(ptr,address)
528*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
5296dd63270SBarry Smith      implicit none
5306dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:,:,:)
5316dd63270SBarry Smith      PetscFortranAddr address
5326dd63270SBarry Smith      PetscInt start1,start2,start3,start4
5336dd63270SBarry Smith
5346dd63270SBarry Smith      start1 = lbound(ptr,1)
5356dd63270SBarry Smith      start2 = lbound(ptr,2)
5366dd63270SBarry Smith      start3 = lbound(ptr,3)
5376dd63270SBarry Smith      start4 = lbound(ptr,4)
5386dd63270SBarry Smith      call F90Array4dGetAddrScalar(ptr(start1,start2,start3,start4),              &
5396dd63270SBarry Smith     &                             address)
5406dd63270SBarry Smith      end subroutine
5416dd63270SBarry Smith
5426dd63270SBarry Smith      subroutine F90Array4dAccessReal(ptr,address)
543*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
5446dd63270SBarry Smith      implicit none
5456dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:,:,:)
5466dd63270SBarry Smith      PetscFortranAddr address
5476dd63270SBarry Smith      PetscInt start1,start2,start3,start4
5486dd63270SBarry Smith
5496dd63270SBarry Smith      start1 = lbound(ptr,1)
5506dd63270SBarry Smith      start2 = lbound(ptr,2)
5516dd63270SBarry Smith      start3 = lbound(ptr,3)
5526dd63270SBarry Smith      start4 = lbound(ptr,4)
5536dd63270SBarry Smith      call F90Array4dGetAddrReal(ptr(start1,start2,start3,start4),                &
5546dd63270SBarry Smith     &                             address)
5556dd63270SBarry Smith      end subroutine
5566dd63270SBarry Smith
5576dd63270SBarry Smith      subroutine F90Array4dAccessInt(ptr,address)
558*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
5596dd63270SBarry Smith      implicit none
5606dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:,:,:)
5616dd63270SBarry Smith      PetscFortranAddr address
5626dd63270SBarry Smith      PetscInt start1,start2,start3,start4
5636dd63270SBarry Smith
5646dd63270SBarry Smith      start1 = lbound(ptr,1)
5656dd63270SBarry Smith      start2 = lbound(ptr,2)
5666dd63270SBarry Smith      start3 = lbound(ptr,3)
5676dd63270SBarry Smith      start4 = lbound(ptr,4)
5686dd63270SBarry Smith      call F90Array4dGetAddrInt(ptr(start1,start2,start3,start4),                 &
5696dd63270SBarry Smith     &                             address)
5706dd63270SBarry Smith      end subroutine
5716dd63270SBarry Smith
5726dd63270SBarry Smith      subroutine F90Array4dAccessFortranAddr(ptr,address)
573*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
5746dd63270SBarry Smith      implicit none
5756dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:,:,:)
5766dd63270SBarry Smith      PetscFortranAddr address
5776dd63270SBarry Smith      PetscFortranAddr start1,start2,start3,start4
5786dd63270SBarry Smith
5796dd63270SBarry Smith      start1 = lbound(ptr,1)
5806dd63270SBarry Smith      start2 = lbound(ptr,2)
5816dd63270SBarry Smith      start3 = lbound(ptr,3)
5826dd63270SBarry Smith      start4 = lbound(ptr,4)
5836dd63270SBarry Smith      call F90Array4dGetAddrFortranAddr(ptr(start1,start2,start3,                 &
5846dd63270SBarry Smith     &                                      start4),address)
5856dd63270SBarry Smith      end subroutine
5866dd63270SBarry Smith
5876dd63270SBarry Smith      subroutine F90Array4dDestroyScalar(ptr)
588*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
5896dd63270SBarry Smith      implicit none
5906dd63270SBarry Smith      PetscScalar, pointer :: ptr(:,:,:,:)
5916dd63270SBarry Smith
5926dd63270SBarry Smith      nullify(ptr)
5936dd63270SBarry Smith      end subroutine
5946dd63270SBarry Smith
5956dd63270SBarry Smith      subroutine F90Array4dDestroyReal(ptr)
596*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
5976dd63270SBarry Smith      implicit none
5986dd63270SBarry Smith      PetscReal, pointer :: ptr(:,:,:,:)
5996dd63270SBarry Smith
6006dd63270SBarry Smith      nullify(ptr)
6016dd63270SBarry Smith      end subroutine
6026dd63270SBarry Smith
6036dd63270SBarry Smith      subroutine F90Array4dDestroyInt(ptr)
604*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
6056dd63270SBarry Smith      implicit none
6066dd63270SBarry Smith      PetscInt, pointer :: ptr(:,:,:,:)
6076dd63270SBarry Smith
6086dd63270SBarry Smith      nullify(ptr)
6096dd63270SBarry Smith      end subroutine
6106dd63270SBarry Smith
6116dd63270SBarry Smith      subroutine F90Array4dDestroyFortranAddr(ptr)
612*fe66ebccSMartin Diehl      use, intrinsic :: ISO_C_binding
6136dd63270SBarry Smith      implicit none
6146dd63270SBarry Smith      PetscFortranAddr, pointer :: ptr(:,:,:,:)
6156dd63270SBarry Smith
6166dd63270SBarry Smith      nullify(ptr)
6176dd63270SBarry Smith      end subroutine
6186dd63270SBarry Smith
6196dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
6206dd63270SBarry Smith!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
6216dd63270SBarry Smith!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
622