polaroid-pp

Schlieren and contour plot tool
git clone https://git.0xfab.ch/polaroid-pp.git
Log | Files | Refs | Submodules | README | LICENSE

commit 930300bb66eaf35a6a8a94b9202afe5a3769a86c
parent c06a4a2a1eb2ba5d6c14258d906bdcb41d274c26
Author: Fabian Wermelinger <fabianw@mavt.ethz.ch>
Date:   Wed, 27 Apr 2016 10:33:39 +0200

added HDF5 photo paper

Diffstat:
Ainclude/PhotoHDF5.h | 49+++++++++++++++++++++++++++++++++++++++++++++++++
Asrc/PhotoHDF5.cpp | 99+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 files changed, 148 insertions(+), 0 deletions(-)

diff --git a/include/PhotoHDF5.h b/include/PhotoHDF5.h @@ -0,0 +1,49 @@ +// File : PhotoHDF5.h +// Date : Tue Apr 26 23:24:00 2016 +// Author : Fabian Wermelinger +// Description: HDF5 Photo Paper +// Copyright 2016 ETH Zurich. All Rights Reserved. +#ifndef PHOTOHDF5_H_VN6WHSEO +#define PHOTOHDF5_H_VN6WHSEO + +#include <hdf5.h> +#include <string> + +#include "common.h" +#include "PhotoPaper.h" + +class PhotoHDF5 : public PhotoPaper +{ +private: + std::string m_fname; + std::string m_description; + bool m_open; + Real* m_hdfraw; + Real m_time; + + // HDF 5 types + hid_t m_hdf5_fileid; + hsize_t m_hdf5_count[4]; + hsize_t m_hdf5_dims[4]; + hsize_t m_hdf5_offset[4]; + + // helper + void _dump_xmf() const; + +public: + PhotoHDF5(const std::string filename="hdf5", const Real t=0) : PhotoPaper(0,0), m_fname(filename), m_description("data"), m_open(false), m_hdfraw(nullptr), m_time(t) {} + PhotoHDF5(const int width, const int height, const std::string filename="hdf5", const Real t=0) : PhotoPaper(0,0), m_fname(filename), m_description("data"), m_open(false), m_hdfraw(nullptr), m_time(t) + { + make_new(width,height); + } + ~PhotoHDF5() { if (m_hdfraw) delete [] m_hdfraw; } + + virtual void make_new(const int width, const int height); + virtual void write(); + virtual void set_pixel(const int x, const int y, const double phi); + virtual std::string suffix() const { return std::string(".h5"); } + virtual void set_description(const char* const desc) { m_description = std::string(desc); } + inline void set_time(const Real t) { m_time = t; } +}; + +#endif /* PHOTOHDF5_H_VN6WHSEO */ diff --git a/src/PhotoHDF5.cpp b/src/PhotoHDF5.cpp @@ -0,0 +1,99 @@ +// File : PhotoHDF5.cpp +// Date : Wed Apr 27 07:52:09 2016 +// Author : Fabian Wermelinger +// Description: HDF5 Photo Implementation +// Copyright 2016 ETH Zurich. All Rights Reserved. +#include <cassert> +#include <cstdio> +#include "PhotoHDF5.h" + +void PhotoHDF5::_dump_xmf() const +{ + FILE *xmf = 0; + xmf = fopen((m_fname+".xmf").c_str(), "w"); + fprintf(xmf, "<?xml version=\"1.0\" ?>\n"); + fprintf(xmf, "<!DOCTYPE Xdmf SYSTEM \"Xdmf.dtd\" []>\n"); + fprintf(xmf, "<Xdmf Version=\"2.0\">\n"); + fprintf(xmf, " <Domain>\n"); + fprintf(xmf, " <Grid GridType=\"Uniform\">\n"); + fprintf(xmf, " <Time Value=\"%e\"/>\n", m_time); + fprintf(xmf, " <Topology TopologyType=\"3DCORECTMesh\" Dimensions=\"%d %d %d\"/>\n", (int)m_hdf5_dims[0], (int)m_hdf5_dims[1], (int)m_hdf5_dims[2]); + fprintf(xmf, " <Geometry GeometryType=\"ORIGIN_DXDYDZ\">\n"); + fprintf(xmf, " <DataItem Name=\"Origin\" Dimensions=\"3\" NumberType=\"Float\" Precision=\"4\" Format=\"XML\">\n"); + fprintf(xmf, " %e %e %e\n", 0.,0.,0.); + fprintf(xmf, " </DataItem>\n"); + fprintf(xmf, " <DataItem Name=\"Spacing\" Dimensions=\"3\" NumberType=\"Float\" Precision=\"4\" Format=\"XML\">\n"); + fprintf(xmf, " %e %e %e\n", 1./(Real)m_hdf5_dims[0],1./(Real)m_hdf5_dims[0],1./(Real)m_hdf5_dims[0]); + fprintf(xmf, " </DataItem>\n"); + fprintf(xmf, " </Geometry>\n"); + fprintf(xmf, " <Attribute Name=\"%s\" AttributeType=\"%s\" Center=\"Node\">\n", m_description.c_str(), "Scalar"); + fprintf(xmf, " <DataItem Dimensions=\"%d %d %d %d\" NumberType=\"Float\" Precision=\"4\" Format=\"HDF\">\n", (int)m_hdf5_dims[0], (int)m_hdf5_dims[1], (int)m_hdf5_dims[2], (int)m_hdf5_dims[3]); + fprintf(xmf, " %s:/%s\n",(m_fname+".h5").c_str(), m_description.c_str()); + fprintf(xmf, " </DataItem>\n"); + fprintf(xmf, " </Attribute>\n"); + fprintf(xmf, " </Grid>\n"); + fprintf(xmf, " </Domain>\n"); + fprintf(xmf, "</Xdmf>\n"); + fclose(xmf); +} + +void PhotoHDF5::make_new(const int width, const int height) +{ + if (m_hdfraw) delete [] m_hdfraw; + m_hdfraw = new Real[width*height]; + + m_width = width; + m_height= height; + + hsize_t tmp[4] = { 1, static_cast<hsize_t>(height), static_cast<hsize_t>(width), 1 }; + for (int i = 0; i < 4; ++i) + { + m_hdf5_count[i] = tmp[i]; + m_hdf5_dims[i] = tmp[i]; + m_hdf5_offset[i]= 0; + } + + H5open(); + hid_t fapl_id = H5Pcreate(H5P_FILE_ACCESS); + m_hdf5_fileid = H5Fcreate((m_fname+".h5").c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, fapl_id); + H5Pclose(fapl_id); + m_open = true; +} + +void PhotoHDF5::write() +{ + if (m_open) + { + hid_t fapl_id = H5Pcreate(H5P_DATASET_XFER); + hid_t fspace_id = H5Screate_simple(4, m_hdf5_dims, NULL); + hid_t dataset_id = H5Dcreate2(m_hdf5_fileid, m_description.c_str(), HDF_PRECISION, fspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); + fspace_id = H5Dget_space(dataset_id); + H5Sselect_hyperslab(fspace_id, H5S_SELECT_SET, m_hdf5_offset, NULL, m_hdf5_count, NULL); + hid_t mspace_id = H5Screate_simple(4, m_hdf5_count, NULL); + H5Dwrite(dataset_id, HDF_PRECISION, mspace_id, fspace_id, fapl_id, m_hdfraw); + + H5Sclose(mspace_id); + H5Sclose(fspace_id); + H5Dclose(dataset_id); + H5Pclose(fapl_id); + + H5Fclose(m_hdf5_fileid); + H5close(); + + _dump_xmf(); + + delete [] m_hdfraw; + m_hdfraw = nullptr; + m_open = false; + } +} + +void PhotoHDF5::set_pixel(const int ix, const int iy, const double phi) +{ + if (m_open) + { + assert(ix >= 0); assert(ix < m_width); + assert(iy >= 0); assert(iy < m_height); + m_hdfraw[iy*m_width + ix] = static_cast<Real>(phi); + } +}