xref: /petsc/src/vec/ftn-mod/petscvec.h90 (revision a4e01a8b1ea6aa127f3e22aee39c726cf71ba8b6)
1#if defined(PETSC_HAVE_FORTRAN_TYPE_STAR)
2      interface
3        subroutine PetscSFBcastBegin(sf, unit, rarray, larray, op, ierr)
4          use, intrinsic :: ISO_C_binding
5#if defined(PETSC_USE_MPI_F08)
6          import tPetscSF, MPI_Datatype, MPI_Op
7#else
8          import tPetscSF
9#endif
10          PetscSF :: sf
11          MPIU_Datatype :: unit
12          MPIU_Op :: op
13          type(*) :: rarray(:)
14          type(*) :: larray(:)
15          PetscErrorCode :: ierr
16        end subroutine PetscSFBcastBegin
17
18        subroutine PetscSFBcastEnd(sf, unit, rarray, larray, op, ierr)
19          use, intrinsic :: ISO_C_binding
20#if defined(PETSC_USE_MPI_F08)
21          import tPetscSF, MPI_Datatype, MPI_Op
22#else
23          import tPetscSF
24#endif
25          PetscSF :: sf
26          MPIU_Datatype :: unit
27          MPIU_Op :: op
28          type(*) :: rarray(:)
29          type(*) :: larray(:)
30          PetscErrorCode :: ierr
31        end subroutine PetscSFBcastEnd
32
33        subroutine PetscSFReduceBegin(sf, unit, larray, rarray, op, ierr)
34          use, intrinsic :: ISO_C_binding
35#if defined(PETSC_USE_MPI_F08)
36          import tPetscSF, MPI_Datatype, MPI_Op
37#else
38          import tPetscSF
39#endif
40          PetscSF :: sf
41          MPIU_Datatype :: unit
42          MPIU_Op :: op
43          type(*) :: larray(:)
44          type(*) :: rarray(:)
45          PetscErrorCode :: ierr
46        end subroutine PetscSFReduceBegin
47
48        subroutine PetscSFReduceEnd(sf, unit, larray, rarray, op, ierr)
49          use, intrinsic :: ISO_C_binding
50#if defined(PETSC_USE_MPI_F08)
51          import tPetscSF, MPI_Datatype, MPI_Op
52#else
53          import tPetscSF
54#endif
55          PetscSF :: sf
56          MPIU_Datatype :: unit
57          MPIU_Op :: op
58          type(*) :: larray(:)
59          type(*) :: rarray(:)
60          PetscErrorCode :: ierr
61        end subroutine PetscSFReduceEnd
62      end interface
63#endif
64
65      interface
66        subroutine VecRestoreOwnershipRanges(v, ptr, ierr)
67          use, intrinsic :: ISO_C_binding
68          import tVec
69          Vec :: v
70          PetscInt, pointer :: ptr(:)
71          PetscErrorCode, intent(out) :: ierr
72        end subroutine VecRestoreOwnershipRanges
73
74        subroutine PetscSFRestoreGraph(sf, nroots, nleaves, ilocal, iremote, ierr)
75          use, intrinsic :: ISO_C_binding
76          import tPetscSF, sPetscSFNode
77          PetscSF :: sf
78          PetscInt :: nroots, nleaves
79          PetscInt, pointer :: ilocal(:)
80          type(sPetscSFNode), pointer :: iremote(:)
81          PetscErrorCode :: ierr
82        end subroutine PetscSFRestoreGraph
83
84        subroutine VecRestoreValuesSection(v, s, p, va, ierr)
85          use, intrinsic :: ISO_C_binding
86          import tVec, tPetscSection
87          PetscScalar, pointer :: va(:)
88          PetscErrorCode ierr
89          Vec v
90          PetscSection s
91          PetscInt p
92        end subroutine
93      end interface
94