Skip to content

Commit

Permalink
Create compound types so that we can avoid using groups to keep
Browse files Browse the repository at this point in the history
attributes together.
  • Loading branch information
jedbrown committed Nov 12, 2009
1 parent bfa55e6 commit 61a3446
Show file tree
Hide file tree
Showing 3 changed files with 148 additions and 13 deletions.
45 changes: 45 additions & 0 deletions src/fs/impls/cont/contviewdhm.c
Original file line number Diff line number Diff line change
Expand Up @@ -126,12 +126,40 @@ dErr dFSView_Cont_DHM(dFS fs,dViewer viewer)
/* Units, \todo not yet stored in the FS so we write defaults here */
err = dViewerDHMWriteDimensions(viewer,fsgrp,"velocity","m s-1",exp(1));dCHK(err);

{
dht_FS fs5;
dht_Field *field5;
hid_t h5t_fs;
dInt i;

err = dViewerDHMGetFSType(viewer,&h5t_fs);dCHK(err);
err = dMeshGetTagName(fs->mesh,fs->degreetag,&fs5.degree);dCHK(err);
err = dMeshGetTagName(fs->mesh,fs->gcoffsetTag,&fs5.global_offset);dCHK(err);
err = PetscStrallocpy("mock_partition",&fs5.partition);dCHK(err);
fs5.fields.len = fs->bs;
err = dMallocA(fs5.fields.len,&field5);dCHK(err);
for (i=0; i<fs->bs; i++) {
field5[i].name = fs->fieldname[i];
field5[i].units.dimensions = (char*)"m s-1"; /* we only use it as \c const */
field5[i].units.scale = exp(1);
}
fs5.fields.p = field5;
att = H5Acreate(fsgrp,"FS_attribute",h5t_fs,dhm->h5s_scalar,H5P_DEFAULT,H5P_DEFAULT);dH5CHK(att,H5Acreate);
herr = H5Awrite(att,h5t_fs,&fs5);dH5CHK(herr,H5Awrite);
herr = H5Aclose(att);dH5CHK(herr,H5Aclose);
err = dFree(field5);dCHK(err);
err = dFree(fs5.partition);dCHK(err);
err = dFree(fs5.global_offset);dCHK(err);
err = dFree(fs5.degree);dCHK(err);
}

herr = H5Gclose(fsgrp); if (herr < 0) dERROR(PETSC_ERR_LIB,"H5Gclose");
dFunctionReturn(0);
}

static dErr VecView_Dohp_FSCont_DHM(Vec X,PetscViewer viewer)
{
dViewer_DHM *dhm = viewer->data;
dFS fs;
hid_t fslink,xlink,fspace,mspace,dset,plist;
hsize_t gdim[2],offset[2],count[2];
Expand Down Expand Up @@ -179,6 +207,23 @@ static dErr VecView_Dohp_FSCont_DHM(Vec X,PetscViewer viewer)
namelen = H5Iget_name(fslink,fsname,sizeof fsname);dH5CHK(namelen,H5Iget_name);
if (!namelen) dERROR(PETSC_ERR_LIB,"Could not get FS path");
herr = H5Lcreate_soft(fsname,xlink,"fs",H5P_DEFAULT,H5P_DEFAULT);dH5CHK(herr,H5Lcreate_soft);
{
dht_Vec vecatt[1];
hsize_t dims[1] = {1};
hid_t space,attr,vectype;

err = dViewerDHMGetVecType(viewer,&vectype);dCHK(err);
/* Initialize data vecatt[0] */
herr = H5Rcreate(&vecatt[0].fs,dhm->file,fsname,H5R_OBJECT,-1);dH5CHK(herr,H5Rcreate);
vecatt[0].time = dhm->time;
err = PetscObjectStateQuery((PetscObject)viewer,&vecatt[0].state);dCHK(err);
/* Write attribute */
space = H5Screate_simple(1,dims,NULL);dH5CHK(space,H5Screate_simple);
attr = H5Acreate(dset,"meta",vectype,space,H5P_DEFAULT,H5P_DEFAULT);dH5CHK(attr,H5Acreate);
herr = H5Awrite(attr,vectype,vecatt);dH5CHK(herr,H5Awrite);
herr = H5Aclose(attr);dH5CHK(herr,H5Aclose);
herr = H5Sclose(space);dH5CHK(herr,H5Aclose);
}
}

