Actual source code: zindexf.c
1: #include <petsc/private/fortranimpl.h>
2: #include <petscis.h>
3: #include <petscviewer.h>
5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
6: #define isgetindices_ ISGETINDICES
7: #define isrestoreindices_ ISRESTOREINDICES
8: #define isgettotalindices_ ISGETTOTALINDICES
9: #define isrestoretotalindices_ ISRESTORETOTALINDICES
10: #define isgetnonlocalindices_ ISGETNONLOCALINDICES
11: #define isrestorenonlocalindices_ ISRESTORENONLOCALINDICES
12: #define islocaltoglobalmappinggetindices_ ISLOCALTOGLOBALMAPPINGGETINDICES
13: #define islocaltoglobalmappingrestoreindices_ ISLOCALTOGLOBALMAPPINGRESTOREINDICES
14: #define islocaltoglobalmappinggetblockindices_ ISLOCALTOGLOBALMAPPINGGETBLOCKINDICES
15: #define islocaltoglobalmappingrestoreblockindices_ ISLOCALTOGLOBALMAPPINGRESTOREBLOCKINDICES
16: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
17: #define isgetindices_ isgetindices
18: #define isrestoreindices_ isrestoreindices
19: #define isgettotalindices_ isgettotalindices
20: #define isrestoretotalindices_ isrestoretotalindices
21: #define isgetnonlocalindices_ isgetnonlocalindices
22: #define isrestorenonlocalindices_ isrestorenonlocalindices
23: #define islocaltoglobalmappinggetindices_ islocaltoglobalmappinggetindices
24: #define islocaltoglobalmappingrestoreindices_ islocaltoglobalmappingrestoreindices
25: #define islocaltoglobalmappinggetblockindices_ islocaltoglobalmappinggetblockindices
26: #define islocaltoglobalmappingrestoreblockindices_ islocaltoglobalmappingrestoreblockindices
27: #endif
29: PETSC_EXTERN void isgetindices_(IS *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
30: {
31: const PetscInt *lx;
33: *ierr = ISGetIndices(*x, &lx);
34: if (*ierr) return;
35: *ia = PetscIntAddressToFortran(fa, (PetscInt *)lx);
36: }
38: PETSC_EXTERN void isrestoreindices_(IS *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
39: {
40: const PetscInt *lx = PetscIntAddressFromFortran(fa, *ia);
41: *ierr = ISRestoreIndices(*x, &lx);
42: }
44: PETSC_EXTERN void isgettotalindices_(IS *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
45: {
46: const PetscInt *lx;
48: *ierr = ISGetTotalIndices(*x, &lx);
49: if (*ierr) return;
50: *ia = PetscIntAddressToFortran(fa, (PetscInt *)lx);
51: }
53: PETSC_EXTERN void isrestoretotalindices_(IS *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
54: {
55: const PetscInt *lx = PetscIntAddressFromFortran(fa, *ia);
56: *ierr = ISRestoreTotalIndices(*x, &lx);
57: }
59: PETSC_EXTERN void isgetnonlocalindices_(IS *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
60: {
61: const PetscInt *lx;
63: *ierr = ISGetNonlocalIndices(*x, &lx);
64: if (*ierr) return;
65: *ia = PetscIntAddressToFortran(fa, (PetscInt *)lx);
66: }
68: PETSC_EXTERN void isrestorenonlocalindices_(IS *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
69: {
70: const PetscInt *lx = PetscIntAddressFromFortran(fa, *ia);
71: *ierr = ISRestoreNonlocalIndices(*x, &lx);
72: }
74: PETSC_EXTERN void islocaltoglobalmappinggetindices_(ISLocalToGlobalMapping *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
75: {
76: const PetscInt *lx;
78: *ierr = ISLocalToGlobalMappingGetIndices(*x, &lx);
79: if (*ierr) return;
80: *ia = PetscIntAddressToFortran(fa, (PetscInt *)lx);
81: }
83: PETSC_EXTERN void islocaltoglobalmappingrestoreindices_(ISLocalToGlobalMapping *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
84: {
85: const PetscInt *lx = PetscIntAddressFromFortran(fa, *ia);
86: *ierr = ISLocalToGlobalMappingRestoreIndices(*x, &lx);
87: }
89: PETSC_EXTERN void islocaltoglobalmappinggetblockindices_(ISLocalToGlobalMapping *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
90: {
91: const PetscInt *lx;
93: *ierr = ISLocalToGlobalMappingGetBlockIndices(*x, &lx);
94: if (*ierr) return;
95: *ia = PetscIntAddressToFortran(fa, (PetscInt *)lx);
96: }
98: PETSC_EXTERN void islocaltoglobalmappingrestoreblockindices_(ISLocalToGlobalMapping *x, PetscInt *fa, size_t *ia, PetscErrorCode *ierr)
99: {
100: const PetscInt *lx = PetscIntAddressFromFortran(fa, *ia);
101: *ierr = ISLocalToGlobalMappingRestoreBlockIndices(*x, &lx);
102: }