/******************************************************************************
* *
* Copyright 2019 Max Planck Institute for Dynamics and Self-Organization *
* *
* This file is part of bfps. *
* *
* bfps 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. *
* *
* bfps 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 bfps. If not, see *
* *
* Contact: Cristian.Lalescu@ds.mpg.de *
* *
******************************************************************************/
#include
#include
#include "NSVE_field_stats.hpp"
#include "fftw_tools.hpp"
#include "scope_timer.hpp"
template
int NSVE_field_stats::initialize(void)
{
TIMEZONE("NSVE_field_stats::initialize");
this->postprocess::read_parameters();
this->vorticity = new field(
nx, ny, nz,
this->comm,
fftw_planner_string_to_flag[this->fftw_plan_rigor]);
this->vorticity->real_space_representation = false;
hid_t parameter_file = H5Fopen(
(this->simname + std::string(".h5")).c_str(),
H5F_ACC_RDONLY,
H5P_DEFAULT);
if (!H5Lexists(parameter_file, "field_dtype", H5P_DEFAULT))
this->bin_IO = NULL;
else
{
hid_t dset = H5Dopen(parameter_file, "field_dtype", H5P_DEFAULT);
hid_t space = H5Dget_space(dset);
hid_t memtype = H5Dget_type(dset);
char *string_data = (char*)malloc(256);
H5Dread(dset, memtype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &string_data);
// check that we're using the correct data type
// field_dtype SHOULD be something like "f4", ">f8"
// first character is ordering, which is machine specific
// for the other two I am checking that they have the correct values
assert(string_data[1] == 'f');
assert(string_data[2] == '0' + sizeof(rnumber));
free(string_data);
H5Sclose(space);
H5Tclose(memtype);
H5Dclose(dset);
this->bin_IO = new field_binary_IO(
this->vorticity->clayout->sizes,
this->vorticity->clayout->subsizes,
this->vorticity->clayout->starts,
this->vorticity->clayout->comm);
}
this->fftw_plan_rigor = hdf5_tools::read_string(parameter_file, "parameters/fftw_plan_rigor");
H5Fclose(parameter_file);
return EXIT_SUCCESS;
}
template
int NSVE_field_stats::read_current_cvorticity(void)
{
TIMEZONE("NSVE_field_stats::read_current_cvorticity");
this->vorticity->real_space_representation = false;
if (this->bin_IO != NULL)
{
char itername[16];
sprintf(itername, "i%.5x", this->iteration);
std::string native_binary_fname = (
this->simname +
std::string("_cvorticity_") +
std::string(itername));
this->bin_IO->read(
native_binary_fname,
this->vorticity->get_cdata());
}
else
{
this->vorticity->io(
this->simname + std::string("_fields.h5"),
"vorticity",
this->iteration,
true);
}
return EXIT_SUCCESS;
}
template
int NSVE_field_stats::read_arbitrary_cvorticity(int arbitrary_iteration)
{
TIMEZONE("NSVE_field_stats::read_arbitrary_cvorticity");
this->vorticity->real_space_representation = false;
if (this->bin_IO != NULL)
{
char itername[16];
sprintf(itername, "i%.5x", arbitrary_iteration);
std::string native_binary_fname = (
this->simname +
std::string("_cvorticity_") +
std::string(itername));
this->bin_IO->read(
native_binary_fname,
this->vorticity->get_cdata());
}
else
{
this->vorticity->io(
this->simname + std::string("_fields.h5"),
"vorticity",
arbitrary_iteration,
true);
}
return EXIT_SUCCESS;
}
template
int NSVE_field_stats::finalize(void)
{
TIMEZONE("NSVE_field_stats::finalize");
if (this->bin_IO != NULL)
delete this->bin_IO;
delete this->vorticity;
return EXIT_SUCCESS;
}
template
int NSVE_field_stats::work_on_current_iteration(void)
{
TIMEZONE("NSVE_field_stats::work_on_current_iteration");
return EXIT_SUCCESS;
}
template class NSVE_field_stats;
template class NSVE_field_stats;