BackupRestoreTest.cpp 3.39 KB
Newer Older
1
#include <amdis/AMDiS.hpp>
2
#include <amdis/Integrate.hpp>
3
4
#include <amdis/ProblemStat.hpp>

5
6
7
8
9
10
11
12
13
14
#if HAVE_DUNE_SPGRID
#include <dune/grid/spgrid.hh>
#endif
#if HAVE_DUNE_ALUGRID
#include <dune/alugrid/grid.hh>
#endif
#if HAVE_DUNE_FOAMGRID
#include <dune/foamgrid/foamgrid.hh>
#endif

15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
#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);

35
36
37
    prob.solution(Dune::Indices::_0) << [](auto const& x) { return x; };
    prob.solution(Dune::Indices::_1) << [](auto const& x) { return two_norm(x); };

38
39
40
41
42
    num_elements = prob.grid()->size(0);
    num_vertices = prob.grid()->size(Grid::dimension);

    AdaptInfo adaptInfo("adapt");
    prob.backup(adaptInfo);
43
    prob.solutionVector()->backup("solution.backup");
44
45
46
47
48
49
50
51
52
  }

  // 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));
53
54
55
56
57
58
59
60
61
62
63
64
65
66

    auto vec = *prob.solutionVector();
    vec.restore("solution.backup");

    prob.solution(Dune::Indices::_0) << [](auto const& x) { return x; };
    prob.solution(Dune::Indices::_1) << [](auto const& x) { return two_norm(x); };

    double error0 = std::sqrt(integrate(
      unary_dot(prob.solution(Dune::Indices::_0) - makeDiscreteFunction(vec, Dune::Indices::_0)), prob.gridView()));
    AMDIS_TEST(error0 < 1.e-10);

    double error1 = std::sqrt(integrate(
      pow<2>(prob.solution(Dune::Indices::_1) - makeDiscreteFunction(vec, Dune::Indices::_1)), prob.gridView()));
    AMDIS_TEST(error1 < 1.e-10);
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
  }
}

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);

92
#if GRID_ID == 0
93
  test_cube<Dune::YaspGrid<2>>();
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
#elif GRID_ID == 1 && HAVE_DUNE_UGGRID
  test_cube<Dune::UGGrid<2>>();
#elif GRID_ID == 2 && HAVE_DUNE_ALUGRID
  test_cube<Dune::ALUGrid<2,2,Dune::cube,Dune::nonconforming>>();
#elif GRID_ID == 3 && HAVE_DUNE_SPGRID
  test<Dune::SPGrid<double,2>>([]() {
    return std::unique_ptr<Dune::SPGrid<double,2>>(
      new Dune::SPGrid<double,2>(Dune::SPDomain<double, 2>({0.0,0.0}, {1.0,1.0}), Dune::SPMultiIndex<2>({2,2}))); });
#elif GRID_ID == 4 && HAVE_DUNE_UGGRID
  test_simplex<Dune::UGGrid<2>>();
#elif GRID_ID == 5 && HAVE_DUNE_ALUGRID
  test_simplex<Dune::ALUGrid<2,2,Dune::simplex,Dune::nonconforming>>();
#elif GRID_ID == 6 && HAVE_DUNE_ALUGRID
  test_simplex<Dune::ALUGrid<2,2,Dune::simplex,Dune::conforming>>();
#endif
  // test_simplex<Dune::AlbertaGrid<2,2>>(); // Segmentation fault
  // test_simplex<Dune::FoamGrid<2,2>>(); // Segmentation fault
111
112

  AMDiS::finalize();
113
  return report_errors();
114
}