Actual source code: zhdf5f.c
1: #include <petsc/private/ftnimpl.h>
2: #include <petscviewerhdf5.h>
4: #if defined(PETSC_HAVE_FORTRAN_CAPS)
5: #define petscviewerhdf5writeattributeint_ PETSCVIEWERHDF5WRITEATTRIBUTEINT
6: #define petscviewerhdf5writeattributescalar_ PETSCVIEWERHDF5WRITEATTRIBUTESCALAR
7: #define petscviewerhdf5writeattributereal_ PETSCVIEWERHDF5WRITEATTRIBUTEREAL
8: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
9: #define petscviewerhdf5writeattributeint_ petscviewerhdf5writeattributeint
10: #define petscviewerhdf5writeattributescalar_ petscviewerhdf5writeattributescalar
11: #define petscviewerhdf5writeattributereal_ petscviewerhdf5writeattributereal
12: #endif
14: 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)
15: {
16: char *c1;
17: char *c2;
19: FIXCHAR(parent, len1, c1);
20: FIXCHAR(name, len2, c2);
21: *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_INT, value);
22: FREECHAR(parent, c1);
23: FREECHAR(name, c2);
24: }
26: 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)
27: {
28: char *c1;
29: char *c2;
31: FIXCHAR(parent, len1, c1);
32: FIXCHAR(name, len2, c2);
33: *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_SCALAR, value);
34: FREECHAR(parent, c1);
35: FREECHAR(name, c2);
36: }
38: 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)
39: {
40: char *c1;
41: char *c2;
43: FIXCHAR(parent, len1, c1);
44: FIXCHAR(name, len2, c2);
45: *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_REAL, value);
46: FREECHAR(parent, c1);
47: FREECHAR(name, c2);
48: }