Actual source code: f90_fwrap.F
1: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
3: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4: subroutine F90Array1dCreateScalar(array,start,len,ptr)
5: implicit none
6: #include "finclude/petscsys.h"
7: PetscInt start,len
8: PetscScalar, target :: array(start+len-1)
9: PetscScalar, pointer :: ptr(:)
11: ptr => array(start:start+len-1)
12: end subroutine
14: subroutine F90Array1dCreateReal(array,start,len,ptr)
15: implicit none
16: #include "finclude/petscsys.h"
17: PetscInt start,len
18: PetscReal, target :: array(start+len-1)
19: PetscReal, pointer :: ptr(:)
21: ptr => array(start:start+len-1)
22: end subroutine
24: subroutine F90Array1dCreateInt(array,start,len,ptr)
25: implicit none
26: #include "finclude/petscsys.h"
27: PetscInt start,len
28: PetscInt, target :: array(start+len-1)
29: PetscInt, pointer :: ptr(:)
31: ptr => array(start:start+len-1)
32: end subroutine
34: subroutine F90Array1dCreateFortranAddr(array,start,len,ptr)
35: implicit none
36: #include "finclude/petscsys.h"
37: PetscInt start,len
38: PetscFortranAddr, target :: array(start+len-1)
39: PetscFortranAddr, pointer :: ptr(:)
41: ptr => array(start:start+len-1)
42: end subroutine
44: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
45: subroutine F90Array1dAccessScalar(ptr,address)
46: implicit none
47: #include "finclude/petscsys.h"
48: PetscScalar, pointer :: ptr(:)
49: PetscFortranAddr address
50: PetscInt start
52: start = lbound(ptr,1)
53: call F90Array1dGetAddrScalar(ptr(start),address)
54: end subroutine
56: subroutine F90Array1dAccessReal(ptr,address)
57: implicit none
58: #include "finclude/petscsys.h"
59: PetscReal, pointer :: ptr(:)
60: PetscFortranAddr address
61: PetscInt start
63: start = lbound(ptr,1)
64: call F90Array1dGetAddrReal(ptr(start),address)
65: end subroutine
67: subroutine F90Array1dAccessInt(ptr,address)
68: implicit none
69: #include "finclude/petscsys.h"
70: PetscInt, pointer :: ptr(:)
71: PetscFortranAddr address
72: PetscInt start
74: start = lbound(ptr,1)
75: call F90Array1dGetAddrInt(ptr(start),address)
76: end subroutine
78: subroutine F90Array1dAccessFortranAddr(ptr,address)
79: implicit none
80: #include "finclude/petscsys.h"
81: PetscFortranAddr, pointer :: ptr(:)
82: PetscFortranAddr address
83: PetscInt start
85: start = lbound(ptr,1)
86: call F90Array1dGetAddrFortranAddr(ptr(start),address)
87: end subroutine
88:
89: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
90: subroutine F90Array1dDestroyScalar(ptr)
91: implicit none
92: #include "finclude/petscsys.h"
93: PetscScalar, pointer :: ptr(:)
95: nullify(ptr)
96: end subroutine
98: subroutine F90Array1dDestroyReal(ptr)
99: implicit none
100: #include "finclude/petscsys.h"
101: PetscReal, pointer :: ptr(:)
103: nullify(ptr)
104: end subroutine
106: subroutine F90Array1dDestroyInt(ptr)
107: implicit none
108: #include "finclude/petscsys.h"
109: PetscInt, pointer :: ptr(:)
111: nullify(ptr)
112: end subroutine
114: subroutine F90Array1dDestroyFortranAddr(ptr)
115: implicit none
116: #include "finclude/petscsys.h"
117: PetscFortranAddr, pointer :: ptr(:)
119: nullify(ptr)
120: end subroutine
121: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
122: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
123: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
124: subroutine F90Array2dCreateScalar(array,start1,len1, &
125: & start2,len2,ptr)
126: implicit none
127: #include "finclude/petscsys.h"
128: PetscInt start1,len1
129: PetscInt start2,len2
130: PetscScalar, target :: &
131: & array(start1+len1-1,start2+len2-1)
132: PetscScalar, pointer :: ptr(:,:)
134: ptr => array(start1:len1,start2:len2)
135: end subroutine
137: subroutine F90Array2dCreateReal(array,start1,len1, &
138: & start2,len2,ptr)
139: implicit none
140: #include "finclude/petscsys.h"
141: PetscInt start1,len1
142: PetscInt start2,len2
143: PetscReal, target :: &
144: & array(start1+len1-1,start2+len2-1)
145: PetscReal, pointer :: ptr(:,:)
147: ptr => array(start1:len1,start2:len2)
148: end subroutine
150: subroutine F90Array2dCreateInt(array,start1,len1, &
151: & start2,len2,ptr)
152: implicit none
153: #include "finclude/petscsys.h"
154: PetscInt start1,len1
155: PetscInt start2,len2
156: PetscInt, target :: &
157: & array(start1+len1-1,start2+len2-1)
158: PetscInt, pointer :: ptr(:,:)
160: ptr => array(start1:len1,start2:len2)
161: end subroutine
163: subroutine F90Array2dCreateFortranAddr(array,start1,len1, &
164: & start2,len2,ptr)
165: implicit none
166: #include "finclude/petscsys.h"
167: PetscInt start1,len1
168: PetscInt start2,len2
169: PetscFortranAddr, target :: &
170: & array(start1+len1-1,start2+len2-1)
171: PetscFortranAddr, pointer :: ptr(:,:)
173: ptr => array(start1:len1,start2:len2)
174: end subroutine
176: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
177: subroutine F90Array2dAccessScalar(ptr,address)
178: implicit none
179: #include "finclude/petscsys.h"
180: PetscScalar, pointer :: ptr(:,:)
181: PetscFortranAddr address
182: PetscInt start1,start2
184: start1 = lbound(ptr,1)
185: start2 = lbound(ptr,2)
186: call F90Array2dGetAddrScalar(ptr(start1,start2),address)
187: end subroutine
189: subroutine F90Array2dAccessReal(ptr,address)
190: implicit none
191: #include "finclude/petscsys.h"
192: PetscReal, pointer :: ptr(:,:)
193: PetscFortranAddr address
194: PetscInt start1,start2
196: start1 = lbound(ptr,1)
197: start2 = lbound(ptr,2)
198: call F90Array2dGetAddrReal(ptr(start1,start2),address)
199: end subroutine
201: subroutine F90Array2dAccessInt(ptr,address)
202: implicit none
203: #include "finclude/petscsys.h"
204: PetscInt, pointer :: ptr(:,:)
205: PetscFortranAddr address
206: PetscInt start1,start2
208: start1 = lbound(ptr,1)
209: start2 = lbound(ptr,2)
210: call F90Array2dGetAddrInt(ptr(start1,start2),address)
211: end subroutine
213: subroutine F90Array2dAccessFortranAddr(ptr,address)
214: implicit none
215: #include "finclude/petscsys.h"
216: PetscFortranAddr, pointer :: ptr(:,:)
217: PetscFortranAddr address
218: PetscInt start1,start2
220: start1 = lbound(ptr,1)
221: start2 = lbound(ptr,2)
222: call F90Array2dGetAddrFortranAddr(ptr(start1,start2),address)
223: end subroutine
224:
225: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
226: subroutine F90Array2dDestroyScalar(ptr)
227: implicit none
228: #include "finclude/petscsys.h"
229: PetscScalar, pointer :: ptr(:,:)
231: nullify(ptr)
232: end subroutine
234: subroutine F90Array2dDestroyReal(ptr)
235: implicit none
236: #include "finclude/petscsys.h"
237: PetscReal, pointer :: ptr(:,:)
239: nullify(ptr)
240: end subroutine
242: subroutine F90Array2dDestroyInt(ptr)
243: implicit none
244: #include "finclude/petscsys.h"
245: PetscInt, pointer :: ptr(:,:)
247: nullify(ptr)
248: end subroutine
250: subroutine F90Array2dDestroyFortranAddr(ptr)
251: implicit none
252: #include "finclude/petscsys.h"
253: PetscFortranAddr, pointer :: ptr(:,:)
255: nullify(ptr)
256: end subroutine
257: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
258: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
259: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!