Actual source code: zisltogf90.c

  1: #include <petscis.h>
  2: #include <petsc/private/f90impl.h>

  4: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  5:   #define islocaltoglobalmappinggetindicesf90_          ISLOCALTOGLOBALMAPPINGGETINDICESF90
  6:   #define islocaltoglobalmappingrestoreindicesf90_      ISLOCALTOGLOBALMAPPINGRESTOREINDICESF90
  7:   #define islocaltoglobalmappinggetblockindicesf90_     ISLOCALTOGLOBALMAPPINGGETBLOCKINDICESF90
  8:   #define islocaltoglobalmappingrestoreblockindicesf90_ ISLOCALTOGLOBALMAPPINGRESTOREBLOCKINDICESF90
  9: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 10:   #define islocaltoglobalmappinggetindicesf90_          islocaltoglobalmappinggetindicesf90
 11:   #define islocaltoglobalmappingrestoreindicesf90_      islocaltoglobalmappingrestoreindicesf90
 12:   #define islocaltoglobalmappinggetblockindicesf90_     islocaltoglobalmappinggetblockindicesf90
 13:   #define islocaltoglobalmappingrestoreblockindicesf90_ islocaltoglobalmappingrestoreblockindicesf90
 14: #endif

 16: PETSC_EXTERN void islocaltoglobalmappinggetindicesf90_(ISLocalToGlobalMapping *da, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 17: {
 18:   const PetscInt *idx;
 19:   PetscInt        n;
 20:   *ierr = ISLocalToGlobalMappingGetIndices(*da, &idx);
 21:   if (*ierr) return;
 22:   *ierr = ISLocalToGlobalMappingGetSize(*da, &n);
 23:   if (*ierr) return;
 24:   *ierr = F90Array1dCreate((void *)idx, MPIU_INT, 1, n, indices PETSC_F90_2PTR_PARAM(ptrd));
 25: }

 27: PETSC_EXTERN void islocaltoglobalmappingrestoreindicesf90_(ISLocalToGlobalMapping *da, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 28: {
 29:   const PetscInt *fa;

 31:   *ierr = F90Array1dAccess(ptr, MPIU_INT, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
 32:   if (*ierr) return;
 33:   *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
 34:   if (*ierr) return;
 35:   *ierr = ISLocalToGlobalMappingRestoreIndices(*da, &fa);
 36:   if (*ierr) return;
 37: }

 39: PETSC_EXTERN void islocaltoglobalmappinggetblockindicesf90_(ISLocalToGlobalMapping *da, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 40: {
 41:   const PetscInt *idx;
 42:   PetscInt        n;
 43:   *ierr = ISLocalToGlobalMappingGetBlockIndices(*da, &idx);
 44:   if (*ierr) return;
 45:   *ierr = ISLocalToGlobalMappingGetSize(*da, &n);
 46:   if (*ierr) return;
 47:   *ierr = F90Array1dCreate((void *)idx, MPIU_INT, 1, n, indices PETSC_F90_2PTR_PARAM(ptrd));
 48: }

 50: PETSC_EXTERN void islocaltoglobalmappingrestoreblockindicesf90_(ISLocalToGlobalMapping *da, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 51: {
 52:   const PetscInt *fa;

 54:   *ierr = F90Array1dAccess(ptr, MPIU_INT, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
 55:   if (*ierr) return;
 56:   *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
 57:   if (*ierr) return;
 58:   *ierr = ISLocalToGlobalMappingRestoreBlockIndices(*da, &fa);
 59:   if (*ierr) return;
 60: }