2#include <amdis/AMDiS.hpp>
3#include <amdis/Integrate.hpp>
4#include <amdis/LocalOperators.hpp>
5#include <amdis/ProblemStat.hpp>
9#include <dune/grid/yaspgrid.hh>
10using Grid = Dune::YaspGrid<GRIDDIM>;
18int main(
int argc,
char** argv)
25 int numLevels = GRIDDIM == 2 ? 6 : 4;
27 numLevels = std::atoi(argv[2]);
29 auto f = [](
auto const& x) {
31 double ux = std::exp(-10.0 * r2);
32 return -(400.0 * r2 - 20.0 * x.size()) * ux;
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};
43 prob.initialize(INIT_ALL);
49 prob.addMatrixOperator(opL);
53 prob.addVectorOperator(opForce);
58 prob.addDirichletBC(1, g);
64 std::vector<double> errL2; errL2.reserve(numLevels);
65 std::vector<double> errH1; errH1.reserve(numLevels);
66 std::vector<double> widths; widths.reserve(numLevels);
70 for (
int l = 0; l < numLevels; ++l) {
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) {
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());
89 prob.assemble(adaptInfo);
90 prob.solve(adaptInfo);
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));
103 prob.writeFiles(adaptInfo);
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");
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]);
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]);
120 msg(
"elapsed time: {} seconds", t.elapsed());
Holds adapt parameters and infos about the problem.
Definition: AdaptInfo.hpp:26
Establishes an environment for sequential and parallel AMDiS programs.
Definition: Environment.hpp:20
Definition: ProblemStat.hpp:55
auto gradientOf(Expr const &expr)
Definition: DerivativeGridFunction.hpp:185
auto two_norm(Vec &&vec)
Applies Operation::TwoNorm to a vector-valued GridFunction.
Definition: Operations.hpp:217
auto dot(Lhs &&lhs, Rhs &&rhs)
Applies Operation::Dot to two vector-valued GridFunctions.
Definition: Operations.hpp:255
auto unary_dot(Vec &&vec)
Applies Operation::UnaryDot to a vector-valued GridFunction.
Definition: Operations.hpp:201
auto sqr(T &&value)
Applies Operation::Sqr to GridFunction.
Definition: Operations.hpp:134
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition: GridFunctionOperator.hpp:235
ProblemStatTraits for a (composite) basis composed of lagrange bases of different degree.
Definition: ProblemStatTraits.hpp:111
Definition: SecondOrderGradTestGradTrial.hpp:20
Definition: ZeroOrderTest.hpp:17