Actual source code: f90_fwrap.F
petsc-3.10.5 2019-03-28
1: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
3: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4: #include <petsc/finclude/petscsys.h>
5: subroutine F90Array1dCreateScalar(array,start,len1,ptr)
6: implicit none
7: PetscInt start,len1
8: PetscScalar, target :: &
9: & array(start:start+len1-1)
10: PetscScalar, pointer :: ptr(:)
12: ptr => array
13: end subroutine
15: subroutine F90Array1dCreateReal(array,start,len1,ptr)
16: implicit none
17: PetscInt start,len1
18: PetscReal, target :: &
19: & array(start:start+len1-1)
20: PetscReal, pointer :: ptr(:)
22: ptr => array
23: end subroutine
25: subroutine F90Array1dCreateInt(array,start,len1,ptr)
26: implicit none
27: PetscInt start,len1
28: PetscInt, target :: &
29: & array(start:start+len1-1)
30: PetscInt, pointer :: ptr(:)
32: ptr => array
33: end subroutine
35: subroutine F90Array1dCreateFortranAddr(array,start,len1,ptr)
36: implicit none
37: PetscInt start,len1
38: PetscFortranAddr, target :: &
39: & array(start:start+len1-1)
40: PetscFortranAddr, pointer :: ptr(:)
42: ptr => array
43: end subroutine
45: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
46: subroutine F90Array1dAccessScalar(ptr,address)
47: implicit none
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: PetscReal, pointer :: ptr(:)
59: PetscFortranAddr address
60: PetscInt start
62: start = lbound(ptr,1)
63: call F90Array1dGetAddrReal(ptr(start),address)
64: end subroutine
66: subroutine F90Array1dAccessInt(ptr,address)
67: implicit none
68: PetscInt, pointer :: ptr(:)
69: PetscFortranAddr address
70: PetscInt start
72: start = lbound(ptr,1)
73: call F90Array1dGetAddrInt(ptr(start),address)
74: end subroutine
76: subroutine F90Array1dAccessFortranAddr(ptr,address)
77: implicit none
78: PetscFortranAddr, pointer :: ptr(:)
79: PetscFortranAddr address
80: PetscInt start
82: start = lbound(ptr,1)
83: call F90Array1dGetAddrFortranAddr(ptr(start),address)
84: end subroutine
86: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
87: subroutine F90Array1dDestroyScalar(ptr)
88: implicit none
89: PetscScalar, pointer :: ptr(:)
91: nullify(ptr)
92: end subroutine
94: subroutine F90Array1dDestroyReal(ptr)
95: implicit none
96: PetscReal, pointer :: ptr(:)
98: nullify(ptr)
99: end subroutine
101: subroutine F90Array1dDestroyInt(ptr)
102: implicit none
103: PetscInt, pointer :: ptr(:)
105: nullify(ptr)
106: end subroutine
108: subroutine F90Array1dDestroyFortranAddr(ptr)
109: implicit none
110: PetscFortranAddr, pointer :: ptr(:)
112: nullify(ptr)
113: end subroutine
114: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
115: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
116: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
117: subroutine F90Array2dCreateScalar(array,start1,len1, &
118: & start2,len2,ptr)
119: implicit none
120: PetscInt start1,len1
121: PetscInt start2,len2
122: PetscScalar, target :: &
123: & array(start1:start1+len1-1,start2:start2+len2-1)
124: PetscScalar, pointer :: ptr(:,:)
126: ptr => array
127: end subroutine
129: subroutine F90Array2dCreateReal(array,start1,len1, &
130: & start2,len2,ptr)
131: implicit none
132: PetscInt start1,len1
133: PetscInt start2,len2
134: PetscReal, target :: &
135: & array(start1:start1+len1-1,start2:start2+len2-1)
136: PetscReal, pointer :: ptr(:,:)
138: ptr => array
139: end subroutine
141: subroutine F90Array2dCreateInt(array,start1,len1, &
142: & start2,len2,ptr)
143: implicit none
144: PetscInt start1,len1
145: PetscInt start2,len2
146: PetscInt, target :: &
147: & array(start1:start1+len1-1,start2:start2+len2-1)
148: PetscInt, pointer :: ptr(:,:)
150: ptr => array
151: end subroutine
153: subroutine F90Array2dCreateFortranAddr(array,start1,len1, &
154: & start2,len2,ptr)
155: implicit none
156: PetscInt start1,len1
157: PetscInt start2,len2
158: PetscFortranAddr, target :: &
159: & array(start1:start1+len1-1,start2:start2+len2-1)
160: PetscFortranAddr, pointer :: ptr(:,:)
162: ptr => array
163: end subroutine
165: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
166: subroutine F90Array2dAccessScalar(ptr,address)
167: implicit none
168: PetscScalar, pointer :: ptr(:,:)
169: PetscFortranAddr address
170: PetscInt start1,start2
172: start1 = lbound(ptr,1)
173: start2 = lbound(ptr,2)
174: call F90Array2dGetAddrScalar(ptr(start1,start2),address)
175: end subroutine
177: subroutine F90Array2dAccessReal(ptr,address)
178: implicit none
179: PetscReal, pointer :: ptr(:,:)
180: PetscFortranAddr address
181: PetscInt start1,start2
183: start1 = lbound(ptr,1)
184: start2 = lbound(ptr,2)
185: call F90Array2dGetAddrReal(ptr(start1,start2),address)
186: end subroutine
188: subroutine F90Array2dAccessInt(ptr,address)
189: implicit none
190: PetscInt, pointer :: ptr(:,:)
191: PetscFortranAddr address
192: PetscInt start1,start2
194: start1 = lbound(ptr,1)
195: start2 = lbound(ptr,2)
196: call F90Array2dGetAddrInt(ptr(start1,start2),address)
197: end subroutine
199: subroutine F90Array2dAccessFortranAddr(ptr,address)
200: implicit none
201: PetscFortranAddr, pointer :: ptr(:,:)
202: PetscFortranAddr address
203: PetscInt start1,start2
205: start1 = lbound(ptr,1)
206: start2 = lbound(ptr,2)
207: call F90Array2dGetAddrFortranAddr(ptr(start1,start2),address)
208: end subroutine
210: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
211: subroutine F90Array2dDestroyScalar(ptr)
212: implicit none
213: PetscScalar, pointer :: ptr(:,:)
215: nullify(ptr)
216: end subroutine
218: subroutine F90Array2dDestroyReal(ptr)
219: implicit none
220: PetscReal, pointer :: ptr(:,:)
222: nullify(ptr)
223: end subroutine
225: subroutine F90Array2dDestroyInt(ptr)
226: implicit none
227: PetscInt, pointer :: ptr(:,:)
229: nullify(ptr)
230: end subroutine
232: subroutine F90Array2dDestroyFortranAddr(ptr)
233: implicit none
234: PetscFortranAddr, pointer :: ptr(:,:)
236: nullify(ptr)
237: end subroutine
238: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
239: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
240: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
241: subroutine F90Array3dCreateScalar(array,start1,len1, &
242: & start2,len2,start3,len3,ptr)
243: implicit none
244: PetscInt start1,len1
245: PetscInt start2,len2
246: PetscInt start3,len3
247: PetscScalar, target :: &
248: & array(start1:start1+len1-1,start2:start2+len2-1, &
249: & start3:start3+len3-1)
250: PetscScalar, pointer :: ptr(:,:,:)
252: ptr => array
253: end subroutine
255: subroutine F90Array3dCreateReal(array,start1,len1, &
256: & start2,len2,start3,len3,ptr)
257: implicit none
258: PetscInt start1,len1
259: PetscInt start2,len2
260: PetscInt start3,len3
261: PetscReal, target :: &
262: & array(start1:start1+len1-1,start2:start2+len2-1, &
263: & start3:start3+len3-1)
264: PetscReal, pointer :: ptr(:,:,:)
266: ptr => array
267: end subroutine
269: subroutine F90Array3dCreateInt(array,start1,len1, &
270: & start2,len2,start3,len3,ptr)
271: implicit none
272: PetscInt start1,len1
273: PetscInt start2,len2
274: PetscInt start3,len3
275: PetscInt, target :: &
276: & array(start1:start1+len1-1,start2:start2+len2-1, &
277: & start3:start3+len3-1)
278: PetscInt, pointer :: ptr(:,:,:)
280: ptr => array
281: end subroutine
283: subroutine F90Array3dCreateFortranAddr(array,start1,len1, &
284: & start2,len2,start3,len3,ptr)
285: implicit none
286: PetscInt start1,len1
287: PetscInt start2,len2
288: PetscInt start3,len3
289: PetscFortranAddr, target :: &
290: & array(start1:start1+len1-1,start2:start2+len2-1, &
291: & start3:start3+len3-1)
292: PetscFortranAddr, pointer :: ptr(:,:,:)
294: ptr => array
295: end subroutine
297: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
298: subroutine F90Array3dAccessScalar(ptr,address)
299: implicit none
300: PetscScalar, pointer :: ptr(:,:,:)
301: PetscFortranAddr address
302: PetscInt start1,start2,start3
304: start1 = lbound(ptr,1)
305: start2 = lbound(ptr,2)
306: start3 = lbound(ptr,3)
307: call F90Array3dGetAddrScalar(ptr(start1,start2,start3),address)
308: end subroutine
310: subroutine F90Array3dAccessReal(ptr,address)
311: implicit none
312: PetscReal, pointer :: ptr(:,:,:)
313: PetscFortranAddr address
314: PetscInt start1,start2,start3
316: start1 = lbound(ptr,1)
317: start2 = lbound(ptr,2)
318: start3 = lbound(ptr,3)
319: call F90Array3dGetAddrReal(ptr(start1,start2,start3),address)
320: end subroutine
322: subroutine F90Array3dAccessInt(ptr,address)
323: implicit none
324: PetscInt, pointer :: ptr(:,:,:)
325: PetscFortranAddr address
326: PetscInt start1,start2,start3
328: start1 = lbound(ptr,1)
329: start2 = lbound(ptr,2)
330: start3 = lbound(ptr,3)
331: call F90Array3dGetAddrInt(ptr(start1,start2,start3),address)
332: end subroutine
334: subroutine F90Array3dAccessFortranAddr(ptr,address)
335: implicit none
336: PetscFortranAddr, pointer :: ptr(:,:,:)
337: PetscFortranAddr address
338: PetscInt start1,start2,start3
340: start1 = lbound(ptr,1)
341: start2 = lbound(ptr,2)
342: start3 = lbound(ptr,3)
343: call F90Array3dGetAddrFortranAddr(ptr(start1,start2,start3), &
344: & address)
345: end subroutine
347: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
348: subroutine F90Array3dDestroyScalar(ptr)
349: implicit none
350: PetscScalar, pointer :: ptr(:,:,:)
352: nullify(ptr)
353: end subroutine
355: subroutine F90Array3dDestroyReal(ptr)
356: implicit none
357: PetscReal, pointer :: ptr(:,:,:)
359: nullify(ptr)
360: end subroutine
362: subroutine F90Array3dDestroyInt(ptr)
363: implicit none
364: PetscInt, pointer :: ptr(:,:,:)
366: nullify(ptr)
367: end subroutine
369: subroutine F90Array3dDestroyFortranAddr(ptr)
370: implicit none
371: PetscFortranAddr, pointer :: ptr(:,:,:)
373: nullify(ptr)
374: end subroutine
376: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
377: subroutine F90Array4dCreateScalar(array,start1,len1, &
378: & start2,len2,start3,len3,start4,len4,ptr)
379: implicit none
380: PetscInt start1,len1
381: PetscInt start2,len2
382: PetscInt start3,len3
383: PetscInt start4,len4
384: PetscScalar, target :: &
385: & array(start1:start1+len1-1,start2:start2+len2-1, &
386: & start3:start3+len3-1,start4:start4+len4-1)
387: PetscScalar, pointer :: ptr(:,:,:,:)
389: ptr => array
390: end subroutine
392: subroutine F90Array4dCreateReal(array,start1,len1, &
393: & start2,len2,start3,len3,start4,len4,ptr)
394: implicit none
395: PetscInt start1,len1
396: PetscInt start2,len2
397: PetscInt start3,len3
398: PetscInt start4,len4
399: PetscReal, target :: &
400: & array(start1:start1+len1-1,start2:start2+len2-1, &
401: & start3:start3+len3-1,start4:start4+len4-1)
402: PetscReal, pointer :: ptr(:,:,:,:)
404: ptr => array
405: end subroutine
407: subroutine F90Array4dCreateInt(array,start1,len1, &
408: & start2,len2,start3,len3,start4,len4,ptr)
409: implicit none
410: PetscInt start1,len1
411: PetscInt start2,len2
412: PetscInt start3,len3
413: PetscInt start4,len4
414: PetscInt, target :: &
415: & array(start1:start1+len1-1,start2:start2+len2-1, &
416: & start3:start3+len3-1,start4:start4+len4-1)
417: PetscInt, pointer :: ptr(:,:,:,:)
419: ptr => array
420: end subroutine
422: subroutine F90Array4dCreateFortranAddr(array,start1,len1, &
423: & start2,len2,start3,len3,start4,len4,ptr)
424: implicit none
425: PetscInt start1,len1
426: PetscInt start2,len2
427: PetscInt start3,len3
428: PetscInt start4,len4
429: PetscFortranAddr, target :: &
430: & array(start1:start1+len1-1,start2:start2+len2-1, &
431: & start3:start3+len3-1,start4:start4+len4-1)
432: PetscFortranAddr, pointer :: ptr(:,:,:,:)
434: ptr => array
435: end subroutine
437: subroutine F90Array4dAccessScalar(ptr,address)
438: implicit none
439: PetscScalar, pointer :: ptr(:,:,:,:)
440: PetscFortranAddr address
441: PetscInt start1,start2,start3,start4
443: start1 = lbound(ptr,1)
444: start2 = lbound(ptr,2)
445: start3 = lbound(ptr,3)
446: start4 = lbound(ptr,4)
447: call F90Array4dGetAddrScalar(ptr(start1,start2,start3,start4), &
448: & address)
449: end subroutine
451: subroutine F90Array4dAccessReal(ptr,address)
452: implicit none
453: PetscReal, pointer :: ptr(:,:,:,:)
454: PetscFortranAddr address
455: PetscInt start1,start2,start3,start4
457: start1 = lbound(ptr,1)
458: start2 = lbound(ptr,2)
459: start3 = lbound(ptr,3)
460: start4 = lbound(ptr,4)
461: call F90Array4dGetAddrReal(ptr(start1,start2,start3,start4), &
462: & address)
463: end subroutine
465: subroutine F90Array4dAccessInt(ptr,address)
466: implicit none
467: PetscInt, pointer :: ptr(:,:,:,:)
468: PetscFortranAddr address
469: PetscInt start1,start2,start3,start4
471: start1 = lbound(ptr,1)
472: start2 = lbound(ptr,2)
473: start3 = lbound(ptr,3)
474: start4 = lbound(ptr,4)
475: call F90Array4dGetAddrInt(ptr(start1,start2,start3,start4), &
476: & address)
477: end subroutine
479: subroutine F90Array4dAccessFortranAddr(ptr,address)
480: implicit none
481: PetscScalar, pointer :: ptr(:,:,:,:)
482: PetscFortranAddr address
483: PetscFortranAddr start1,start2,start3,start4
485: start1 = lbound(ptr,1)
486: start2 = lbound(ptr,2)
487: start3 = lbound(ptr,3)
488: start4 = lbound(ptr,4)
489: call F90Array4dGetAddrFortranAddr(ptr(start1,start2,start3, &
490: & start4),address)
491: end subroutine
493: subroutine F90Array4dDestroyScalar(ptr)
494: implicit none
495: PetscScalar, pointer :: ptr(:,:,:,:)
497: nullify(ptr)
498: end subroutine
500: subroutine F90Array4dDestroyReal(ptr)
501: implicit none
502: PetscReal, pointer :: ptr(:,:,:,:)
504: nullify(ptr)
505: end subroutine
507: subroutine F90Array4dDestroyInt(ptr)
508: implicit none
509: PetscInt, pointer :: ptr(:,:,:,:)
511: nullify(ptr)
512: end subroutine
514: subroutine F90Array4dDestroyFortranAddr(ptr)
515: implicit none
516: PetscFortranAddr, pointer :: ptr(:,:,:,:)
518: nullify(ptr)
519: end subroutine
521: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
522: !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
523: !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!