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