ellipt.cc 2.95 KB
Newer Older
1
2
3
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif
Praetorius, Simon's avatar
Praetorius, Simon committed
4
5
#include <iostream>

6
#include <amdis/AMDiS.hpp>
7
#include <amdis/Integrate.hpp>
8
#include <amdis/LocalOperators.hpp>
9
10
#include <amdis/ProblemStat.hpp>
#include <amdis/common/Literals.hpp>
Praetorius, Simon's avatar
Praetorius, Simon committed
11

12
using namespace AMDiS;
13
using namespace Dune::Indices;
14

Praetorius, Simon's avatar
Praetorius, Simon committed
15
// 1 component with polynomial degree 1
16
using Param   = YaspGridBasis<GRIDDIM, 2>;
17
using ElliptProblem = ProblemStat<Param>;
Praetorius, Simon's avatar
Praetorius, Simon committed
18
19
20

int main(int argc, char** argv)
{
21
  Environment env(argc, argv);
22

23
  int numLevels = GRIDDIM == 2 ? 8 : 5;
24
25
26
27
28
29
30
31
32
  if (argc > 2)
    numLevels = std::atoi(argv[2]);

  auto f = [](auto const& x) {
    double r2 = dot(x,x);
    double ux = std::exp(-10.0 * r2);
    return -(400.0 * r2 - 20.0 * x.size()) * ux;
  };
  auto g = [](auto const& x){ return std::exp(-10.0 * dot(x,x)); };
33
  auto grad_g = [](auto const& x) -> FieldMatrix<double,1,GRIDDIM> {
34
35
    return {-20.0 * std::exp(-10.0 * dot(x,x)) * x};
  };
36

37
  ElliptProblem prob("ellipt");
Praetorius, Simon's avatar
Praetorius, Simon committed
38
  prob.initialize(INIT_ALL);
39

40
  auto opL = makeOperator(tag::gradtest_gradtrial{}, 1.0);
Praetorius, Simon's avatar
Praetorius, Simon committed
41
  prob.addMatrixOperator(opL, 0, 0);
42

43
  auto opForce = makeOperator(tag::test{}, f, 6);
Praetorius, Simon's avatar
Praetorius, Simon committed
44
  prob.addVectorOperator(opForce, 0);
45

Praetorius, Simon's avatar
Praetorius, Simon committed
46
  // set boundary condition
47
  prob.addDirichletBC(BoundaryType{1}, 0, 0, g);
48
49

  AdaptInfo adaptInfo("adapt");
50

51
52
53
54
  std::vector<double> errL2; errL2.reserve(numLevels);
  std::vector<double> errH1; errH1.reserve(numLevels);
  std::vector<double> widths; widths.reserve(numLevels);
  for (int i = 0; i < numLevels; ++i) {
55
    prob.grid()->globalRefine(1);
56
    auto gridView = prob.gridView();
57
58
59
60
61
62

    double h = 0;
    for (auto const& e : edges(gridView))
      h = std::max(h, e.geometry().volume());
    widths.push_back(h);

63
64
    prob.globalBasis()->update(gridView);
    prob.solutionVector()->compress();
Praetorius, Simon's avatar
Praetorius, Simon committed
65
    prob.assemble(adaptInfo);
66
67
    prob.solve(adaptInfo);

Praetorius, Simon's avatar
Praetorius, Simon committed
68
    double errorL2 = integrate(sqr(g - prob.solution(0)), gridView, 6);
69
    errL2.push_back(std::sqrt(errorL2));
Praetorius, Simon's avatar
Praetorius, Simon committed
70
    double errorH1 = errorL2 + integrate(unary_dot(grad_g - gradientAtQP(prob.solution(0))), gridView, 6);
71
72
73
    errH1.push_back(std::sqrt(errorH1));

#if WRITE_FILES
74
    Dune::VTKWriter<typename ElliptProblem::GridView> vtkWriter(gridView);
Praetorius, Simon's avatar
Praetorius, Simon committed
75
    vtkWriter.addVertexData(prob.solution(0), Dune::VTK::FieldInfo("u", Dune::VTK::FieldInfo::Type::scalar, 1));
76
77
78
    vtkWriter.write("u_" + std::to_string(i));
#endif
  }
79

80
  std::cout << "\n";
81
  msg("{:5} | {:12} | {:12} | {:12} | {:12} | {:12}",
82
      "level", "h", "|u - u*|_L2","|u - u*|_H1","eoc_L2","eoc_H1");
83
  msg_("{0:->7}{0:->15}{0:->15}{0:->15}{0:->15}{1:->14}","+","\n");
84

85
86
87
88
  std::vector<double> eocL2(numLevels, 0.0), eocH1(numLevels, 0.0);
  for (int i = 1; i < numLevels; ++i) {
    eocL2[i] = std::log(errL2[i]/errL2[i-1]) / std::log(widths[i]/widths[i-1]);
    eocH1[i] = std::log(errH1[i]/errH1[i-1]) / std::log(widths[i]/widths[i-1]);
Praetorius, Simon's avatar
Praetorius, Simon committed
89
  }
90

91
  for (int i = 0; i < numLevels; ++i)
92
    msg("{:<5} | {:<12} | {:<12} | {:<12} | {:<12} | {:<12}",
93
        i+1, widths[i], errL2[i], errH1[i], eocL2[i], eocH1[i]);
94

Praetorius, Simon's avatar
Praetorius, Simon committed
95
  return 0;
96
}