xref: /petsc/src/sys/classes/viewer/impls/hdf5/ftn-custom/zhdf5f.c (revision a2d6be1bbe027d53c625ca74a1867ef3b55d1170)
1af0996ceSBarry Smith #include <petsc/private/fortranimpl.h>
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
55c6c1daeSBarry Smith #define petscviewerhdf5open_            PETSCVIEWERHDF5OPEN
63ef9c667SSatish Balay #define petscviewerhdf5pushgroup_       PETSCVIEWERHDF5PUSHGROUP
73ef9c667SSatish Balay #define petscviewerhdf5getgroup_        PETSCVIEWERHDF5GETGROUP
889e0ef10SVaclav Hapla #define petscviewerhdf5hasdataset_      PETSCVIEWERHDF5HASDATASET
906bbad5aSVaclav Hapla #define petscviewerhdf5hasattribute_    PETSCVIEWERHDF5HASATTRIBUTE
10df863907SAlex Fikl #define petscviewerhdf5writeattribute_  PETSCVIEWERHDF5WRITEATTRIBUTE
11df863907SAlex Fikl #define petscviewerhdf5readattribute_   PETSCVIEWERHDF5READATTRIBUTE
124302210dSVaclav Hapla #define petscviewerhdf5hasgroup_        PETSCVIEWERHDF5HASGROUP
135c6c1daeSBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
145c6c1daeSBarry Smith #define petscviewerhdf5open_            petscviewerhdf5open
153ef9c667SSatish Balay #define petscviewerhdf5pushgroup_       petscviewerhdf5pushgroup
163ef9c667SSatish Balay #define petscviewerhdf5getgroup_        petscviewerhdf5getgroup
1789e0ef10SVaclav Hapla #define petscviewerhdf5hasdataset_      petscviewerhdf5hasdataset
1806bbad5aSVaclav Hapla #define petscviewerhdf5hasattribute_    petscviewerhdf5hasattribute
19df863907SAlex Fikl #define petscviewerhdf5writeattribute_  petscviewerhdf5writeattribute
20df863907SAlex Fikl #define petscviewerhdf5readattribute_   petscviewerhdf5readattribute
214302210dSVaclav Hapla #define petscviewerhdf5hasgroup_        petscviewerhdf5hasgroup
225c6c1daeSBarry Smith #endif
235c6c1daeSBarry Smith 
2419caf8f3SSatish Balay PETSC_EXTERN void petscviewerhdf5open_(MPI_Comm *comm, char* name, PetscFileMode *type,
2519caf8f3SSatish Balay     PetscViewer *binv, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
265c6c1daeSBarry Smith {
275c6c1daeSBarry Smith   char *c1;
28df863907SAlex Fikl 
295c6c1daeSBarry Smith   FIXCHAR(name, len, c1);
30d49bb8f9SBarry Smith   *ierr = PetscViewerHDF5Open(MPI_Comm_f2c(*(MPI_Fint*)&*comm), c1, *type, binv);if (*ierr) return;
315c6c1daeSBarry Smith   FREECHAR(name, c1);
325c6c1daeSBarry Smith }
335c6c1daeSBarry Smith 
3419caf8f3SSatish Balay PETSC_EXTERN void petscviewerhdf5pushgroup_(PetscViewer *viewer, char* name,
3519caf8f3SSatish Balay     PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
363ef9c667SSatish Balay {
373ef9c667SSatish Balay   char *c1;
38df863907SAlex Fikl 
393ef9c667SSatish Balay   FIXCHAR(name, len, c1);
40d49bb8f9SBarry Smith   *ierr = PetscViewerHDF5PushGroup(*viewer, c1);if (*ierr) return;
413ef9c667SSatish Balay   FREECHAR(name, c1);
423ef9c667SSatish Balay }
433ef9c667SSatish Balay 
4419caf8f3SSatish Balay PETSC_EXTERN void petscviewerhdf5getgroup_(PetscViewer *viewer, char* name,
4519caf8f3SSatish Balay     PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
463ef9c667SSatish Balay {
473ef9c667SSatish Balay   const char *c1;
48df863907SAlex Fikl 
49d49bb8f9SBarry Smith   *ierr = PetscViewerHDF5GetGroup(*viewer, &c1);if (*ierr) return;
503ef9c667SSatish Balay   *ierr = PetscStrncpy(name, c1, len);
513ef9c667SSatish Balay   FIXRETURNCHAR(PETSC_TRUE,name,len);
523ef9c667SSatish Balay }
53df863907SAlex Fikl 
5419caf8f3SSatish Balay PETSC_EXTERN void petscviewerhdf5hasattribute_(PetscViewer *viewer, char* parent,
5519caf8f3SSatish Balay     char* name, PetscBool *has, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T plen,PETSC_FORTRAN_CHARLEN_T nlen)
56df863907SAlex Fikl {
57df863907SAlex Fikl    char *c1, *c2;
58df863907SAlex Fikl 
59df863907SAlex Fikl    FIXCHAR(parent, plen, c1);
60df863907SAlex Fikl    FIXCHAR(name, nlen, c2);
61d49bb8f9SBarry Smith    *ierr = PetscViewerHDF5HasAttribute(*viewer, c1, c2, has);if (*ierr) return;
62df863907SAlex Fikl    FREECHAR(parent, c1);
63df863907SAlex Fikl    FREECHAR(name, c2);
64df863907SAlex Fikl }
65df863907SAlex Fikl 
6619caf8f3SSatish Balay PETSC_EXTERN void petscviewerhdf5writeattribute_(PetscViewer *viewer, char* parent,
6719caf8f3SSatish Balay     char* name, PetscDataType *datatype, const void *value, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T plen,PETSC_FORTRAN_CHARLEN_T nlen)
68df863907SAlex Fikl {
69df863907SAlex Fikl    char *c1, *c2;
70df863907SAlex Fikl 
71df863907SAlex Fikl    FIXCHAR(parent, plen, c1);
72df863907SAlex Fikl    FIXCHAR(name, nlen, c2);
73*a2d6be1bSVaclav Hapla    *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, *datatype, value);if (*ierr) return;
74df863907SAlex Fikl    FREECHAR(parent, c1);
75df863907SAlex Fikl    FREECHAR(name, c2);
76df863907SAlex Fikl }
77df863907SAlex Fikl 
7819caf8f3SSatish Balay PETSC_EXTERN void petscviewerhdf5readattribute_(PetscViewer *viewer, char* parent,
79*a2d6be1bSVaclav Hapla     char* name, PetscDataType *datatype, void *defaultValue, void *value, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T plen,PETSC_FORTRAN_CHARLEN_T nlen)
80df863907SAlex Fikl {
81df863907SAlex Fikl    char *c1, *c2;
82df863907SAlex Fikl 
83df863907SAlex Fikl    FIXCHAR(parent, plen, c1);
84df863907SAlex Fikl    FIXCHAR(name, nlen, c2);
85*a2d6be1bSVaclav Hapla    *ierr = PetscViewerHDF5ReadAttribute(*viewer, c1, c2, *datatype, defaultValue, value);if (*ierr) return;
86df863907SAlex Fikl    FREECHAR(parent, c1);
87df863907SAlex Fikl    FREECHAR(name, c2);
88df863907SAlex Fikl }
894302210dSVaclav Hapla 
9089e0ef10SVaclav Hapla PETSC_EXTERN void petscviewerhdf5hasdataset_(PetscViewer *viewer, char* path, PetscBool *has, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
9189e0ef10SVaclav Hapla {
9289e0ef10SVaclav Hapla    char *c1;
9389e0ef10SVaclav Hapla 
9489e0ef10SVaclav Hapla    FIXCHAR(path, len, c1);
9589e0ef10SVaclav Hapla    *ierr = PetscViewerHDF5HasDataset(*viewer, c1, has);if (*ierr) return;
9689e0ef10SVaclav Hapla    FREECHAR(path, c1);
9789e0ef10SVaclav Hapla }
9889e0ef10SVaclav Hapla 
994302210dSVaclav Hapla PETSC_EXTERN void petscviewerhdf5hasgroup_(PetscViewer *viewer, char* path, PetscBool *has, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
1004302210dSVaclav Hapla {
1014302210dSVaclav Hapla    char *c1;
1024302210dSVaclav Hapla 
1034302210dSVaclav Hapla    FIXCHAR(path, len, c1);
1044302210dSVaclav Hapla    *ierr = PetscViewerHDF5HasGroup(*viewer, c1, has);if (*ierr) return;
1054302210dSVaclav Hapla    FREECHAR(path, c1);
1064302210dSVaclav Hapla }
107