AMDiS  2.10
The Adaptive Multi-Dimensional Simulation Toolbox
ellipt.cc
1 #include <amdis/AMDiS.hpp>
3 #include <amdis/Integrate.hpp>
4 #include <amdis/LocalOperators.hpp>
5 #include <amdis/ProblemStat.hpp>
7 
9 #include <dune/grid/yaspgrid.hh>
10 using Grid = Dune::YaspGrid<GRIDDIM>;
12 
14 using namespace AMDiS;
16 
18 int main(int argc, char** argv)
19 {
22  Environment env(argc, argv);
23  Dune::Timer t;
24 
25  int numLevels = GRIDDIM == 2 ? 6 : 4;
26  if (argc > 2)
27  numLevels = std::atoi(argv[2]);
28 
29  auto f = [](auto const& x) {
30  double r2 = dot(x,x);
31  double ux = std::exp(-10.0 * r2);
32  return -(400.0 * r2 - 20.0 * x.size()) * ux;
33  };
34  auto g = [](auto const& x){ return std::exp(-10.0 * dot(x,x)); };
35  auto grad_g = [](auto const& x) -> FieldMatrix<double,1,GRIDDIM> {
36  return {-20.0 * std::exp(-10.0 * dot(x,x)) * x};
37  };
39 
41  using Param = LagrangeBasis<Grid, 2>;
42  ProblemStat<Param> prob("ellipt");
43  prob.initialize(INIT_ALL);
45 
47  // ⟨∇u,∇v⟩
48  auto opL = makeOperator(tag::gradtest_gradtrial{}, 1.0);
49  prob.addMatrixOperator(opL);
50 
51  // (f,v)
52  auto opForce = makeOperator(tag::test{}, f, 6);
53  prob.addVectorOperator(opForce);
55 
57  // set boundary condition
58  prob.addDirichletBC(1, g);
60 
62  AdaptInfo adaptInfo("adapt");
63 
64  std::vector<double> errL2; errL2.reserve(numLevels);
65  std::vector<double> errH1; errH1.reserve(numLevels);
66  std::vector<double> widths; widths.reserve(numLevels);
68 
70  for (int l = 0; l < numLevels; ++l) {
71  prob.globalRefine(1);
73 
75  double h = 0;
76  for (auto const& e : elements(prob.gridView(), Dune::Partitions::interior)) {
77  auto refElem = Dune::referenceElement<double,GRIDDIM>(e.type());
78  auto geo = e.geometry();
79  for (int i = 0; i < refElem.size(GRIDDIM-1); ++i) { // edges
80  auto v0 = geo.global(refElem.position(refElem.subEntity(i,GRIDDIM-1,0,GRIDDIM), GRIDDIM));
81  auto v1 = geo.global(refElem.position(refElem.subEntity(i,GRIDDIM-1,1,GRIDDIM), GRIDDIM));
82  h = std::max(h, (v0 - v1).two_norm());
83  }
84  }
85  widths.push_back(h);
87 
89  prob.assemble(adaptInfo);
90  prob.solve(adaptInfo);
92 
94  double errorL2 = integrate(sqr(g - prob.solution()), prob.gridView(), 6);
95  errL2.push_back(std::sqrt(errorL2));
96  double errorH1 = errorL2 + integrate(unary_dot(grad_g - gradientOf(prob.solution())), prob.gridView(), 6);
97  errH1.push_back(std::sqrt(errorH1));
98  }
100 
102  // write last solution to file
103  prob.writeFiles(adaptInfo);
104 
105  msg("");
106  msg("{:5} | {:12} | {:12} | {:12} | {:12} | {:12}",
107  "level", "h", "|u - u*|_L2","|u - u*|_H1","eoc_L2","eoc_H1");
108  msg_("{0:->7}{0:->15}{0:->15}{0:->15}{0:->15}{1:->14}","+","\n");
109 
110  std::vector<double> eocL2(numLevels, 0.0), eocH1(numLevels, 0.0);
111  for (int i = 1; i < numLevels; ++i) {
112  eocL2[i] = std::log(errL2[i]/errL2[i-1]) / std::log(widths[i]/widths[i-1]);
113  eocH1[i] = std::log(errH1[i]/errH1[i-1]) / std::log(widths[i]/widths[i-1]);
114  }
115 
116  for (int i = 0; i < numLevels; ++i)
117  msg("{:<5} | {:<12} | {:<12} | {:<12} | {:<12} | {:<12}",
118  i+1, widths[i], errL2[i], errH1[i], eocL2[i], eocH1[i]);
119 
120  msg("elapsed time: {} seconds", t.elapsed());
122 
124  return 0;
125 }
Definition: SecondOrderGradTestGradTrial.hpp:20
Definition: AdaptBase.hpp:6
Establishes an environment for sequential and parallel AMDiS programs.
Definition: Environment.hpp:19
ProblemStatTraits for a (composite) basis composed of lagrange bases of different degree...
Definition: ProblemStatTraits.hpp:109
Definition: ProblemStat.hpp:52
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition: GridFunctionOperator.hpp:235
auto gradientOf(Expr const &expr)
Definition: DerivativeGridFunction.hpp:185
Definition: ZeroOrderTest.hpp:17
Holds adapt parameters and infos about the problem.
Definition: AdaptInfo.hpp:25