Skip to content
Snippets Groups Projects
Commit a5b484da authored by Chichi Lalescu's avatar Chichi Lalescu Committed by Cristian Lalescu
Browse files

WIP: change base class hierarchy

I want to properly accomodate postprocessing codes, so I need something
more general than `direct_numerical_simulation`.
parent df0f18fc
No related branches found
No related tags found
1 merge request!21Bugfix/nansampling
#include "code_base.hpp"
#include "scope_timer.hpp"
code_base::code_base(
const MPI_Comm COMMUNICATOR,
const std::string &simulation_name):
comm(COMMUNICATOR),
simname(simulation_name)
{
MPI_Comm_rank(this->comm, &this->myrank);
MPI_Comm_size(this->comm, &this->nprocs);
this->stop_code_now = false;
}
int code_base::check_stopping_condition(void)
{
if (myrank == 0)
{
std::string fname = (
std::string("stop_") +
std::string(this->simname));
{
struct stat file_buffer;
this->stop_code_now = (
stat(fname.c_str(), &file_buffer) == 0);
}
}
MPI_Bcast(
&this->stop_code_now,
1,
MPI_C_BOOL,
0,
MPI_COMM_WORLD);
return EXIT_SUCCESS;
}
/**********************************************************************
* *
* Copyright 2017 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 <http://www.gnu.org/licenses/> *
* *
* Contact: Cristian.Lalescu@ds.mpg.de *
* *
**********************************************************************/
#ifndef CODE_BASE_HPP
#define CODE_BASE_HPP
#include <cstdlib>
#include <sys/types.h>
#include <sys/stat.h>
#include "base.hpp"
class code_base
{
private:
clock_t time0, time1;
public:
int myrank, nprocs;
MPI_Comm comm;
std::string simname;
bool stop_code_now;
int nx;
int ny;
int nz;
int dealias_type;
double dkx;
double dky;
double dkz;
code_base(
const MPI_Comm COMMUNICATOR,
const std::string &simulation_name);
virtual ~code_base(){}
int check_stopping_condition(void);
int start_simple_timer(void)
{
this->time0 = clock();
return EXIT_SUCCESS;
}
int print_simple_timer(void)
{
this->time1 = clock();
double local_time_difference = ((
(unsigned int)(this->time1 - this->time0)) /
((double)CLOCKS_PER_SEC));
double time_difference = 0.0;
MPI_Allreduce(
&local_time_difference,
&time_difference,
1,
MPI_DOUBLE,
MPI_SUM,
MPI_COMM_WORLD);
if (this->myrank == 0)
std::cout << "iteration " << this->iteration <<
" took " << time_difference/this->nprocs <<
" seconds" << std::endl;
if (this->myrank == 0)
std::cerr << "iteration " << this->iteration <<
" took " << time_difference/this->nprocs <<
" seconds" << std::endl;
this->time0 = this->time1;
return EXIT_SUCCESS;
}
};
#endif//CODE_BASE_HPP
...@@ -3,64 +3,15 @@ ...@@ -3,64 +3,15 @@
#include <sys/stat.h> #include <sys/stat.h>
#include "direct_numerical_simulation.hpp" #include "direct_numerical_simulation.hpp"
#include "scope_timer.hpp" #include "scope_timer.hpp"
#include "hdf5_tools.hpp"
int grow_single_dataset(hid_t dset, int tincrement)
{
int ndims;
hsize_t space;
space = H5Dget_space(dset);
ndims = H5Sget_simple_extent_ndims(space);
hsize_t *dims = new hsize_t[ndims];
H5Sget_simple_extent_dims(space, dims, NULL);
dims[0] += tincrement;
H5Dset_extent(dset, dims);
H5Sclose(space);
delete[] dims;
return EXIT_SUCCESS;
}
herr_t grow_dataset_visitor(
hid_t o_id,
const char *name,
const H5O_info_t *info,
void *op_data)
{
if (info->type == H5O_TYPE_DATASET)
{
hsize_t dset = H5Dopen(o_id, name, H5P_DEFAULT);
grow_single_dataset(dset, *((int*)(op_data)));
H5Dclose(dset);
}
return EXIT_SUCCESS;
}
direct_numerical_simulation::direct_numerical_simulation(
const MPI_Comm COMMUNICATOR,
const std::string &simulation_name):
comm(COMMUNICATOR),
simname(simulation_name)
{
MPI_Comm_rank(this->comm, &this->myrank);
MPI_Comm_size(this->comm, &this->nprocs);
this->stop_code_now = false;
}
int direct_numerical_simulation::grow_file_datasets() int direct_numerical_simulation::grow_file_datasets()
{ {
int file_problems = 0; return hdf5_tools::grow_file_datasets(
this->stat_file,
hid_t group; "statistics",
group = H5Gopen(this->stat_file, "/statistics", H5P_DEFAULT); this->niter_todo / this->niter_stat);
int tincrement = this->niter_todo / this->niter_stat;
H5Ovisit(
group,
H5_INDEX_NAME,
H5_ITER_NATIVE,
grow_dataset_visitor,
&tincrement);
H5Gclose(group);
return file_problems;
} }
int direct_numerical_simulation::read_iteration(void) int direct_numerical_simulation::read_iteration(void)
...@@ -164,25 +115,3 @@ int direct_numerical_simulation::main_loop(void) ...@@ -164,25 +115,3 @@ int direct_numerical_simulation::main_loop(void)
return EXIT_SUCCESS; return EXIT_SUCCESS;
} }
int direct_numerical_simulation::check_stopping_condition(void)
{
if (myrank == 0)
{
std::string fname = (
std::string("stop_") +
std::string(this->simname));
{
struct stat file_buffer;
this->stop_code_now = (
stat(fname.c_str(), &file_buffer) == 0);
}
}
MPI_Bcast(
&this->stop_code_now,
1,
MPI_C_BOOL,
0,
MPI_COMM_WORLD);
return EXIT_SUCCESS;
}
...@@ -31,45 +31,24 @@ ...@@ -31,45 +31,24 @@
#include <sys/types.h> #include <sys/types.h>
#include <sys/stat.h> #include <sys/stat.h>
#include "base.hpp" #include "base.hpp"
#include "full_code/code_base.hpp"
int grow_single_dataset(hid_t dset, int tincrement); class direct_numerical_simulation: public code_base
herr_t grow_dataset_visitor(
hid_t o_id,
const char *name,
const H5O_info_t *info,
void *op_data);
class direct_numerical_simulation
{ {
private:
clock_t time0, time1;
public: public:
int myrank, nprocs;
MPI_Comm comm;
std::string simname;
int iteration, checkpoint; int iteration, checkpoint;
int checkpoints_per_file; int checkpoints_per_file;
int niter_out; int niter_out;
int niter_stat; int niter_stat;
int niter_todo; int niter_todo;
hid_t stat_file; hid_t stat_file;
bool stop_code_now;
int nx;
int ny;
int nz;
int dealias_type;
double dkx;
double dky;
double dkz;
direct_numerical_simulation( direct_numerical_simulation(
const MPI_Comm COMMUNICATOR, const MPI_Comm COMMUNICATOR,
const std::string &simulation_name); const std::string &simulation_name):
code_base(
COMMUNICATOR,
simulation_name){}
virtual ~direct_numerical_simulation(){} virtual ~direct_numerical_simulation(){}
virtual int write_checkpoint(void) = 0; virtual int write_checkpoint(void) = 0;
...@@ -82,39 +61,6 @@ class direct_numerical_simulation ...@@ -82,39 +61,6 @@ class direct_numerical_simulation
int read_iteration(void); int read_iteration(void);
int write_iteration(void); int write_iteration(void);
int grow_file_datasets(void); int grow_file_datasets(void);
int check_stopping_condition(void);
int start_simple_timer(void)
{
this->time0 = clock();
return EXIT_SUCCESS;
}
int print_simple_timer(void)
{
this->time1 = clock();
double local_time_difference = ((
(unsigned int)(this->time1 - this->time0)) /
((double)CLOCKS_PER_SEC));
double time_difference = 0.0;
MPI_Allreduce(
&local_time_difference,
&time_difference,
1,
MPI_DOUBLE,
MPI_SUM,
MPI_COMM_WORLD);
if (this->myrank == 0)
std::cout << "iteration " << this->iteration <<
" took " << time_difference/this->nprocs <<
" seconds" << std::endl;
if (this->myrank == 0)
std::cerr << "iteration " << this->iteration <<
" took " << time_difference/this->nprocs <<
" seconds" << std::endl;
this->time0 = this->time1;
return EXIT_SUCCESS;
}
}; };
#endif//DIRECT_NUMERICAL_SIMULATION_HPP #endif//DIRECT_NUMERICAL_SIMULATION_HPP
......
#include "hdf5_tools.hpp"
int hdf5_tools::grow_single_dataset(hid_t dset, int tincrement)
{
int ndims;
hsize_t space;
space = H5Dget_space(dset);
ndims = H5Sget_simple_extent_ndims(space);
hsize_t *dims = new hsize_t[ndims];
H5Sget_simple_extent_dims(space, dims, NULL);
dims[0] += tincrement;
H5Dset_extent(dset, dims);
H5Sclose(space);
delete[] dims;
return EXIT_SUCCESS;
}
herr_t hdf5_tools::grow_dataset_visitor(
hid_t o_id,
const char *name,
const H5O_info_t *info,
void *op_data)
{
if (info->type == H5O_TYPE_DATASET)
{
hsize_t dset = H5Dopen(o_id, name, H5P_DEFAULT);
grow_single_dataset(dset, *((int*)(op_data)));
H5Dclose(dset);
}
return EXIT_SUCCESS;
}
int hdf5_tools::grow_file_datasets(
const hid_t stat_file,
const std::string group_name,
const int tincrement)
{
int file_problems = 0;
hid_t group;
group = H5Gopen(stat_file, group_name.c_str(), H5P_DEFAULT);
H5Ovisit(
group,
H5_INDEX_NAME,
H5_ITER_NATIVE,
grow_dataset_visitor,
&tincrement);
H5Gclose(group);
return file_problems;
}
/**********************************************************************
* *
* Copyright 2017 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 <http://www.gnu.org/licenses/> *
* *
* Contact: Cristian.Lalescu@ds.mpg.de *
* *
**********************************************************************/
#ifndef HDF5_TOOLS_HPP
#define HDF5_TOOLS_HPP
#include "base.hpp"
namespace hdf5_tools
{
int grow_single_dataset(
hid_t dset,
int tincrement);
herr_t grow_dataset_visitor(
hid_t o_id,
const char *name,
const H5O_info_t *info,
void *op_data);
int grow_file_datasets(
const hid_t stat_file,
const std::string group_name,
const int tincrement);
}
#endif//HDF5_TOOLS_HPP
...@@ -35,6 +35,36 @@ int NSVEparticles<rnumber>::do_stats() ...@@ -35,6 +35,36 @@ int NSVEparticles<rnumber>::do_stats()
return EXIT_SUCCESS; return EXIT_SUCCESS;
} }
template <typename rnumber>
int direct_numerical_simulation::main_loop(void)
{
this->start_simple_timer();
int max_iter = (this->iteration + this->niter_todo -
(this->iteration % this->niter_todo));
for (; this->iteration < max_iter;)
{
#ifdef USE_TIMINGOUTPUT
const std::string loopLabel = ("code::main_start::loop-" +
std::to_string(this->iteration));
TIMEZONE(loopLabel.c_str());
#endif
this->do_stats();
this->step();
if (this->iteration % this->niter_out == 0)
this->write_checkpoint();
this->check_stopping_condition();
if (this->stop_code_now)
break;
this->print_simple_timer();
}
this->do_stats();
this->print_simple_timer();
if (this->iteration % this->niter_out != 0)
this->write_checkpoint();
return EXIT_SUCCESS;
}
template class NSVEparticles<float>; template class NSVEparticles<float>;
template class NSVEparticles<double>; template class NSVEparticles<double>;
...@@ -39,8 +39,6 @@ template <typename rnumber> ...@@ -39,8 +39,6 @@ template <typename rnumber>
class ppNSVE: public NSVE<rnumber> class ppNSVE: public NSVE<rnumber>
{ {
public: public:
ppNSVE( ppNSVE(
const MPI_Comm COMMUNICATOR, const MPI_Comm COMMUNICATOR,
const std::string &simulation_name): const std::string &simulation_name):
......
...@@ -90,6 +90,8 @@ print('This is bfps version ' + VERSION) ...@@ -90,6 +90,8 @@ print('This is bfps version ' + VERSION)
### lists of files and MANIFEST.in ### lists of files and MANIFEST.in
src_file_list = ['full_code/ppNSVE', src_file_list = ['full_code/ppNSVE',
'field_binary_IO', 'field_binary_IO',
'full_code/hdf5_tools',
'full_code/code_base',
'full_code/direct_numerical_simulation', 'full_code/direct_numerical_simulation',
'full_code/NSVE', 'full_code/NSVE',
'full_code/NSVEparticles', 'full_code/NSVEparticles',
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment