Skip to content
Snippets Groups Projects
Commit ee7732b3 authored by Praetorius, Simon's avatar Praetorius, Simon
Browse files

add functions backup() and restore() to ProblemStat

parent 73513f33
Branches
Tags
1 merge request!38Backup and restore facilities
...@@ -47,6 +47,12 @@ namespace AMDiS ...@@ -47,6 +47,12 @@ namespace AMDiS
value = pt().get(key, value); value = pt().get(key, value);
} }
template <class T>
static void set(std::string const& key, T const& value)
{
pt()[key] = value;
}
/// Returns specified info level /// Returns specified info level
static int getMsgInfo() static int getMsgInfo()
{ {
......
...@@ -274,6 +274,12 @@ namespace AMDiS ...@@ -274,6 +274,12 @@ namespace AMDiS
/// Writes output files. /// Writes output files.
void writeFiles(AdaptInfo& adaptInfo, bool force = false); void writeFiles(AdaptInfo& adaptInfo, bool force = false);
/// Implements \ref ProblemStatBase::backup
void backup(AdaptInfo& adaptInfo) override;
/// Implements \ref ProblemStatBase::restore
void restore(Flag initFlag) override;
public: // get-methods public: // get-methods
......
...@@ -484,4 +484,56 @@ writeFiles(AdaptInfo& adaptInfo, bool force) ...@@ -484,4 +484,56 @@ writeFiles(AdaptInfo& adaptInfo, bool force)
} }
template <class Traits>
void ProblemStat<Traits>::
backup(AdaptInfo& adaptInfo)
{
auto param = Parameters::get<std::string>(name_ + "->backup filename");
std::string filename = param ? *param : name_ + "_" + std::to_string(adaptInfo.timestepNumber()) + ".backup";
Dune::BackupRestoreFacility<Grid>::backup(*grid_, filename);
// TODO: backup data
msg("Problem backuped to file '{}'.", filename);
}
template <class Traits>
void ProblemStat<Traits>::
restore(Flag initFlag)
{
std::string filename = Parameters::get<std::string>(name_ + "->restore filename").value();
test_exit(filesystem::exists(filename), "Restore file '{}' not found.", filename);
// restore grid from file
std::shared_ptr<Grid> grid(Dune::BackupRestoreFacility<Grid>::restore(filename));
adoptGrid(grid);
// create fespace
if (initFlag.isSet(INIT_FE_SPACE) || initFlag.isSet(INIT_SYSTEM))
createGlobalBasis();
// create system
if (initFlag.isSet(INIT_SYSTEM))
createMatricesAndVectors();
// create solver
if (linearSolver_)
warning("solver already created\n");
else if (initFlag.isSet(INIT_SOLVER))
createSolver();
// create marker
if (initFlag.isSet(INIT_MARKER))
createMarker();
// create file writer
if (initFlag.isSet(INIT_FILEWRITER))
createFileWriter();
solution_->compress();
// TODO: read data
}
} // end namespace AMDiS } // end namespace AMDiS
...@@ -111,6 +111,12 @@ namespace AMDiS ...@@ -111,6 +111,12 @@ namespace AMDiS
*/ */
virtual void estimate(AdaptInfo& adaptInfo) = 0; virtual void estimate(AdaptInfo& adaptInfo) = 0;
/// \brief backup the grid and the solution to a file
virtual void backup(AdaptInfo& adaptInfo) = 0;
/// \brief restore the grid and the solution from a file
virtual void restore(Flag initFlag) = 0;
/// Returns the name of the problem. /// Returns the name of the problem.
virtual std::string const& name() const = 0; virtual std::string const& name() const = 0;
}; };
......
#include <amdis/AMDiS.hpp>
#include <amdis/ProblemStat.hpp>
#include "Tests.hpp"
using namespace AMDiS;
template <class Grid, class Factory>
void test(Factory factory)
{
using Param = TaylorHoodBasis<Grid>;
using Problem = ProblemStat<Param>;
std::size_t num_elements = 0;
std::size_t num_vertices = 0;
// backup
{
std::unique_ptr<Grid> grid(factory());
Problem prob("test", *grid);
prob.initialize(INIT_ALL);
prob.globalRefine(2);
num_elements = prob.grid()->size(0);
num_vertices = prob.grid()->size(Grid::dimension);
AdaptInfo adaptInfo("adapt");
prob.backup(adaptInfo);
}
// restore
{
Problem prob("test");
prob.restore(INIT_ALL);
AMDIS_TEST_EQ(num_elements, prob.grid()->size(0));
AMDIS_TEST_EQ(num_vertices, prob.grid()->size(Grid::dimension));
}
}
template <class Grid>
void test_cube()
{
using Factory = Dune::StructuredGridFactory<Grid>;
test<Grid>([]() { return Factory::createCubeGrid({0.0,0.0}, {1.0,1.0}, std::array<unsigned int,2>{2,2}); });
}
template <class Grid>
void test_simplex()
{
using Factory = Dune::StructuredGridFactory<Grid>;
test<Grid>([]() { return Factory::createSimplexGrid({0.0,0.0}, {1.0,1.0}, std::array<unsigned int,2>{2,2}); });
}
int main(int argc, char** argv)
{
AMDiS::init(argc, argv);
std::string filename = "test.backup";
Parameters::set("test->backup filename", filename);
Parameters::set("test->restore filename", filename);
test_cube<Dune::YaspGrid<2>>();
// test_cube<Dune::UGGrid<2>>();
// test_simplex<Dune::AlbertaGrid<2,2>>();
// test_simplex<Dune::UGGrid<2>>();
AMDiS::finalize();
return 0;
}
...@@ -2,6 +2,10 @@ ...@@ -2,6 +2,10 @@
dune_add_test(SOURCES AdaptInfoTest.cpp dune_add_test(SOURCES AdaptInfoTest.cpp
LINK_LIBRARIES amdis) LINK_LIBRARIES amdis)
dune_add_test(SOURCES BackupRestoreTest.cpp
LINK_LIBRARIES amdis)
add_dune_alberta_flags(GRIDDIM 2 WORLDDIM 2 BackupRestoreTest)
dune_add_test(SOURCES ConceptsTest.cpp dune_add_test(SOURCES ConceptsTest.cpp
LINK_LIBRARIES amdis) LINK_LIBRARIES amdis)
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment