#ifdef HAVE_CONFIG_H #include "config.h" #endif #include <iostream> #include <amdis/AMDiS.hpp> #include <amdis/Integrate.hpp> #include <amdis/LocalOperators.hpp> #include <amdis/ProblemStat.hpp> #include <amdis/common/Literals.hpp> using namespace AMDiS; using namespace Dune::Indices; // 1 component with polynomial degree 1 using Param = YaspGridBasis<GRIDDIM, 2>; using ElliptProblem = ProblemStat<Param>; int main(int argc, char** argv) { AMDiS::init(argc, argv); int numLevels = GRIDDIM == 2 ? 8 : 5; 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<double,1,GRIDDIM> { return {-20.0 * std::exp(-10.0 * dot(x,x)) * x}; }; ElliptProblem prob("ellipt"); prob.initialize(INIT_ALL); auto opL = makeOperator(tag::gradtest_gradtrial{}, 1.0); prob.addMatrixOperator(opL, 0, 0); auto opForce = makeOperator(tag::test{}, f, 6); prob.addVectorOperator(opForce, 0); // set boundary condition prob.addDirichletBC(BoundaryType{1}, 0, 0, g); AdaptInfo adaptInfo("adapt"); 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) { prob.grid()->globalRefine(1); auto gridView = prob.gridView(); double h = 0; for (auto const& e : edges(gridView)) h = std::max(h, e.geometry().volume()); widths.push_back(h); prob.globalBasis()->update(gridView); prob.solutionVector()->compress(); prob.assemble(adaptInfo); prob.solve(adaptInfo); double errorL2 = integrate(sqr(g - prob.solution(0)), gridView, 6); errL2.push_back(std::sqrt(errorL2)); double errorH1 = errorL2 + integrate(unary_dot(grad_g - gradientAtQP(prob.solution(0))), gridView, 6); errH1.push_back(std::sqrt(errorH1)); #if WRITE_FILES Dune::VTKWriter<typename ElliptProblem::GridView> vtkWriter(gridView); vtkWriter.addVertexData(prob.solution(0), Dune::VTK::FieldInfo("u", Dune::VTK::FieldInfo::Type::scalar, 1)); vtkWriter.write("u_" + std::to_string(i)); #endif } std::cout << "\n"; 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<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]); } 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]); AMDiS::finalize(); return 0; }