processfields: reorganized hdf field & mesh dumps

pull/1/head
Thorsten Liebig 2010-12-27 13:54:09 +01:00
parent 2254e01fd6
commit fbab37e0ad
5 changed files with 103 additions and 57 deletions

View File

@ -46,53 +46,6 @@ ProcessFields::~ProcessFields()
} }
} }
void ProcessFields::InitProcess()
{
if (Enabled==false) return;
//get the correct direction names for all coordinate systems
string names[] = {Op->GetDirName(0),Op->GetDirName(1),Op->GetDirName(2)};
if (m_fileType==HDF5_FILETYPE)
{
m_filename+= ".h5";
H5::H5File* file = new H5::H5File( m_filename , H5F_ACC_TRUNC );
H5::Group* group = new H5::Group( file->createGroup( "/Mesh" ));
for (int n=0; n<3; ++n)
{
hsize_t dimsf[1]; // dataset dimensions
dimsf[0] = numLines[n];
H5::DataSpace dataspace( 1, dimsf );
H5::FloatType datatype( H5::PredType::NATIVE_FLOAT );
H5::DataSet dataset = group->createDataSet( names[n].c_str(), datatype, dataspace );
//convert to float...
float* array = new float[numLines[n]];
for (unsigned int i=0; i<numLines[n]; ++i)
{
#ifdef OUTPUT_IN_DRAWINGUNITS
array[i] = Lines[n][i];
#else
if ((m_Mesh_Type==CYLINDRICAL_MESH) && (n==1)) //check for alpha-direction
array[i] = discLines[n][i];
else
array[i] = discLines[n][i] * Op->GetGridDelta();
#endif
}
//write to dataset
dataset.write( array, H5::PredType::NATIVE_FLOAT );
}
delete group;
group = new H5::Group( file->createGroup( "/FieldData" ));
delete group;
group = new H5::Group( file->createGroup( "/FieldData/FD" ));
delete group;
group = new H5::Group( file->createGroup( "/FieldData/TD" ));
delete group;
delete file;
}
}
string ProcessFields::GetFieldNameByType(DumpType type) string ProcessFields::GetFieldNameByType(DumpType type)
{ {
switch (type) switch (type)
@ -418,14 +371,57 @@ bool ProcessFields::DumpMultiScalarArray2VTK(ofstream &file, string names[], FDT
return true; return true;
} }
bool ProcessFields::DumpVectorArray2HDF5(string filename, string name, FDTD_FLOAT const* const* const* const* array, unsigned int const* numLines, float time)
bool ProcessFields::WriteMesh2HDF5(string filename, string groupName, unsigned int const* numLines, double const* const* discLines, MeshType meshT, double discLines_scaling)
{
H5::H5File file( filename, H5F_ACC_RDWR );
H5::Group hdf_group( file.openGroup( groupName ));
string names[] = {"x","y","z"};
if (meshT==CYLINDRICAL_MESH)
{
names[0]="rho";
names[1]="alpha";
}
H5::Group* group = new H5::Group( hdf_group.createGroup( "/Mesh" ));
for (int n=0; n<3; ++n)
{
hsize_t dimsf[1]; // dataset dimensions
dimsf[0] = numLines[n];
H5::DataSpace dataspace( 1, dimsf );
H5::FloatType datatype( H5::PredType::NATIVE_FLOAT );
H5::DataSet dataset = group->createDataSet( names[n].c_str(), datatype, dataspace );
//convert to float...
float* array = new float[numLines[n]];
for (unsigned int i=0; i<numLines[n]; ++i)
{
#ifdef OUTPUT_IN_DRAWINGUNITS
array[i] = Lines[n][i];
#else
if ((meshT==CYLINDRICAL_MESH) && (n==1)) //check for alpha-direction
array[i] = discLines[n][i];
else
array[i] = discLines[n][i] * discLines_scaling;
#endif
}
//write to dataset
dataset.write( array, H5::PredType::NATIVE_FLOAT );
}
delete group;
return true;
}
bool ProcessFields::DumpVectorArray2HDF5(string filename, string groupName, string name, FDTD_FLOAT const* const* const* const* array, unsigned int const* numLines, float time)
{ {
const H5std_string FILE_NAME(filename); const H5std_string FILE_NAME(filename);
const H5std_string DATASET_NAME( name ); const H5std_string DATASET_NAME( name );
H5::H5File file( FILE_NAME, H5F_ACC_RDWR ); H5::H5File file( FILE_NAME, H5F_ACC_RDWR );
H5::Group group( file.openGroup( "/FieldData/TD" )); H5::Group group( file.openGroup( groupName ));
hsize_t dimsf[4]; // dataset dimensions hsize_t dimsf[4]; // dataset dimensions
@ -466,7 +462,7 @@ bool ProcessFields::DumpVectorArray2HDF5(string filename, string name, FDTD_FLOA
return true; return true;
} }
bool ProcessFields::DumpVectorArray2HDF5(string filename, string name, std::complex<float> const* const* const* const* array, unsigned int const* numLines, float weight, float frequency) bool ProcessFields::DumpVectorArray2HDF5(string filename, string groupName, string name, std::complex<float> const* const* const* const* array, unsigned int const* numLines, float weight, float frequency)
{ {
const H5std_string FILE_NAME(filename); const H5std_string FILE_NAME(filename);
const H5std_string DATASET_NAME_RE( name + "_real"); const H5std_string DATASET_NAME_RE( name + "_real");
@ -474,7 +470,7 @@ bool ProcessFields::DumpVectorArray2HDF5(string filename, string name, std::comp
H5::H5File file( FILE_NAME, H5F_ACC_RDWR ); H5::H5File file( FILE_NAME, H5F_ACC_RDWR );
H5::Group group( file.openGroup( "/FieldData/FD" )); H5::Group group( file.openGroup( groupName ));
hsize_t t_dimsf[] = {1}; hsize_t t_dimsf[] = {1};
H5::DataSpace t_dataspace( 1, t_dimsf ); H5::DataSpace t_dataspace( 1, t_dimsf );

View File

@ -32,8 +32,6 @@ public:
enum FileType { VTK_FILETYPE, HDF5_FILETYPE}; enum FileType { VTK_FILETYPE, HDF5_FILETYPE};
enum DumpType { E_FIELD_DUMP, H_FIELD_DUMP}; enum DumpType { E_FIELD_DUMP, H_FIELD_DUMP};
virtual void InitProcess();
virtual void DefineStartStopCoord(double* dstart, double* dstop); virtual void DefineStartStopCoord(double* dstart, double* dstop);
//! Define a field dump sub sampling rate for a given direction (default: \a dir = -1 means all directions) //! Define a field dump sub sampling rate for a given direction (default: \a dir = -1 means all directions)
@ -71,11 +69,14 @@ public:
static bool DumpScalarArray2VTK(ofstream &file, string name, FDTD_FLOAT const* const* const* array, double const* const* discLines, unsigned int const* numLines, unsigned int precision=12, string header_info = string(), MeshType meshT = CARTESIAN_MESH, double discLines_scaling = 1); static bool DumpScalarArray2VTK(ofstream &file, string name, FDTD_FLOAT const* const* const* array, double const* const* discLines, unsigned int const* numLines, unsigned int precision=12, string header_info = string(), MeshType meshT = CARTESIAN_MESH, double discLines_scaling = 1);
static bool DumpMultiScalarArray2VTK(ofstream &file, string names[], FDTD_FLOAT const* const* const* const* array, unsigned int numFields, double const* const* discLines, unsigned int const* numLines, unsigned int precision=12, string header_info = string(), MeshType meshT = CARTESIAN_MESH, double discLines_scaling = 1); static bool DumpMultiScalarArray2VTK(ofstream &file, string names[], FDTD_FLOAT const* const* const* const* array, unsigned int numFields, double const* const* discLines, unsigned int const* numLines, unsigned int precision=12, string header_info = string(), MeshType meshT = CARTESIAN_MESH, double discLines_scaling = 1);
//! Write a mesh information to the given hdf5-group
static bool WriteMesh2HDF5(string filename, string groupName, unsigned int const* numLines, double const* const* discLines, MeshType meshT = CARTESIAN_MESH, double discLines_scaling = 1);
//! Dump a time-domain vector dump to an HDF5 file //! Dump a time-domain vector dump to an HDF5 file
static bool DumpVectorArray2HDF5(string filename, string name, FDTD_FLOAT const* const* const* const* array, unsigned int const* numLines, float time=0); static bool DumpVectorArray2HDF5(string filename, string groupName, string name, FDTD_FLOAT const* const* const* const* array, unsigned int const* numLines, float time=0);
//! Dump a frequency-domain complex-vector dump to an HDF5 file //! Dump a frequency-domain complex-vector dump to an HDF5 file
static bool DumpVectorArray2HDF5(string filename, string name, std::complex<float> const* const* const* const* array, unsigned int const* numLines, float weight, float frequency); static bool DumpVectorArray2HDF5(string filename, string groupName, string name, std::complex<float> const* const* const* const* array, unsigned int const* numLines, float weight, float frequency);
double CalcTotalEnergy() const; double CalcTotalEnergy() const;

View File

@ -17,6 +17,7 @@
#include "processfields_fd.h" #include "processfields_fd.h"
#include "Common/operator_base.h" #include "Common/operator_base.h"
#include <H5Cpp.h>
#include <iomanip> #include <iomanip>
#include <sstream> #include <sstream>
#include <string> #include <string>
@ -48,6 +49,26 @@ void ProcessFieldsFD::InitProcess()
//setup the hdf5 file //setup the hdf5 file
ProcessFields::InitProcess(); ProcessFields::InitProcess();
if (m_fileType==HDF5_FILETYPE)
{
//create hdf5 file & necessary groups
m_filename+= ".h5";
H5::H5File* file = new H5::H5File( m_filename , H5F_ACC_TRUNC );
H5::Group* group = new H5::Group( file->createGroup( "/FieldData" ));
delete group;
group = new H5::Group( file->createGroup( "/FieldData/FD" ));
delete group;
delete file;
//write mesh information in main root-group
#ifdef OUTPUT_IN_DRAWINGUNITS
double discScaling = 1;
#else
double discScaling = Op->GetGridDelta();
#endif
ProcessFields::WriteMesh2HDF5(m_filename,"/",numLines,discLines,m_Mesh_Type, discScaling);
}
//create data structures... //create data structures...
for (size_t n = 0; n<m_FD_Samples.size(); ++n) for (size_t n = 0; n<m_FD_Samples.size(); ++n)
{ {
@ -196,7 +217,7 @@ void ProcessFieldsFD::DumpFDData()
{ {
stringstream ss; stringstream ss;
ss << "f" << n; ss << "f" << n;
DumpVectorArray2HDF5(m_filename.c_str(), ss.str(), m_FD_Fields.at(n),numLines,1.0/(float)m_FD_SampleCount,m_FD_Samples.at(n)); DumpVectorArray2HDF5(m_filename.c_str(), "/FieldData/FD", ss.str(), m_FD_Fields.at(n),numLines,1.0/(float)m_FD_SampleCount,m_FD_Samples.at(n));
} }
return; return;
} }

View File

@ -17,6 +17,7 @@
#include "processfields_td.h" #include "processfields_td.h"
#include "Common/operator_base.h" #include "Common/operator_base.h"
#include <H5Cpp.h>
#include <iomanip> #include <iomanip>
#include <sstream> #include <sstream>
#include <string> #include <string>
@ -30,6 +31,31 @@ ProcessFieldsTD::~ProcessFieldsTD()
{ {
} }
void ProcessFieldsTD::InitProcess()
{
ProcessFields::InitProcess();
if (m_fileType==HDF5_FILETYPE)
{
//create hdf5 file & necessary groups
m_filename+= ".h5";
H5::H5File* file = new H5::H5File( m_filename, H5F_ACC_TRUNC );
H5::Group* group = new H5::Group( file->createGroup( "/FieldData" ));
delete group;
group = new H5::Group( file->createGroup( "/FieldData/TD" ));
delete group;
delete file;
//write mesh information in main root-group
#ifdef OUTPUT_IN_DRAWINGUNITS
double discScaling = 1;
#else
double discScaling = Op->GetGridDelta();
#endif
ProcessFields::WriteMesh2HDF5(m_filename,"/",numLines,discLines,m_Mesh_Type, discScaling);
}
}
int ProcessFieldsTD::Process() int ProcessFieldsTD::Process()
{ {
if (Enabled==false) return -1; if (Enabled==false) return -1;
@ -69,7 +95,7 @@ int ProcessFieldsTD::Process()
{ {
stringstream ss; stringstream ss;
ss << std::setw( pad_length ) << std::setfill( '0' ) << m_Eng_Interface->GetNumberOfTimesteps(); ss << std::setw( pad_length ) << std::setfill( '0' ) << m_Eng_Interface->GetNumberOfTimesteps();
DumpVectorArray2HDF5(filename.c_str(),string( ss.str() ),field,numLines,m_Eng_Interface->GetTime(m_dualTime)); DumpVectorArray2HDF5(filename.c_str(), "/FieldData/TD", string( ss.str() ), field, numLines, m_Eng_Interface->GetTime(m_dualTime));
} }
else else
cerr << "ProcessFieldsTD::Process: unknown File-Type" << endl; cerr << "ProcessFieldsTD::Process: unknown File-Type" << endl;

View File

@ -26,6 +26,8 @@ public:
ProcessFieldsTD(Engine_Interface_Base* eng_if); ProcessFieldsTD(Engine_Interface_Base* eng_if);
virtual ~ProcessFieldsTD(); virtual ~ProcessFieldsTD();
virtual void InitProcess();
virtual int Process(); virtual int Process();
//! Set the length of the filename timestep pad filled with zeros (default is 8) //! Set the length of the filename timestep pad filled with zeros (default is 8)