Actual source code: ex21f90.F
petsc-3.7.3 2016-08-01
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: ! -----------------------------------------------------------------------
7: #include <petsc/finclude/petscsysdef.h>
8: #include <petsc/finclude/petscvecdef.h>
10: module mymodule
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 mymodule
24: implicit none
25: #include <petsc/finclude/petscsys.h>
26: PetscInt start,len
27: type(MyStruct), target :: &
28: & array(start:start+len-1)
29: type(MyStruct), pointer :: ptr(:)
31: ptr => array
32: end subroutine
34: subroutine F90Array1dAccessMyStruct(ptr,address)
35: use mymodule
36: implicit none
37: #include <petsc/finclude/petscsys.h>
38: type(MyStruct), pointer :: ptr(:)
39: PetscFortranAddr address
40: PetscInt start
42: start = lbound(ptr,1)
43: call F90Array1dGetAddrMyStruct(ptr(start),address)
44: end subroutine
46: subroutine F90Array1dDestroyMyStruct(ptr)
47: use mymodule
48: implicit none
49: #include <petsc/finclude/petscsys.h>
50: type(MyStruct), pointer :: ptr(:)
52: nullify(ptr)
53: end subroutine
56: program main
57: use mymodule
58: implicit none
60: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: ! Include files
62: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: !
64: ! The following include statements are required for Fortran programs
65: ! that use PETSc vectors:
66: ! petscsys.h - base PETSc routines
67: ! petscvec.h - vectors
68: ! petscvec.h90 - to allow access to Fortran90 features of vectors
69: !
70: ! Additional include statements may be needed if using additional
71: ! PETSc routines in a Fortran program, e.g.,
72: ! petscviewer.h - viewers
73: ! petscis.h - index sets
74: !
75: #include <petsc/finclude/petscsys.h>
76: #include <petsc/finclude/petscviewer.h>
77: #include <petsc/finclude/petscvec.h>
78: #include <petsc/finclude/petscvec.h90>
80: !
81: ! These two routines are defined in ex21.c they create the Fortran pointer to the derived type
82: !
83: Interface
84: Subroutine VecGetArrayMyStruct(v,array,ierr)
85: use mymodule
86: type(MyStruct), pointer :: array(:)
87: PetscErrorCode ierr
88: Vec v
89: End Subroutine
90: End Interface
92: Interface
93: Subroutine VecRestoreArrayMyStruct(v,array,ierr)
94: use mymodule
95: type(MyStruct), pointer :: array(:)
96: PetscErrorCode ierr
97: Vec v
98: End Subroutine
99: End Interface
101: !
102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: ! Variable declarations
104: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: !
106: ! Variables:
107: ! x, y, w - vectors
108: ! z - array of vectors
109: !
110: Vec x,y
111: type(MyStruct), pointer :: xarray(:)
112: PetscInt n
113: PetscErrorCode ierr
114: PetscBool flg
115: integer i
117: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: ! Beginning of program
119: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
122: n = 30
124: call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER, &
125: & '-n',n,flg,ierr)
126: call VecCreate(PETSC_COMM_WORLD,x,ierr)
127: call VecSetSizes(x,PETSC_DECIDE,n,ierr)
128: call VecSetFromOptions(x,ierr)
129: call VecDuplicate(x,y,ierr)
131: call VecGetArrayMyStruct(x,xarray,ierr)
132: do i=1,10
133: xarray(i)%a = i
134: xarray(i)%b = 100*i
135: xarray(i)%c = 10000*i
136: enddo
138: call VecRestoreArrayMyStruct(x,xarray,ierr)
139: call VecView(x,PETSC_VIEWER_STDOUT_SELF,ierr)
140: call VecGetArrayMyStruct(x,xarray,ierr)
141: do i = 1 , 10
142: write(*,*) abs(xarray(i)%a),abs(xarray(i)%b),abs(xarray(i)%c)
143: end do
144: call VecRestoreArrayMyStruct(x,xarray,ierr)
147: call VecDestroy(x,ierr)
148: call VecDestroy(y,ierr)
149: call PetscFinalize(ierr)
151: end