#include #include #include #include #include #include using Grid = Dune::YaspGrid; using namespace AMDiS; int main(int argc, char** argv) { Environment env(argc, argv); Dune::Timer t; int numLevels = GRIDDIM == 2 ? 6 : 4; 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)); }; auto grad_g = [](auto const& x) -> FieldMatrix { return {-20.0 * std::exp(-10.0 * dot(x,x)) * x}; }; using Param = LagrangeBasis; ProblemStat prob("ellipt"); prob.initialize(INIT_ALL); auto opL = makeOperator(tag::gradtest_gradtrial{}, 1.0); prob.addMatrixOperator(opL); auto opForce = makeOperator(tag::test{}, f, 6); prob.addVectorOperator(opForce); // set boundary condition prob.addDirichletBC(1, g); AdaptInfo adaptInfo("adapt"); std::vector errL2; errL2.reserve(numLevels); std::vector errH1; errH1.reserve(numLevels); std::vector widths; widths.reserve(numLevels); for (int l = 0; l < numLevels; ++l) { prob.globalRefine(1); double h = 0; for (auto const& e : elements(prob.gridView(), Dune::Partitions::interior)) { auto refElem = Dune::referenceElement(e.type()); auto geo = e.geometry(); for (int i = 0; i < refElem.size(GRIDDIM-1); ++i) { // edges auto v0 = geo.global(refElem.position(refElem.subEntity(i,GRIDDIM-1,0,GRIDDIM), GRIDDIM)); auto v1 = geo.global(refElem.position(refElem.subEntity(i,GRIDDIM-1,1,GRIDDIM), GRIDDIM)); h = std::max(h, (v0 - v1).two_norm()); } } widths.push_back(h); prob.assemble(adaptInfo); prob.solve(adaptInfo); double errorL2 = integrate(sqr(g - prob.solution()), prob.gridView(), 6); errL2.push_back(std::sqrt(errorL2)); double errorH1 = errorL2 + integrate(unary_dot(grad_g - gradientAtQP(prob.solution())), prob.gridView(), 6); errH1.push_back(std::sqrt(errorH1)); } // write last solution to file prob.writeFiles(adaptInfo); msg(""); msg("{:5} | {:12} | {:12} | {:12} | {:12} | {:12}", "level", "h", "|u - u*|_L2","|u - u*|_H1","eoc_L2","eoc_H1"); msg_("{0:->7}{0:->15}{0:->15}{0:->15}{0:->15}{1:->14}","+","\n"); std::vector 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]); } for (int i = 0; i < numLevels; ++i) msg("{:<5} | {:<12} | {:<12} | {:<12} | {:<12} | {:<12}", i+1, widths[i], errL2[i], errH1[i], eocL2[i], eocH1[i]); msg("elapsed time: {} seconds", t.elapsed()); return 0; }