herr = H5Dclose(dset);dH5CHK(herr,H5Dclose);
Expand Down
88 changes: 76 additions & 12 deletions src/viewer/dhm.c
Original file line number Diff line number Diff line change
Expand Up @@ -9,7 +9,7 @@ dErr dViewerDHMGetStringTypes(PetscViewer viewer,hid_t *fstring,hid_t *mstring,h

dFunctionBegin;
if (dhm->h5t_fstring < 0) {
dhm->h5t_fstring = H5Tcopy(H5T_FORTRAN_S1);dH5CHK(dhm->h5t_fstring,H5Tcopy);
dhm->h5t_fstring = H5Tcopy(H5T_C_S1);dH5CHK(dhm->h5t_fstring,H5Tcopy);
herr = H5Tset_size(dhm->h5t_fstring,H5T_VARIABLE);dH5CHK(herr,H5Tset_size);
herr = H5Tcommit(dhm->typeroot,"string",dhm->h5t_fstring,H5P_DEFAULT,H5P_DEFAULT,H5P_DEFAULT);dH5CHK(herr,H5Tcommit);
}
Expand Down Expand Up @@ -55,6 +55,66 @@ dErr dViewerDHMWriteDimensions(PetscViewer viewer,hid_t grp,const char *name,con
dFunctionReturn(0);
}

dErr dViewerDHMGetFSType(PetscViewer viewer,hid_t *intype)
{
dViewer_DHM *dhm = viewer->data;
dErr err;

dFunctionBegin;
if (dhm->h5t_fs < 0) {
hid_t strtype,fstype,unittype,fieldtype,fieldstype;
herr_t herr;

err = dViewerDHMGetStringTypes(viewer,&strtype,NULL,NULL);dCHK(err);
unittype = H5Tcreate(H5T_COMPOUND,sizeof(dht_Units));dH5CHK(unittype,H5Tcreate);
herr = H5Tinsert(unittype,"dimensions",offsetof(dht_Units,dimensions),strtype);dH5CHK(herr,H5Tinsert);
herr = H5Tinsert(unittype,"scale",offsetof(dht_Units,scale),dH5T_SCALAR);dH5CHK(herr,H5Tinsert);

herr = H5Tcommit(dhm->typeroot,"Units_type",unittype,H5P_DEFAULT,H5P_DEFAULT,H5P_DEFAULT);dH5CHK(herr,H5Tcommit); /* skip this? */

fieldtype = H5Tcreate(H5T_COMPOUND,sizeof(dht_Field));dH5CHK(fieldtype,H5Tcreate);
herr = H5Tinsert(fieldtype,"name",offsetof(dht_Field,name),strtype);dH5CHK(herr,H5Tinsert);
herr = H5Tinsert(fieldtype,"units",offsetof(dht_Field,units),unittype);dH5CHK(herr,H5Tinsert);

herr = H5Tcommit(dhm->typeroot,"Field_type",fieldtype,H5P_DEFAULT,H5P_DEFAULT,H5P_DEFAULT);dH5CHK(herr,H5Tcommit); /* skip this? */

fieldstype = H5Tvlen_create(fieldtype);dH5CHK(fieldstype,H5Tvlen_create);

fstype = H5Tcreate(H5T_COMPOUND,sizeof(dht_FS));dH5CHK(fstype,H5Tcreate);
herr = H5Tinsert(fstype,"degree",offsetof(dht_FS,degree),strtype);dH5CHK(herr,H5Tinsert);
herr = H5Tinsert(fstype,"global_offset",offsetof(dht_FS,global_offset),strtype);dH5CHK(herr,H5Tinsert);
herr = H5Tinsert(fstype,"partition",offsetof(dht_FS,partition),strtype);dH5CHK(herr,H5Tinsert);
herr = H5Tinsert(fstype,"fields",offsetof(dht_FS,fields),fieldstype);dH5CHK(herr,H5Tinsert);

herr = H5Tcommit(dhm->typeroot,"FS_type",fstype,H5P_DEFAULT,H5P_DEFAULT,H5P_DEFAULT);dH5CHK(herr,H5Tcommit);
herr = H5Tclose(fieldstype);dH5CHK(herr,H5Tclose);
herr = H5Tclose(fieldtype);dH5CHK(herr,H5Tclose);
herr = H5Tclose(unittype);dH5CHK(herr,H5Tclose);
dhm->h5t_fs = fstype;
}
*intype = dhm->h5t_fs;
dFunctionReturn(0);
}

dErr dViewerDHMGetVecType(PetscViewer viewer,hid_t *intype)
{
dViewer_DHM *dhm = viewer->data;

dFunctionBegin;
if (dhm->h5t_vec < 0) {
hid_t vectype;
herr_t herr;
vectype = H5Tcreate(H5T_COMPOUND,sizeof(dht_Vec));dH5CHK(vectype,H5Tcreate);
herr = H5Tinsert(vectype,"FS",offsetof(dht_Vec,fs),H5T_STD_REF_OBJ);dH5CHK(herr,H5Tinsert);
herr = H5Tinsert(vectype,"Time",offsetof(dht_Vec,time),dH5T_REAL);dH5CHK(herr,H5Tinsert);
herr = H5Tinsert(vectype,"State",offsetof(dht_Vec,state),dH5T_INT);dH5CHK(herr,H5Tinsert);
herr = H5Tcommit(dhm->typeroot,"Vec_type",vectype,H5P_DEFAULT,H5P_DEFAULT,H5P_DEFAULT);dH5CHK(herr,H5Tinsert);
dhm->h5t_vec = vectype;
}
*intype = dhm->h5t_vec;
dFunctionReturn(0);
}

dErr dViewerDHMSetUp(PetscViewer viewer)
{
dViewer_DHM *dhm = viewer->data;
Expand All @@ -63,7 +123,7 @@ dErr dViewerDHMSetUp(PetscViewer viewer)
dErr err;

dFunctionBegin;
if (dhm->file > 0) dFunctionReturn(0);
if (dhm->file >= 0) dFunctionReturn(0);
/* Set attributes for opening the file */
plist_id = H5Pcreate(H5P_FILE_ACCESS);
herr = H5Pset_fapl_mpio(plist_id,((PetscObject)viewer)->comm,MPI_INFO_NULL);dH5CHK(herr,H5Pset_fapl_mpio);
Expand Down Expand Up @@ -110,16 +170,18 @@ static dErr PetscViewerDestroy_DHM(PetscViewer v)

dFunctionBegin;
err = dFree(dhm->filename);dCHK(err);
if (dhm->h5t_fstring > 0) {herr = H5Tclose(dhm->h5t_fstring);dH5CHK(herr,H5Tclose);}
if (dhm->h5t_mstring > 0) {herr = H5Tclose(dhm->h5t_mstring);dH5CHK(herr,H5Tclose);}
if (dhm->h5s_scalar > 0) {herr = H5Sclose(dhm->h5s_scalar);dH5CHK(herr,H5Sclose);}
if (dhm->curstep > 0) {herr = H5Gclose(dhm->curstep);dH5CHK(herr,H5Gclose);}
if (dhm->meshroot > 0) {herr = H5Gclose(dhm->meshroot);dH5CHK(herr,H5Gclose);}
if (dhm->fsroot > 0) {herr = H5Gclose(dhm->fsroot);dH5CHK(herr,H5Gclose);}
if (dhm->steproot > 0) {herr = H5Gclose(dhm->steproot);dH5CHK(herr,H5Gclose);}
if (dhm->typeroot > 0) {herr = H5Gclose(dhm->typeroot);dH5CHK(herr,H5Gclose);}
if (dhm->dohproot > 0) {herr = H5Gclose(dhm->dohproot);dH5CHK(herr,H5Gclose);}
if (dhm->file > 0) {herr = H5Fclose(dhm->file);dH5CHK(herr,H5Fclose);}
if (dhm->h5t_vec >= 0) {herr = H5Tclose(dhm->h5t_vec);dH5CHK(herr,H5Tclose);}
if (dhm->h5t_fs >= 0) {herr = H5Tclose(dhm->h5t_fs);dH5CHK(herr,H5Tclose);}
if (dhm->h5t_fstring >= 0) {herr = H5Tclose(dhm->h5t_fstring);dH5CHK(herr,H5Tclose);}
if (dhm->h5t_mstring >= 0) {herr = H5Tclose(dhm->h5t_mstring);dH5CHK(herr,H5Tclose);}
if (dhm->h5s_scalar >= 0) {herr = H5Sclose(dhm->h5s_scalar);dH5CHK(herr,H5Sclose);}
if (dhm->curstep >= 0) {herr = H5Gclose(dhm->curstep);dH5CHK(herr,H5Gclose);}
if (dhm->meshroot >= 0) {herr = H5Gclose(dhm->meshroot);dH5CHK(herr,H5Gclose);}
if (dhm->fsroot >= 0) {herr = H5Gclose(dhm->fsroot);dH5CHK(herr,H5Gclose);}
if (dhm->steproot >= 0) {herr = H5Gclose(dhm->steproot);dH5CHK(herr,H5Gclose);}
if (dhm->typeroot >= 0) {herr = H5Gclose(dhm->typeroot);dH5CHK(herr,H5Gclose);}
if (dhm->dohproot >= 0) {herr = H5Gclose(dhm->dohproot);dH5CHK(herr,H5Gclose);}
if (dhm->file >= 0) {herr = H5Fclose(dhm->file);dH5CHK(herr,H5Fclose);}
err = dFree(dhm);dCHK(err);
dFunctionReturn(0);
}
Expand Down Expand Up @@ -216,6 +278,8 @@ PetscErrorCode PetscViewerCreate_DHM(PetscViewer v)
dhm->curstep = -1;
dhm->h5t_fstring = -1;
dhm->h5t_mstring = -1;
dhm->h5t_vec = -1;
dhm->h5t_fs = -1;
dhm->h5s_scalar = -1;

err = PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscViewerFileSetName_C","PetscViewerFileSetName_DHM",
Expand Down
28 changes: 27 additions & 1 deletion src/viewer/dhm.h
Original file line number Diff line number Diff line change
Expand Up @@ -49,13 +49,39 @@ typedef struct {
char *timeunits;
dReal timescale;
dInt stepnumber;
hid_t h5t_mstring,h5t_fstring,h5s_scalar;
hid_t h5t_mstring,h5t_fstring,h5s_scalar,h5t_fs,h5t_vec;
} dViewer_DHM;

extern dErr dViewerDHMSetUp(dViewer);
extern dErr dViewerDHMGetStringTypes(PetscViewer,hid_t *fstring,hid_t *mstring,hid_t *sspace);
extern dErr dViewerDHMGetStep(PetscViewer viewer,hid_t *step);
extern dErr dViewerDHMAttributeStringWrite(PetscViewer viewer,hid_t grp,const char *attname,const char *str);
extern dErr dViewerDHMWriteDimensions(PetscViewer viewer,hid_t grp,const char *name,const char *units,dReal scale);
extern dErr dViewerDHMGetFSType(PetscViewer viewer,hid_t *type);
extern dErr dViewerDHMGetVecType(PetscViewer viewer,hid_t *type);

/* Compound structures for HDF5 attributes */
typedef struct {
char *dimensions;
dScalar scale;
} dht_Units;

typedef struct {
char *name;
dht_Units units;
} dht_Field;

typedef struct {
char *degree;
char *global_offset;
char *partition;
hvl_t fields;
} dht_FS;

typedef struct {
hobj_ref_t fs;
dReal time;
dInt state;
} dht_Vec;

#endif

0 comments on commit 61a3446

Please sign in to comment.