#ifdef HAVE_CONFIG_H #include "config.h" #endif #include #include #include #include #include #include #include #include #include #include using namespace AMDiS; struct NavierStokesBasis { using Grid = Dune::YaspGrid; using GlobalBasis = typename TaylorHoodBasis::GlobalBasis; }; using StokesProblem = ProblemStat; using StokesProblemInstat = ProblemInstat; int main(int argc, char** argv) { AMDiS::init(argc, argv); StokesProblem prob("stokes"); prob.initialize(INIT_ALL); StokesProblemInstat probInstat("stokes", prob); probInstat.initialize(INIT_UH_OLD); double viscosity = 1.0; double density = 1.0; double vel = 1.0; Parameters::get("stokes->viscosity", viscosity); Parameters::get("stokes->density", density); Parameters::get("stokes->boundary velocity", vel); AdaptInfo adaptInfo("adapt"); // tree-paths for components auto _v = Dune::Indices::_0; auto _p = Dune::Indices::_1; auto invTau = std::ref(probInstat.invTau()); // <1/tau * u, v> auto opTime = makeOperator(tag::testvec_trialvec{}, density * invTau); prob.addMatrixOperator(opTime, _v, _v); // <1/tau * u^old, v> auto opTimeOld = makeOperator(tag::testvec{}, density * invTau * prob.solution(_v)); prob.addVectorOperator(opTimeOld, _v); // sum_i + + auto opStokes = makeOperator(tag::stokes{}, viscosity); prob.addMatrixOperator(opStokes, treepath(), treepath()); // <(u * nabla)u_i^old, v_i> auto opNonlin1 = makeOperator(tag::testvec_trialvec{}, density * trans(gradientAtQP(prob.solution(_v)))); prob.addMatrixOperator(opNonlin1, _v, _v); for (std::size_t i = 0; i < AMDIS_DOW; ++i) { // <(u^old * nabla)u_i, v_i> auto opNonlin2 = makeOperator(tag::test_gradtrial{}, density * prob.solution(_v)); prob.addMatrixOperator(opNonlin2, treepath(_v,i), treepath(_v,i)); } // <(u^old * grad(u_i^old)), v_i> auto opNonlin3 = makeOperator(tag::testvec{}, trans(gradientAtQP(prob.solution(_v))) * prob.solution(_v)); prob.addVectorOperator(opNonlin3, _v); // define boundary regions auto left = [](auto const& x) { return x[0] < 1.e-8; }; auto not_left = [](auto const& x) { return x[0] > 1.0 - 1.e-8 || x[1] < 1.e-8 || x[1] > 1.0 - 1.e-8; }; // define boundary values auto parabolic_y = [](auto const& x) -> Dune::FieldVector { return {0.0, x[1]*(1.0 - x[1])}; }; auto zero = [](auto const& x) -> Dune::FieldVector { return {0.0, 0.0}; }; // set boundary conditions for velocity prob.addDirichletBC(left, _v, _v, parabolic_y); prob.addDirichletBC(not_left, _v, _v, zero); // set point constraint for pressure prob.addDirichletBC([](auto const& x) { return x[0] < 1.e-8 && x[1] < 1.e-8; }, _p, _p, 0.0); // set initial conditions prob.solution(_v).interpolate(parabolic_y); prob.solution(_p).interpolate(0.0); AdaptInstationary adapt("adapt", prob, adaptInfo, probInstat, adaptInfo); adapt.adapt(); // output solution prob.writeFiles(adaptInfo); AMDiS::finalize(); return 0; }