Actual source code: zvsectionisf90.c

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

  4: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  5:   #define petscsectiongetconstraintindicesf90_          PETSCSECTIONGETCONSTRAINTINDICESF90
  6:   #define petscsectionrestoreconstraintindicesf90_      PETSCSECTIONRESTORECONSTRAINTINDICESF90
  7:   #define petscsectionsetconstraintindicesf90_          PETSCSECTIONSETCONSTRAINTINDICESF90
  8:   #define petscsectiongetfieldconstraintindicesf90_     PETSCSECTIONGETFIELDINDICESF90
  9:   #define petscsectionrestorefieldconstraintindicesf90_ PETSCSECTIONRESTOREFIELDCONSTRAINTINDICESF90
 10:   #define petscsectionsetfieldconstraintindicesf90_     PETSCSECTIONSETFIELDCONSTRAINTINDICESF90
 11: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 12:   #define petscsectiongetconstraintindicesf90_          petscsectiongetconstraintindicesf90
 13:   #define petscsectionrestoreconstraintindicesf90_      petscsectionrestoreconstraintindicesf90
 14:   #define petscsectionsetconstraintindicesf90_          petscsectionsetconstraintindicesf90
 15:   #define petscsectiongetfieldconstraintindicesf90_     petscsectiongetfieldconstraintindicesf90
 16:   #define petscsectionrestorefieldconstraintindicesf90_ petscsectionrestorefieldconstraintindicesf90
 17:   #define petscsectionsetfieldconstraintindicesf90_     petscsectionsetfieldconstraintindicesf90
 18: #endif

 20: PETSC_EXTERN void petscsectiongetconstraintindicesf90_(PetscSection *s, PetscInt *point, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 21: {
 22:   const PetscInt *idx;
 23:   PetscInt        n;

 25:   *ierr = PetscSectionGetConstraintIndices(*s, *point, &idx);
 26:   if (*ierr) return;
 27:   *ierr = PetscSectionGetConstraintDof(*s, *point, &n);
 28:   if (*ierr) return;
 29:   *ierr = F90Array1dCreate((void *)idx, MPIU_INT, 1, n, indices PETSC_F90_2PTR_PARAM(ptrd));
 30: }

 32: PETSC_EXTERN void petscsectionrestoreconstraintindicesf90_(PetscSection *s, PetscInt *point, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 33: {
 34:   *ierr = F90Array1dDestroy(indices, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
 35:   if (*ierr) return;
 36: }

 38: PETSC_EXTERN void petscsectionsetconstraintindicesf90_(PetscSection *s, PetscInt *point, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 39: {
 40:   const PetscInt *idx;

 42:   *ierr = F90Array1dAccess(indices, MPIU_INT, (void **)&idx PETSC_F90_2PTR_PARAM(ptrd));
 43:   if (*ierr) return;
 44:   *ierr = PetscSectionSetConstraintIndices(*s, *point, idx);
 45:   if (*ierr) return;
 46: }

 48: PETSC_EXTERN void petscsectiongetfieldconstraintindicesf90_(PetscSection *s, PetscInt *point, PetscInt *field, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 49: {
 50:   const PetscInt *idx;
 51:   PetscInt        n;

 53:   *ierr = PetscSectionGetFieldConstraintIndices(*s, *point, *field, &idx);
 54:   if (*ierr) return;
 55:   *ierr = PetscSectionGetFieldConstraintDof(*s, *point, *field, &n);
 56:   if (*ierr) return;
 57:   *ierr = F90Array1dCreate((void *)idx, MPIU_INT, 1, n, indices PETSC_F90_2PTR_PARAM(ptrd));
 58: }

 60: PETSC_EXTERN void petscsectionrestorefieldconstraintindicesf90_(PetscSection *s, PetscInt *point, PetscInt *field, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 61: {
 62:   *ierr = F90Array1dDestroy(indices, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
 63:   if (*ierr) return;
 64: }

 66: PETSC_EXTERN void petscsectionsetfieldconstraintindicesf90_(PetscSection *s, PetscInt *point, PetscInt *field, F90Array1d *indices, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
 67: {
 68:   const PetscInt *idx;

 70:   *ierr = F90Array1dAccess(indices, MPIU_INT, (void **)&idx PETSC_F90_2PTR_PARAM(ptrd));
 71:   if (*ierr) return;
 72:   *ierr = PetscSectionSetFieldConstraintIndices(*s, *point, *field, idx);
 73:   if (*ierr) return;
 74: }