Actual source code: zvsectionf90.c

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

  5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  6:   #define vecsetvaluessectionf90_     VECSETVALUESSECTIONF90
  7:   #define vecgetvaluessectionf90_     VECGETVALUESSECTIONF90
  8:   #define vecrestorevaluessectionf90_ VECRESTOREVALUESSECTIONF90
  9: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) && !defined(FORTRANDOUBLEUNDERSCORE)
 10:   #define vecsetvaluessectionf90_     vecsetvaluessectionf90
 11:   #define vecgetvaluessectionf90_     vecgetvaluessectionf90
 12:   #define vecrestorevaluessectionf90_ vecrestorevaluessectionf90
 13: #endif

 15: /* Definitions of Fortran Wrapper routines */

 17: PETSC_EXTERN void vecsetvaluessectionf90_(Vec *v, PetscSection *section, PetscInt *point, F90Array1d *ptr, InsertMode *mode, int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
 18: {
 19:   PetscScalar *array;

 21:   *__ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&array PETSC_F90_2PTR_PARAM(ptrd));
 22:   if (*__ierr) return;
 23:   *__ierr = VecSetValuesSection(*v, *section, *point, array, *mode);
 24: }

 26: PETSC_EXTERN void vecgetvaluessectionf90_(Vec *v, PetscSection *section, PetscInt *point, F90Array1d *ptr, int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
 27: {
 28:   PetscScalar *fa;
 29:   PetscInt     len;

 31:   *__ierr = VecGetValuesSection(*v, *section, *point, &fa);
 32:   if (*__ierr) return;
 33:   *__ierr = PetscSectionGetDof(*section, *point, &len);
 34:   if (*__ierr) return;
 35:   *__ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, len, ptr PETSC_F90_2PTR_PARAM(ptrd));
 36: }

 38: PETSC_EXTERN void vecrestorevaluessectionf90_(Vec *v, PetscSection *section, PetscInt *point, F90Array1d *ptr, int *__ierr PETSC_F90_2PTR_PROTO(ptrd))
 39: {
 40:   PetscScalar *fa;

 42:   *__ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
 43:   if (*__ierr) return;
 44:   *__ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
 45:   if (*__ierr) return;
 46: }