Actual source code: ex21f90.F90

  1: !
  2: !
  3: !    Demonstrates how one may access entries of a PETSc Vec as if it was an array of Fortran derived types
  4: !
  5: !
  6: ! -----------------------------------------------------------------------

  8:       module ex21f90module
  9: #include <petsc/finclude/petscsys.h>
 10:       use petscsys
 11:       type MyStruct
 12:         sequence
 13:         PetscScalar :: a,b,c
 14:       end type MyStruct
 15:       end module

 17: !
 18: !  These routines are used internally by the C functions VecGetArrayMyStruct() and VecRestoreArrayMyStruct()
 19: !  Because Fortran requires "knowing" exactly what derived types the pointers to point too, these have to be
 20: !  customized for exactly the derived type in question
 21: !
 22:       subroutine F90Array1dCreateMyStruct(array,start,len,ptr)
 23:       use ex21f90module
 24:       implicit none
 25:       PetscInt start,len
 26:       type(MyStruct), target :: array(start:start+len-1)
 27:       type(MyStruct), pointer :: ptr(:)

 29:       ptr => array
 30:       end subroutine

 32:       subroutine F90Array1dAccessMyStruct(ptr,address)
 33:       use ex21f90module
 34:       implicit none
 35:       type(MyStruct), pointer :: ptr(:)
 36:       PetscFortranAddr address
 37:       PetscInt start

 39:       start = lbound(ptr,1)
 40:       call F90Array1dGetAddrMyStruct(ptr(start),address)
 41:       end subroutine

 43:       subroutine F90Array1dDestroyMyStruct(ptr)
 44:       use ex21f90module
 45:       implicit none
 46:       type(MyStruct), pointer :: ptr(:)

 48:       nullify(ptr)
 49:       end subroutine

 51:       program main
 52: #include <petsc/finclude/petscvec.h>
 53:           use petscvec
 54:       use ex21f90module
 55:       implicit none

 57: !
 58: !
 59: !   These two routines are defined in ex21.c they create the Fortran pointer to the derived type
 60: !
 61:       Interface
 62:         Subroutine VecGetArrayMyStruct(v,array,ierr)
 63: #include <petsc/finclude/petscvec.h>
 64:           use petscvec
 65:           use ex21f90module
 66:           type(MyStruct), pointer :: array(:)
 67:           PetscErrorCode ierr
 68:           Vec     v
 69:         End Subroutine
 70:       End Interface

 72:       Interface
 73:         Subroutine VecRestoreArrayMyStruct(v,array,ierr)
 74: #include <petsc/finclude/petscvec.h>
 75:           use petscvec
 76:           use ex21f90module
 77:           type(MyStruct), pointer :: array(:)
 78:           PetscErrorCode ierr
 79:           Vec     v
 80:         End Subroutine
 81:       End Interface

 83: !
 84: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85: !                   Variable declarations
 86: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87: !
 88: !  Variables:
 89: !     x, y, w - vectors
 90: !     z       - array of vectors
 91: !
 92:       Vec              x,y
 93:       type(MyStruct),  pointer :: xarray(:)
 94:       PetscInt         n
 95:       PetscErrorCode   ierr
 96:       PetscBool        flg
 97:       integer          i

 99: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: !                 Beginning of program
101: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

103:       PetscCallA(PetscInitialize(ierr))
104:       n     = 30

106:       PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))
107:       PetscCallA(VecCreate(PETSC_COMM_WORLD,x,ierr))
108:       PetscCallA(VecSetSizes(x,PETSC_DECIDE,n,ierr))
109:       PetscCallA(VecSetFromOptions(x,ierr))
110:       PetscCallA(VecDuplicate(x,y,ierr))

112:       PetscCallA(VecGetArrayMyStruct(x,xarray,ierr))
113:       do i=1,10
114:       xarray(i)%a = i
115:       xarray(i)%b = 100*i
116:       xarray(i)%c = 10000*i
117:       enddo

119:       PetscCallA(VecRestoreArrayMyStruct(x,xarray,ierr))
120:       PetscCallA(VecView(x,PETSC_VIEWER_STDOUT_SELF,ierr))
121:       PetscCallA(VecGetArrayMyStruct(x,xarray,ierr))
122:       do i = 1 , 10
123:         write(*,*) abs(xarray(i)%a),abs(xarray(i)%b),abs(xarray(i)%c)
124:       end do
125:       PetscCallA(VecRestoreArrayMyStruct(x,xarray,ierr))

127:       PetscCallA(VecDestroy(x,ierr))
128:       PetscCallA(VecDestroy(y,ierr))
129:       PetscCallA(PetscFinalize(ierr))

131:       end

133: !/*TEST
134: !   build:
135: !     depends: ex21.c
136: !
137: !   test:
138: !
139: !TEST*/