xref: /petsc/src/sys/ftn-src/f90_fwrap.F90 (revision 61c8d4ed94f37c13af459aa1b8bd670fce1ab5da)
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
3!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4#include <petsc/finclude/petscsys.h>
5      subroutine F90Array1dCreateScalar(array,start,len1,ptr)
6      use, intrinsic :: ISO_C_binding
7      implicit none
8      PetscInt start,len1
9      PetscScalar, target ::                                                      &
10     &     array(start:start+len1-1)
11      PetscScalar, pointer :: ptr(:)
12
13      ptr => array
14      end subroutine
15
16      subroutine F90Array1dCreateReal(array,start,len1,ptr)
17      use, intrinsic :: ISO_C_binding
18      implicit none
19      PetscInt start,len1
20      PetscReal, target ::                                                        &
21     &     array(start:start+len1-1)
22      PetscReal, pointer :: ptr(:)
23
24      ptr => array
25      end subroutine
26
27      subroutine F90Array1dCreateInt(array,start,len1,ptr)
28      use, intrinsic :: ISO_C_binding
29      implicit none
30      PetscInt start,len1
31      PetscInt, target ::                                                         &
32     &     array(start:start+len1-1)
33      PetscInt, pointer :: ptr(:)
34
35      ptr => array
36      end subroutine
37
38      subroutine F90Array1dCreateMPIInt(array,start,len1,ptr)
39      use, intrinsic :: ISO_C_binding
40      implicit none
41      PetscInt start,len1
42      PetscMPIInt, target ::                                                      &
43      &     array(start:start+len1-1)
44      PetscMPIInt, pointer :: ptr(:)
45
46      ptr => array
47      end subroutine
48
49      subroutine F90Array1dCreateFortranAddr(array,start,len1,ptr)
50      use, intrinsic :: ISO_C_binding
51      implicit none
52      PetscInt start,len1
53      PetscFortranAddr, target ::                                                 &
54     &     array(start:start+len1-1)
55      PetscFortranAddr, pointer :: ptr(:)
56
57      ptr => array
58      end subroutine
59
60!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
61      subroutine F90Array1dAccessScalar(ptr,address)
62      use, intrinsic :: ISO_C_binding
63      implicit none
64      PetscScalar, pointer :: ptr(:)
65      PetscFortranAddr address
66      PetscInt start
67
68      if (associated(ptr) .eqv. .false.) then
69        address = 0
70      else
71        start = lbound(ptr,1)
72        call F90Array1dGetAddrScalar(ptr(start),address)
73      endif
74      end subroutine
75
76      subroutine F90Array1dAccessReal(ptr,address)
77      use, intrinsic :: ISO_C_binding
78      implicit none
79      PetscReal, pointer :: ptr(:)
80      PetscFortranAddr address
81      PetscInt start
82
83      if (associated(ptr) .eqv. .false.) then
84        address = 0
85      else
86        start = lbound(ptr,1)
87        call F90Array1dGetAddrReal(ptr(start),address)
88      endif
89      end subroutine
90
91      subroutine F90Array1dAccessInt(ptr,address)
92      use, intrinsic :: ISO_C_binding
93      implicit none
94      PetscInt, pointer :: ptr(:)
95      PetscFortranAddr address
96      PetscInt start
97
98      if (associated(ptr) .eqv. .false.) then
99        address = 0
100      else
101        start = lbound(ptr,1)
102        call F90Array1dGetAddrInt(ptr(start),address)
103      endif
104      end subroutine
105
106      subroutine F90Array1dAccessMPIInt(ptr,address)
107      use, intrinsic :: ISO_C_binding
108        implicit none
109        PetscMPIInt, pointer :: ptr(:)
110        PetscFortranAddr address
111        PetscInt start
112
113        if (associated(ptr) .eqv. .false.) then
114          address = 0
115        else
116          start = lbound(ptr,1)
117          call F90Array1dGetAddrMPIInt(ptr(start),address)
118        endif
119        end subroutine
120
121      subroutine F90Array1dAccessFortranAddr(ptr,address)
122      use, intrinsic :: ISO_C_binding
123      implicit none
124      PetscFortranAddr, pointer :: ptr(:)
125      PetscFortranAddr address
126      PetscInt start
127
128      if (associated(ptr) .eqv. .false.) then
129        address = 0
130      else
131        start = lbound(ptr,1)
132        call F90Array1dGetAddrFortranAddr(ptr(start),address)
133      endif
134      end subroutine
135
136!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
137      subroutine F90Array1dDestroyScalar(ptr)
138      use, intrinsic :: ISO_C_binding
139      implicit none
140      PetscScalar, pointer :: ptr(:)
141
142      nullify(ptr)
143      end subroutine
144
145      subroutine F90Array1dDestroyReal(ptr)
146      use, intrinsic :: ISO_C_binding
147      implicit none
148      PetscReal, pointer :: ptr(:)
149
150      nullify(ptr)
151      end subroutine
152
153      subroutine F90Array1dDestroyInt(ptr)
154      use, intrinsic :: ISO_C_binding
155      implicit none
156      PetscInt, pointer :: ptr(:)
157
158      nullify(ptr)
159      end subroutine
160
161      subroutine F90Array1dDestroyMPIInt(ptr)
162      use, intrinsic :: ISO_C_binding
163      implicit none
164      PetscMPIInt, pointer :: ptr(:)
165
166      nullify(ptr)
167      end subroutine
168
169      subroutine F90Array1dDestroyFortranAddr(ptr)
170      use, intrinsic :: ISO_C_binding
171      implicit none
172      PetscFortranAddr, pointer :: ptr(:)
173
174      nullify(ptr)
175      end subroutine
176!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
177!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
178!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
179      subroutine F90Array2dCreateScalar(array,start1,len1,                        &
180     &     start2,len2,ptr)
181      use, intrinsic :: ISO_C_binding
182      implicit none
183      PetscInt start1,len1
184      PetscInt start2,len2
185      PetscScalar, target ::                                                      &
186     &     array(start1:start1+len1-1,start2:start2+len2-1)
187      PetscScalar, pointer :: ptr(:,:)
188
189      ptr => array
190      end subroutine
191
192      subroutine F90Array2dCreateReal(array,start1,len1,                          &
193     &     start2,len2,ptr)
194      use, intrinsic :: ISO_C_binding
195      implicit none
196      PetscInt start1,len1
197      PetscInt start2,len2
198      PetscReal, target ::                                                        &
199     &     array(start1:start1+len1-1,start2:start2+len2-1)
200      PetscReal, pointer :: ptr(:,:)
201
202      ptr => array
203      end subroutine
204
205      subroutine F90Array2dCreateInt(array,start1,len1,                           &
206     &     start2,len2,ptr)
207      use, intrinsic :: ISO_C_binding
208      implicit none
209      PetscInt start1,len1
210      PetscInt start2,len2
211      PetscInt, target ::                                                         &
212     &     array(start1:start1+len1-1,start2:start2+len2-1)
213      PetscInt, pointer :: ptr(:,:)
214
215      ptr => array
216      end subroutine
217
218      subroutine F90Array2dCreateFortranAddr(array,start1,len1,                   &
219     &     start2,len2,ptr)
220      use, intrinsic :: ISO_C_binding
221      implicit none
222      PetscInt start1,len1
223      PetscInt start2,len2
224      PetscFortranAddr, target ::                                                 &
225     &     array(start1:start1+len1-1,start2:start2+len2-1)
226      PetscFortranAddr, pointer :: ptr(:,:)
227
228      ptr => array
229      end subroutine
230
231!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
232      subroutine F90Array2dAccessScalar(ptr,address)
233      use, intrinsic :: ISO_C_binding
234      implicit none
235      PetscScalar, pointer :: ptr(:,:)
236      PetscFortranAddr address
237      PetscInt start1,start2
238
239      start1 = lbound(ptr,1)
240      start2 = lbound(ptr,2)
241      call F90Array2dGetAddrScalar(ptr(start1,start2),address)
242      end subroutine
243
244      subroutine F90Array2dAccessReal(ptr,address)
245      use, intrinsic :: ISO_C_binding
246      implicit none
247      PetscReal, pointer :: ptr(:,:)
248      PetscFortranAddr address
249      PetscInt start1,start2
250
251      start1 = lbound(ptr,1)
252      start2 = lbound(ptr,2)
253      call F90Array2dGetAddrReal(ptr(start1,start2),address)
254      end subroutine
255
256      subroutine F90Array2dAccessInt(ptr,address)
257      use, intrinsic :: ISO_C_binding
258      implicit none
259      PetscInt, pointer :: ptr(:,:)
260      PetscFortranAddr address
261      PetscInt start1,start2
262
263      start1 = lbound(ptr,1)
264      start2 = lbound(ptr,2)
265      call F90Array2dGetAddrInt(ptr(start1,start2),address)
266      end subroutine
267
268      subroutine F90Array2dAccessFortranAddr(ptr,address)
269      use, intrinsic :: ISO_C_binding
270      implicit none
271      PetscFortranAddr, pointer :: ptr(:,:)
272      PetscFortranAddr address
273      PetscInt start1,start2
274
275      start1 = lbound(ptr,1)
276      start2 = lbound(ptr,2)
277      call F90Array2dGetAddrFortranAddr(ptr(start1,start2),address)
278      end subroutine
279
280!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
281      subroutine F90Array2dDestroyScalar(ptr)
282      use, intrinsic :: ISO_C_binding
283      implicit none
284      PetscScalar, pointer :: ptr(:,:)
285
286      nullify(ptr)
287      end subroutine
288
289      subroutine F90Array2dDestroyReal(ptr)
290      use, intrinsic :: ISO_C_binding
291      implicit none
292      PetscReal, pointer :: ptr(:,:)
293
294      nullify(ptr)
295      end subroutine
296
297      subroutine F90Array2dDestroyInt(ptr)
298      use, intrinsic :: ISO_C_binding
299      implicit none
300      PetscInt, pointer :: ptr(:,:)
301
302      nullify(ptr)
303      end subroutine
304
305      subroutine F90Array2dDestroyFortranAddr(ptr)
306      use, intrinsic :: ISO_C_binding
307      implicit none
308      PetscFortranAddr, pointer :: ptr(:,:)
309
310      nullify(ptr)
311      end subroutine
312!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
313!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
314!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
315      subroutine F90Array3dCreateScalar(array,start1,len1,                        &
316     &     start2,len2,start3,len3,ptr)
317      use, intrinsic :: ISO_C_binding
318      implicit none
319      PetscInt start1,len1
320      PetscInt start2,len2
321      PetscInt start3,len3
322      PetscScalar, target ::                                                      &
323     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
324     &           start3:start3+len3-1)
325      PetscScalar, pointer :: ptr(:,:,:)
326
327      ptr => array
328      end subroutine
329
330      subroutine F90Array3dCreateReal(array,start1,len1,                          &
331     &     start2,len2,start3,len3,ptr)
332      use, intrinsic :: ISO_C_binding
333      implicit none
334      PetscInt start1,len1
335      PetscInt start2,len2
336      PetscInt start3,len3
337      PetscReal, target ::                                                        &
338     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
339     &           start3:start3+len3-1)
340      PetscReal, pointer :: ptr(:,:,:)
341
342      ptr => array
343      end subroutine
344
345      subroutine F90Array3dCreateInt(array,start1,len1,                           &
346     &     start2,len2,start3,len3,ptr)
347      use, intrinsic :: ISO_C_binding
348      implicit none
349      PetscInt start1,len1
350      PetscInt start2,len2
351      PetscInt start3,len3
352      PetscInt, target ::                                                         &
353     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
354     &           start3:start3+len3-1)
355      PetscInt, pointer :: ptr(:,:,:)
356
357      ptr => array
358      end subroutine
359
360      subroutine F90Array3dCreateFortranAddr(array,start1,len1,                   &
361     &     start2,len2,start3,len3,ptr)
362      use, intrinsic :: ISO_C_binding
363      implicit none
364      PetscInt start1,len1
365      PetscInt start2,len2
366      PetscInt start3,len3
367      PetscFortranAddr, target ::                                                 &
368     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
369     &           start3:start3+len3-1)
370      PetscFortranAddr, pointer :: ptr(:,:,:)
371
372      ptr => array
373      end subroutine
374
375!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
376      subroutine F90Array3dAccessScalar(ptr,address)
377      use, intrinsic :: ISO_C_binding
378      implicit none
379      PetscScalar, pointer :: ptr(:,:,:)
380      PetscFortranAddr address
381      PetscInt start1,start2,start3
382
383      start1 = lbound(ptr,1)
384      start2 = lbound(ptr,2)
385      start3 = lbound(ptr,3)
386      call F90Array3dGetAddrScalar(ptr(start1,start2,start3),address)
387      end subroutine
388
389      subroutine F90Array3dAccessReal(ptr,address)
390      use, intrinsic :: ISO_C_binding
391      implicit none
392      PetscReal, pointer :: ptr(:,:,:)
393      PetscFortranAddr address
394      PetscInt start1,start2,start3
395
396      start1 = lbound(ptr,1)
397      start2 = lbound(ptr,2)
398      start3 = lbound(ptr,3)
399      call F90Array3dGetAddrReal(ptr(start1,start2,start3),address)
400      end subroutine
401
402      subroutine F90Array3dAccessInt(ptr,address)
403      use, intrinsic :: ISO_C_binding
404      implicit none
405      PetscInt, pointer :: ptr(:,:,:)
406      PetscFortranAddr address
407      PetscInt start1,start2,start3
408
409      start1 = lbound(ptr,1)
410      start2 = lbound(ptr,2)
411      start3 = lbound(ptr,3)
412      call F90Array3dGetAddrInt(ptr(start1,start2,start3),address)
413      end subroutine
414
415      subroutine F90Array3dAccessFortranAddr(ptr,address)
416      use, intrinsic :: ISO_C_binding
417      implicit none
418      PetscFortranAddr, pointer :: ptr(:,:,:)
419      PetscFortranAddr address
420      PetscInt start1,start2,start3
421
422      start1 = lbound(ptr,1)
423      start2 = lbound(ptr,2)
424      start3 = lbound(ptr,3)
425      call F90Array3dGetAddrFortranAddr(ptr(start1,start2,start3),        &
426     &                                  address)
427      end subroutine
428
429!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
430      subroutine F90Array3dDestroyScalar(ptr)
431      use, intrinsic :: ISO_C_binding
432      implicit none
433      PetscScalar, pointer :: ptr(:,:,:)
434
435      nullify(ptr)
436      end subroutine
437
438      subroutine F90Array3dDestroyReal(ptr)
439      use, intrinsic :: ISO_C_binding
440      implicit none
441      PetscReal, pointer :: ptr(:,:,:)
442
443      nullify(ptr)
444      end subroutine
445
446      subroutine F90Array3dDestroyInt(ptr)
447      use, intrinsic :: ISO_C_binding
448      implicit none
449      PetscInt, pointer :: ptr(:,:,:)
450
451      nullify(ptr)
452      end subroutine
453
454      subroutine F90Array3dDestroyFortranAddr(ptr)
455      use, intrinsic :: ISO_C_binding
456      implicit none
457      PetscFortranAddr, pointer :: ptr(:,:,:)
458
459      nullify(ptr)
460      end subroutine
461
462!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
463      subroutine F90Array4dCreateScalar(array,start1,len1,                        &
464     &     start2,len2,start3,len3,start4,len4,ptr)
465      use, intrinsic :: ISO_C_binding
466      implicit none
467      PetscInt start1,len1
468      PetscInt start2,len2
469      PetscInt start3,len3
470      PetscInt start4,len4
471      PetscScalar, target ::                                                      &
472     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
473     &           start3:start3+len3-1,start4:start4+len4-1)
474      PetscScalar, pointer :: ptr(:,:,:,:)
475
476      ptr => array
477      end subroutine
478
479      subroutine F90Array4dCreateReal(array,start1,len1,                          &
480     &     start2,len2,start3,len3,start4,len4,ptr)
481      use, intrinsic :: ISO_C_binding
482      implicit none
483      PetscInt start1,len1
484      PetscInt start2,len2
485      PetscInt start3,len3
486      PetscInt start4,len4
487      PetscReal, target ::                                                        &
488     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
489     &           start3:start3+len3-1,start4:start4+len4-1)
490      PetscReal, pointer :: ptr(:,:,:,:)
491
492      ptr => array
493      end subroutine
494
495      subroutine F90Array4dCreateInt(array,start1,len1,                           &
496     &     start2,len2,start3,len3,start4,len4,ptr)
497      use, intrinsic :: ISO_C_binding
498      implicit none
499      PetscInt start1,len1
500      PetscInt start2,len2
501      PetscInt start3,len3
502      PetscInt start4,len4
503      PetscInt, target ::                                                         &
504     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
505     &           start3:start3+len3-1,start4:start4+len4-1)
506      PetscInt, pointer :: ptr(:,:,:,:)
507
508      ptr => array
509      end subroutine
510
511      subroutine F90Array4dCreateFortranAddr(array,start1,len1,                   &
512     &     start2,len2,start3,len3,start4,len4,ptr)
513      use, intrinsic :: ISO_C_binding
514      implicit none
515      PetscInt start1,len1
516      PetscInt start2,len2
517      PetscInt start3,len3
518      PetscInt start4,len4
519      PetscFortranAddr, target ::                                                 &
520     &     array(start1:start1+len1-1,start2:start2+len2-1,                       &
521     &           start3:start3+len3-1,start4:start4+len4-1)
522      PetscFortranAddr, pointer :: ptr(:,:,:,:)
523
524      ptr => array
525      end subroutine
526
527      subroutine F90Array4dAccessScalar(ptr,address)
528      use, intrinsic :: ISO_C_binding
529      implicit none
530      PetscScalar, pointer :: ptr(:,:,:,:)
531      PetscFortranAddr address
532      PetscInt start1,start2,start3,start4
533
534      start1 = lbound(ptr,1)
535      start2 = lbound(ptr,2)
536      start3 = lbound(ptr,3)
537      start4 = lbound(ptr,4)
538      call F90Array4dGetAddrScalar(ptr(start1,start2,start3,start4),              &
539     &                             address)
540      end subroutine
541
542      subroutine F90Array4dAccessReal(ptr,address)
543      use, intrinsic :: ISO_C_binding
544      implicit none
545      PetscReal, pointer :: ptr(:,:,:,:)
546      PetscFortranAddr address
547      PetscInt start1,start2,start3,start4
548
549      start1 = lbound(ptr,1)
550      start2 = lbound(ptr,2)
551      start3 = lbound(ptr,3)
552      start4 = lbound(ptr,4)
553      call F90Array4dGetAddrReal(ptr(start1,start2,start3,start4),                &
554     &                             address)
555      end subroutine
556
557      subroutine F90Array4dAccessInt(ptr,address)
558      use, intrinsic :: ISO_C_binding
559      implicit none
560      PetscInt, pointer :: ptr(:,:,:,:)
561      PetscFortranAddr address
562      PetscInt start1,start2,start3,start4
563
564      start1 = lbound(ptr,1)
565      start2 = lbound(ptr,2)
566      start3 = lbound(ptr,3)
567      start4 = lbound(ptr,4)
568      call F90Array4dGetAddrInt(ptr(start1,start2,start3,start4),                 &
569     &                             address)
570      end subroutine
571
572      subroutine F90Array4dAccessFortranAddr(ptr,address)
573      use, intrinsic :: ISO_C_binding
574      implicit none
575      PetscScalar, pointer :: ptr(:,:,:,:)
576      PetscFortranAddr address
577      PetscFortranAddr start1,start2,start3,start4
578
579      start1 = lbound(ptr,1)
580      start2 = lbound(ptr,2)
581      start3 = lbound(ptr,3)
582      start4 = lbound(ptr,4)
583      call F90Array4dGetAddrFortranAddr(ptr(start1,start2,start3,                 &
584     &                                      start4),address)
585      end subroutine
586
587      subroutine F90Array4dDestroyScalar(ptr)
588      use, intrinsic :: ISO_C_binding
589      implicit none
590      PetscScalar, pointer :: ptr(:,:,:,:)
591
592      nullify(ptr)
593      end subroutine
594
595      subroutine F90Array4dDestroyReal(ptr)
596      use, intrinsic :: ISO_C_binding
597      implicit none
598      PetscReal, pointer :: ptr(:,:,:,:)
599
600      nullify(ptr)
601      end subroutine
602
603      subroutine F90Array4dDestroyInt(ptr)
604      use, intrinsic :: ISO_C_binding
605      implicit none
606      PetscInt, pointer :: ptr(:,:,:,:)
607
608      nullify(ptr)
609      end subroutine
610
611      subroutine F90Array4dDestroyFortranAddr(ptr)
612      use, intrinsic :: ISO_C_binding
613      implicit none
614      PetscFortranAddr, pointer :: ptr(:,:,:,:)
615
616      nullify(ptr)
617      end subroutine
618
619!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
620!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX!
621!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
622