Actual source code: zsfutilsf90.c

  1: #include <petsc/private/f90impl.h>
  2: #include <petscsf.h>
  3: #include <petscsection.h>

  5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  6:   #define petscsfdistributesectionf90_   PETSCSFDISTRIBUTESECTIONF90
  7:   #define petscsfcreatesectionsff90_     PETSCSFCREATESECTIONSFF90
  8:   #define petscsfcreateremoteoffsetsf90_ PETSCSFCREATEREMOTEOFFSETSF90
  9: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 10:   #define petscsfdistributesectionf90_   petscsfdistributesectionf90
 11:   #define petscsfcreatesectionsff90_     petscsfcreatesectionsff90
 12:   #define petscsfcreateremoteoffsetsf90_ petscsfcreateremoteoffsetsf90
 13: #endif

 15: PETSC_EXTERN void petscsfdistributesectionf90_(PetscSF *sf, PetscSection *rootSection, F90Array1d *ptr, PetscSection *leafSection, int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
 16: {
 17:   PetscInt *fa;

 19:   *__ierr = F90Array1dAccess(ptr, MPIU_INT, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
 20:   if (*__ierr) return;
 21:   *__ierr = PetscSFDistributeSection(*sf, *rootSection, &fa, *leafSection);
 22: }

 24: PETSC_EXTERN void petscsfcreatesectionsff90_(PetscSF *pointSF, PetscSection *rootSection, F90Array1d *ptr, PetscSection *leafSection, PetscSF *sf, int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
 25: {
 26:   PetscInt *fa;

 28:   *__ierr = F90Array1dAccess(ptr, MPIU_INT, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
 29:   if (*__ierr) return;
 30:   *__ierr = PetscSFCreateSectionSF(*pointSF, *rootSection, fa, *leafSection, sf);
 31: }

 33: PETSC_EXTERN void petscsfcreateremoteoffsetsf90_(PetscSF *pointSF, PetscSection *rootSection, PetscSection *leafSection, F90Array1d *ptr, int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
 34: {
 35:   PetscInt *fa;
 36:   PetscInt  lpStart, lpEnd;

 38:   *__ierr = PetscSFCreateRemoteOffsets(*pointSF, *rootSection, *leafSection, &fa);
 39:   *__ierr = PetscSectionGetChart(*leafSection, &lpStart, &lpEnd);
 40:   if (*__ierr) return;
 41:   *__ierr = F90Array1dCreate((void *)fa, MPIU_INT, 1, lpEnd - lpStart, ptr PETSC_F90_2PTR_PARAM(ptrd));
 42: }