xref: /petsc/src/sys/classes/viewer/impls/hdf5/ftn-custom/zhdf5f.c (revision c6a13b3964cb660c2492050346160b55652290a0)
16dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
53014b61aSVaclav Hapla   #define petscviewerhdf5opengroup_            PETSCVIEWERHDF5OPENGROUP
6*a968899dSTapashree Pradhan   #define petscviewerhdf5writeattributeint_    PETSCVIEWERHDF5WRITEATTRIBUTEINT
7*a968899dSTapashree Pradhan   #define petscviewerhdf5writeattributescalar_ PETSCVIEWERHDF5WRITEATTRIBUTESCALAR
8*a968899dSTapashree Pradhan   #define petscviewerhdf5writeattributereal_   PETSCVIEWERHDF5WRITEATTRIBUTEREAL
95c6c1daeSBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
103014b61aSVaclav Hapla   #define petscviewerhdf5opengroup_            petscviewerhdf5opengroup
11*a968899dSTapashree Pradhan   #define petscviewerhdf5writeattributeint_    petscviewerhdf5writeattributeint
12*a968899dSTapashree Pradhan   #define petscviewerhdf5writeattributescalar_ petscviewerhdf5writeattributescalar
13*a968899dSTapashree Pradhan   #define petscviewerhdf5writeattributereal_   petscviewerhdf5writeattributereal
145c6c1daeSBarry Smith #endif
155c6c1daeSBarry Smith 
petscviewerhdf5opengroup_(PetscViewer * viewer,char path[],hid_t * fileId,hid_t * groupId,int * ierr,PETSC_FORTRAN_CHARLEN_T len)163014b61aSVaclav Hapla PETSC_EXTERN void petscviewerhdf5opengroup_(PetscViewer *viewer, char path[], hid_t *fileId, hid_t *groupId, int *ierr, PETSC_FORTRAN_CHARLEN_T len)
173014b61aSVaclav Hapla {
183014b61aSVaclav Hapla   char *c1;
193014b61aSVaclav Hapla 
203014b61aSVaclav Hapla   FIXCHAR(path, len, c1);
213014b61aSVaclav Hapla   *ierr = PetscViewerHDF5OpenGroup(*viewer, c1, fileId, groupId);
223014b61aSVaclav Hapla   FREECHAR(path, c1);
233014b61aSVaclav Hapla }
24*a968899dSTapashree Pradhan 
petscviewerhdf5writeattributeint_(PetscViewer * viewer,const char parent[],const char name[],PetscInt * value,int * ierr,PETSC_FORTRAN_CHARLEN_T len1,PETSC_FORTRAN_CHARLEN_T len2)25*a968899dSTapashree Pradhan 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*a968899dSTapashree Pradhan {
27*a968899dSTapashree Pradhan   char *c1;
28*a968899dSTapashree Pradhan   char *c2;
29*a968899dSTapashree Pradhan 
30*a968899dSTapashree Pradhan   FIXCHAR(parent, len1, c1);
31*a968899dSTapashree Pradhan   FIXCHAR(name, len2, c2);
32*a968899dSTapashree Pradhan   *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_INT, value);
33*a968899dSTapashree Pradhan   FREECHAR(parent, c1);
34*a968899dSTapashree Pradhan   FREECHAR(name, c2);
35*a968899dSTapashree Pradhan }
36*a968899dSTapashree Pradhan 
petscviewerhdf5writeattributescalar_(PetscViewer * viewer,const char parent[],const char name[],PetscScalar * value,int * ierr,PETSC_FORTRAN_CHARLEN_T len1,PETSC_FORTRAN_CHARLEN_T len2)37*a968899dSTapashree Pradhan 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*a968899dSTapashree Pradhan {
39*a968899dSTapashree Pradhan   char *c1;
40*a968899dSTapashree Pradhan   char *c2;
41*a968899dSTapashree Pradhan 
42*a968899dSTapashree Pradhan   FIXCHAR(parent, len1, c1);
43*a968899dSTapashree Pradhan   FIXCHAR(name, len2, c2);
44*a968899dSTapashree Pradhan   *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_SCALAR, value);
45*a968899dSTapashree Pradhan   FREECHAR(parent, c1);
46*a968899dSTapashree Pradhan   FREECHAR(name, c2);
47*a968899dSTapashree Pradhan }
48*a968899dSTapashree Pradhan 
petscviewerhdf5writeattributereal_(PetscViewer * viewer,const char parent[],const char name[],PetscReal * value,int * ierr,PETSC_FORTRAN_CHARLEN_T len1,PETSC_FORTRAN_CHARLEN_T len2)49*a968899dSTapashree Pradhan 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*a968899dSTapashree Pradhan {
51*a968899dSTapashree Pradhan   char *c1;
52*a968899dSTapashree Pradhan   char *c2;
53*a968899dSTapashree Pradhan 
54*a968899dSTapashree Pradhan   FIXCHAR(parent, len1, c1);
55*a968899dSTapashree Pradhan   FIXCHAR(name, len2, c2);
56*a968899dSTapashree Pradhan   *ierr = PetscViewerHDF5WriteAttribute(*viewer, c1, c2, PETSC_REAL, value);
57*a968899dSTapashree Pradhan   FREECHAR(parent, c1);
58*a968899dSTapashree Pradhan   FREECHAR(name, c2);
59*a968899dSTapashree Pradhan }
60