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 :: &
9: & array(start:start+len-1)
10: PetscScalar, pointer :: ptr(:)
12: ptr => array
13: end subroutine
15: subroutine F90Array1dCreateReal(array,start,len,ptr)
16: implicit none
17: #include <finclude/petscsys.h>
18: PetscInt start,len
19: PetscReal, target :: &
20: & array(start:start+len-1)
21: PetscReal, pointer :: ptr(:)
23: ptr => array
24: end subroutine
26: subroutine F90Array1dCreateInt(array,start,len,ptr)
27: implicit none
28: #include <finclude/petscsys.h>
29: PetscInt start,len
30: PetscInt, target :: &
31: & array(start:start+len-1)
32: PetscInt, pointer :: ptr(:)
34: ptr => array
35: end subroutine
37: subroutine F90Array1dCreateFortranAddr(array,start,len,ptr)
38: implicit none
39: #include <finclude/petscsys.h>
40: PetscInt start,len
41: PetscFortranAddr, target :: &
42: & array(start:start+len-1)
43: PetscFortranAddr, pointer :: ptr(:)
45: ptr => array
46: end subroutine
48: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
49: subroutine F90Array1dAccessScalar(ptr,address)
50: implicit none
51: #include <finclude/petscsys.h>
52: PetscScalar, pointer :: ptr(:)
53: PetscFortranAddr address
54: PetscInt start
56: start = lbound(ptr,1)
57: call F90Array1dGetAddrScalar(ptr(start),address)
58: end subroutine
60: subroutine F90Array1dAccessReal(ptr,address)
61: implicit none
62: #include <finclude/petscsys.h>
63: PetscReal, pointer :: ptr(:)
64: PetscFortranAddr address
65: PetscInt start
67: start = lbound(ptr,1)
68: call F90Array1dGetAddrReal(ptr(start),address)
69: end subroutine
71: subroutine F90Array1dAccessInt(ptr,address)
72: implicit none
73: #include <finclude/petscsys.h>
74: PetscInt, pointer :: ptr(:)
75: PetscFortranAddr address
76: PetscInt start
78: start = lbound(ptr,1)
79: call F90Array1dGetAddrInt(ptr(start),address)
80: end subroutine
82: subroutine F90Array1dAccessFortranAddr(ptr,address)
83: implicit none
84: #include <finclude/petscsys.h>
85: PetscFortranAddr, pointer :: ptr(:)
86: PetscFortranAddr address
87: PetscInt start
89: start = lbound(ptr,1)
90: call F90Array1dGetAddrFortranAddr(ptr(start),address)
91: end subroutine
92:
93: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
94: subroutine F90Array1dDestroyScalar(ptr)
95: implicit none
96: #include <finclude/petscsys.h>
97: PetscScalar, pointer :: ptr(:)
99: nullify(ptr)
100: end subroutine
102: subroutine F90Array1dDestroyReal(ptr)
103: implicit none
104: #include <finclude/petscsys.h>
105: PetscReal, pointer :: ptr(:)
107: nullify(ptr)
108: end subroutine
110: subroutine F90Array1dDestroyInt(ptr)
111: implicit none
112: #include <finclude/petscsys.h>
113: PetscInt, pointer :: ptr(:)
115: nullify(ptr)
116: end subroutine
118: subroutine F90Array1dDestroyFortranAddr(ptr)
119: implicit none
120: #include <finclude/petscsys.h>
121: PetscFortranAddr, pointer :: ptr(:)
123: nullify(ptr)
124: end subroutine
125: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
126: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
127: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
128: subroutine F90Array2dCreateScalar(array,start1,len1, &
129: & start2,len2,ptr)
130: implicit none
131: #include <finclude/petscsys.h>
132: PetscInt start1,len1
133: PetscInt start2,len2
134: PetscScalar, target :: &
135: & array(start1:start1+len1-1,start2:start2+len2-1)
136: PetscScalar, pointer :: ptr(:,:)
138: ptr => array
139: end subroutine
141: subroutine F90Array2dCreateReal(array,start1,len1, &
142: & start2,len2,ptr)
143: implicit none
144: #include <finclude/petscsys.h>
145: PetscInt start1,len1
146: PetscInt start2,len2
147: PetscReal, target :: &
148: & array(start1:start1+len1-1,start2:start2+len2-1)
149: PetscReal, pointer :: ptr(:,:)
151: ptr => array
152: end subroutine
154: subroutine F90Array2dCreateInt(array,start1,len1, &
155: & start2,len2,ptr)
156: implicit none
157: #include <finclude/petscsys.h>
158: PetscInt start1,len1
159: PetscInt start2,len2
160: PetscInt, target :: &
161: & array(start1:start1+len1-1,start2:start2+len2-1)
162: PetscInt, pointer :: ptr(:,:)
164: ptr => array
165: end subroutine
167: subroutine F90Array2dCreateFortranAddr(array,start1,len1, &
168: & start2,len2,ptr)
169: implicit none
170: #include <finclude/petscsys.h>
171: PetscInt start1,len1
172: PetscInt start2,len2
173: PetscFortranAddr, target :: &
174: & array(start1:start1+len1-1,start2:start2+len2-1)
175: PetscFortranAddr, pointer :: ptr(:,:)
177: ptr => array
178: end subroutine
180: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
181: subroutine F90Array2dAccessScalar(ptr,address)
182: implicit none
183: #include <finclude/petscsys.h>
184: PetscScalar, pointer :: ptr(:,:)
185: PetscFortranAddr address
186: PetscInt start1,start2
188: start1 = lbound(ptr,1)
189: start2 = lbound(ptr,2)
190: call F90Array2dGetAddrScalar(ptr(start1,start2),address)
191: end subroutine
193: subroutine F90Array2dAccessReal(ptr,address)
194: implicit none
195: #include <finclude/petscsys.h>
196: PetscReal, pointer :: ptr(:,:)
197: PetscFortranAddr address
198: PetscInt start1,start2
200: start1 = lbound(ptr,1)
201: start2 = lbound(ptr,2)
202: call F90Array2dGetAddrReal(ptr(start1,start2),address)
203: end subroutine
205: subroutine F90Array2dAccessInt(ptr,address)
206: implicit none
207: #include <finclude/petscsys.h>
208: PetscInt, pointer :: ptr(:,:)
209: PetscFortranAddr address
210: PetscInt start1,start2
212: start1 = lbound(ptr,1)
213: start2 = lbound(ptr,2)
214: call F90Array2dGetAddrInt(ptr(start1,start2),address)
215: end subroutine
217: subroutine F90Array2dAccessFortranAddr(ptr,address)
218: implicit none
219: #include <finclude/petscsys.h>
220: PetscFortranAddr, pointer :: ptr(:,:)
221: PetscFortranAddr address
222: PetscInt start1,start2
224: start1 = lbound(ptr,1)
225: start2 = lbound(ptr,2)
226: call F90Array2dGetAddrFortranAddr(ptr(start1,start2),address)
227: end subroutine
228:
229: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
230: subroutine F90Array2dDestroyScalar(ptr)
231: implicit none
232: #include <finclude/petscsys.h>
233: PetscScalar, pointer :: ptr(:,:)
235: nullify(ptr)
236: end subroutine
238: subroutine F90Array2dDestroyReal(ptr)
239: implicit none
240: #include <finclude/petscsys.h>
241: PetscReal, pointer :: ptr(:,:)
243: nullify(ptr)
244: end subroutine
246: subroutine F90Array2dDestroyInt(ptr)
247: implicit none
248: #include <finclude/petscsys.h>
249: PetscInt, pointer :: ptr(:,:)
251: nullify(ptr)
252: end subroutine
254: subroutine F90Array2dDestroyFortranAddr(ptr)
255: implicit none
256: #include <finclude/petscsys.h>
257: PetscFortranAddr, pointer :: ptr(:,:)
259: nullify(ptr)
260: end subroutine
261: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
262: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
263: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
264: subroutine F90Array3dCreateScalar(array,start1,len1, &
265: & start2,len2,start3,len3,ptr)
266: implicit none
267: #include <finclude/petscsys.h>
268: PetscInt start1,len1
269: PetscInt start2,len2
270: PetscInt start3,len3
271: PetscScalar, target :: &
272: & array(start1:start1+len1-1,start2:start2+len2-1, &
273: & start3:start3+len3-1)
274: PetscScalar, pointer :: ptr(:,:,:)
276: ptr => array
277: end subroutine
279: subroutine F90Array3dCreateReal(array,start1,len1, &
280: & start2,len2,start3,len3,ptr)
281: implicit none
282: #include <finclude/petscsys.h>
283: PetscInt start1,len1
284: PetscInt start2,len2
285: PetscInt start3,len3
286: PetscReal, target :: &
287: & array(start1:start1+len1-1,start2:start2+len2-1, &
288: & start3:start3+len3-1)
289: PetscReal, pointer :: ptr(:,:,:)
291: ptr => array
292: end subroutine
294: subroutine F90Array3dCreateInt(array,start1,len1, &
295: & start2,len2,start3,len3,ptr)
296: implicit none
297: #include <finclude/petscsys.h>
298: PetscInt start1,len1
299: PetscInt start2,len2
300: PetscInt start3,len3
301: PetscInt, target :: &
302: & array(start1:start1+len1-1,start2:start2+len2-1, &
303: & start3:start3+len3-1)
304: PetscInt, pointer :: ptr(:,:,:)
306: ptr => array
307: end subroutine
309: subroutine F90Array3dCreateFortranAddr(array,start1,len1, &
310: & start2,len2,start3,len3,ptr)
311: implicit none
312: #include <finclude/petscsys.h>
313: PetscInt start1,len1
314: PetscInt start2,len2
315: PetscInt start3,len3
316: PetscFortranAddr, target :: &
317: & array(start1:start1+len1-1,start2:start2+len2-1, &
318: & start3:start3+len3-1)
319: PetscFortranAddr, pointer :: ptr(:,:,:)
321: ptr => array
322: end subroutine
324: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
325: subroutine F90Array4dCreateScalar(array,start1,len1, &
326: & start2,len2,start3,len3,start4,len4,ptr)
327: implicit none
328: #include <finclude/petscsys.h>
329: PetscInt start1,len1
330: PetscInt start2,len2
331: PetscInt start3,len3
332: PetscInt start4,len4
333: PetscScalar, target :: &
334: & array(start1:start1+len1-1,start2:start2+len2-1, &
335: & start3:start3+len3-1,start4:start4+len4-1)
336: PetscScalar, pointer :: ptr(:,:,:,:)
338: ptr => array
339: end subroutine
341: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
342: subroutine F90Array3dAccessScalar(ptr,address)
343: implicit none
344: #include <finclude/petscsys.h>
345: PetscScalar, pointer :: ptr(:,:,:)
346: PetscFortranAddr address
347: PetscInt start1,start2,start3
349: start1 = lbound(ptr,1)
350: start2 = lbound(ptr,2)
351: start3 = lbound(ptr,3)
352: call F90Array3dGetAddrScalar(ptr(start1,start2,start3),address)
353: end subroutine
355: subroutine F90Array3dAccessReal(ptr,address)
356: implicit none
357: #include <finclude/petscsys.h>
358: PetscReal, pointer :: ptr(:,:,:)
359: PetscFortranAddr address
360: PetscInt start1,start2,start3
362: start1 = lbound(ptr,1)
363: start2 = lbound(ptr,2)
364: start3 = lbound(ptr,3)
365: call F90Array3dGetAddrReal(ptr(start1,start2,start3),address)
366: end subroutine
368: subroutine F90Array3dAccessInt(ptr,address)
369: implicit none
370: #include <finclude/petscsys.h>
371: PetscInt, pointer :: ptr(:,:,:)
372: PetscFortranAddr address
373: PetscInt start1,start2,start3
375: start1 = lbound(ptr,1)
376: start2 = lbound(ptr,2)
377: start3 = lbound(ptr,3)
378: call F90Array3dGetAddrInt(ptr(start1,start2,start3),address)
379: end subroutine
381: subroutine F90Array3dAccessFortranAddr(ptr,address)
382: implicit none
383: #include <finclude/petscsys.h>
384: PetscFortranAddr, pointer :: ptr(:,:,:)
385: PetscFortranAddr address
386: PetscInt start1,start2,start3
388: start1 = lbound(ptr,1)
389: start2 = lbound(ptr,2)
390: start3 = lbound(ptr,3)
391: call F90Array3dGetAddrFortranAddr(ptr(start1,start2,start3), &
392: & address)
393: end subroutine
394:
395: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
396: subroutine F90Array3dDestroyScalar(ptr)
397: implicit none
398: #include <finclude/petscsys.h>
399: PetscScalar, pointer :: ptr(:,:,:)
401: nullify(ptr)
402: end subroutine
404: subroutine F90Array3dDestroyReal(ptr)
405: implicit none
406: #include <finclude/petscsys.h>
407: PetscReal, pointer :: ptr(:,:,:)
409: nullify(ptr)
410: end subroutine
412: subroutine F90Array3dDestroyInt(ptr)
413: implicit none
414: #include <finclude/petscsys.h>
415: PetscInt, pointer :: ptr(:,:,:)
417: nullify(ptr)
418: end subroutine
420: subroutine F90Array3dDestroyFortranAddr(ptr)
421: implicit none
422: #include <finclude/petscsys.h>
423: PetscFortranAddr, pointer :: ptr(:,:,:)
425: nullify(ptr)
426: end subroutine
428: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
429: subroutine F90Array4dDestroyScalar(ptr)
430: implicit none
431: #include <finclude/petscsys.h>
432: PetscScalar, pointer :: ptr(:,:,:)
434: nullify(ptr)
435: end subroutine
437: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
438: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
439: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!