Actual source code: zhdf5f.c

  1: #include <petsc/private/ftnimpl.h>
  2: #include <petscviewerhdf5.h>

  4: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  5:   #define petscviewerhdf5opengroup_            PETSCVIEWERHDF5OPENGROUP
  6:   #define petscviewerhdf5writeattributeint_    PETSCVIEWERHDF5WRITEATTRIBUTEINT
  7:   #define petscviewerhdf5writeattributescalar_ PETSCVIEWERHDF5WRITEATTRIBUTESCALAR
  8:   #define petscviewerhdf5writeattributereal_   PETSCVIEWERHDF5WRITEATTRIBUTEREAL
  9: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 10:   #define petscviewerhdf5opengroup_            petscviewerhdf5opengroup
 11:   #define petscviewerhdf5writeattributeint_    petscviewerhdf5writeattributeint
 12:   #define petscviewerhdf5writeattributescalar_ petscviewerhdf5writeattributescalar
 13:   #define petscviewerhdf5writeattributereal_   petscviewerhdf5writeattributereal
 14: #endif

 16: PETSC_EXTERN void petscviewerhdf5opengroup_(PetscViewer *viewer, char path[], hid_t *fileId, hid_t *groupId, int *ierr, PETSC_FORTRAN_CHARLEN_T len)
 17: {
 18:   char *c1;

 20:   FIXCHAR(path, len, c1);
 21:   *ierr = PetscViewerHDF5OpenGroup(*viewer, c1, fileId, groupId);
 22:   FREECHAR(path, c1);
 23: }

 25: PETSC_EXTERN void petscviewerhdf5writeattributeint_(PetscViewer *viewer, const char parent[], const char name[], PetscInt *value, int *ierr, PETSC_FORTRAN_CHARLEN_T len1, PETSC_FORTRAN_CHARLEN_T len2)
 26: {
 27:   char *c1;
 28:   char *c2;

 30:   FIXCHAR(parent, len1, c1);
 31:   FIXCHAR(name, len2, c2);
 32:   *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_INT, value);
 33:   FREECHAR(parent, c1);
 34:   FREECHAR(name, c2);
 35: }

 37: PETSC_EXTERN void petscviewerhdf5writeattributescalar_(PetscViewer *viewer, const char parent[], const char name[], PetscScalar *value, int *ierr, PETSC_FORTRAN_CHARLEN_T len1, PETSC_FORTRAN_CHARLEN_T len2)
 38: {
 39:   char *c1;
 40:   char *c2;

 42:   FIXCHAR(parent, len1, c1);
 43:   FIXCHAR(name, len2, c2);
 44:   *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_SCALAR, value);
 45:   FREECHAR(parent, c1);
 46:   FREECHAR(name, c2);
 47: }

 49: PETSC_EXTERN void petscviewerhdf5writeattributereal_(PetscViewer *viewer, const char parent[], const char name[], PetscReal *value, int *ierr, PETSC_FORTRAN_CHARLEN_T len1, PETSC_FORTRAN_CHARLEN_T len2)
 50: {
 51:   char *c1;
 52:   char *c2;

 54:   FIXCHAR(parent, len1, c1);
 55:   FIXCHAR(name, len2, c2);
 56:   *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_REAL, value);
 57:   FREECHAR(parent, c1);
 58:   FREECHAR(name, c2);
 59: }