/* * Copyright (C) 2011 Thorsten Liebig (Thorsten.Liebig@gmx.de) * * This program is free software: you can redistribute it and/or modify * it under the terms of the GNU General Public License as published by * the Free Software Foundation, either version 3 of the License, or * (at your option) any later version. * * This program is distributed in the hope that it will be useful, * but WITHOUT ANY WARRANTY; without even the implied warranty of * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * GNU General Public License for more details. * * You should have received a copy of the GNU General Public License * along with this program. If not, see . */ using namespace std; #include "hdf5_file_writer.h" #include #include #include #include #include HDF5_File_Writer::HDF5_File_Writer(string filename) { m_filename = filename; m_Group = "/"; hid_t hdf5_file = H5Fcreate(m_filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT); if (hdf5_file<0) { cerr << "HDF5_File_Writer::HDF5_File_Writer: Error, creating the given file """ << m_filename << """ failed" << endl; } H5Fclose(hdf5_file); } HDF5_File_Writer::~HDF5_File_Writer() { } hid_t HDF5_File_Writer::OpenGroup(hid_t hdf5_file, string group) { if (hdf5_file<0) { cerr << "HDF5_File_Writer::CreateGroup: Error, invalid file id" << endl; return -1; } vector results; boost::split(results, group, boost::is_any_of("/")); hid_t grp=H5Gopen(hdf5_file,"/", H5P_DEFAULT); if (grp<0) { cerr << "HDF5_File_Writer::OpenGroup: Error, opening root group " << endl; return -1; } for (size_t n=0;n0)) //check for theta/phi-direction array[i] = discLines[n][i]; else array[i] = discLines[n][i] * scaling; } if (H5Dwrite(dataset, H5T_NATIVE_FLOAT, space, H5P_DEFAULT, H5P_DEFAULT, array)) { cerr << "HDF5_File_Writer::WriteRectMesh: Error, writing to dataset failed" << endl; delete[] array; H5Dclose(dataset); H5Sclose(space); H5Gclose(mesh_grp); H5Fclose(hdf5_file); return false; } delete[] array; H5Dclose(dataset); H5Sclose(space); } H5Gclose(mesh_grp); H5Fclose(hdf5_file); return true; } bool HDF5_File_Writer::WriteScalarField(std::string dataSetName, float const* const* const* field, size_t datasize[3]) { size_t pos = 0; size_t size = datasize[0]*datasize[1]*datasize[2]; size_t n_size[3]={datasize[2],datasize[1],datasize[0]}; float* buffer = new float[size]; for (size_t k=0;k const* const* const* field, size_t datasize[3]) { size_t pos = 0; size_t size = datasize[0]*datasize[1]*datasize[2]; size_t n_size[3]={datasize[2],datasize[1],datasize[0]}; float* buffer = new float[size]; for (size_t k=0;k const* const* const* field, size_t datasize[3]) { size_t pos = 0; size_t size = datasize[0]*datasize[1]*datasize[2]; size_t n_size[3]={datasize[2],datasize[1],datasize[0]}; double* buffer = new double[size]; for (size_t k=0;k const* const* const* const* field, size_t datasize[3]) { size_t pos = 0; size_t size = datasize[0]*datasize[1]*datasize[2]*3; size_t n_size[4]={3,datasize[2],datasize[1],datasize[0]}; float* buffer = new float[size]; for (int n=0;n<3;++n) for (size_t k=0;k const* const* const* const* field, size_t datasize[3]) { size_t pos = 0; size_t size = datasize[0]*datasize[1]*datasize[2]*3; size_t n_size[4]={3,datasize[2],datasize[1],datasize[0]}; double* buffer = new double[size]; for (int n=0;n<3;++n) for (size_t k=0;k values) { float val[values.size()]; for (size_t n=0;n values) { double val[values.size()]; for (size_t n=0;n