Skip to content
Snippets Groups Projects
Commit 05eec5cd authored by Praetorius, Simon's avatar Praetorius, Simon
Browse files

Merge branch 'feature/integrate_documentation' into 'develop'

Feature/integrate documentation

See merge request spraetor/dune-amdis!47
parents c27b2e73 5e247bf0
No related branches found
No related tags found
No related merge requests found
......@@ -8,12 +8,13 @@ namespace AMDiS
{
namespace Impl
{
template <class GridFct, class GridView, class QuadProvider>
auto integrateImpl(GridFct&& gridFct, GridView const& gridView, QuadProvider makeQuad)
template <class GF, class GridView, class QuadProvider>
auto integrateImpl(GF&& gf, GridView const& gridView, QuadProvider makeQuad)
{
auto localFct = localFunction(gridFct);
auto localFct = localFunction(std::forward<GF>(gf));
using Range = typename std::decay_t<GridFct>::Range;
using GridFct = std::decay_t<GF>;
using Range = typename GridFct::Range;
Range result(0);
for (auto const& element : elements(gridView)) {
......@@ -33,7 +34,15 @@ namespace AMDiS
} // end namespace Impl
/// Integrate expression with quadrature rule determined by polynomial order of expression
/// \brief Integrate expression with quadrature rule determined by polynomial order of expression
/**
* **Example:**
* ```
* double i1 = integrate(prob.solution(0), prob.gridView());
* double i2 = integrate(X(0) + X(1), prob.gridView());
* double i3 = integrate([](auto const& x) { return x[0] + x[1]; }, prob.gridView()); // ERROR
* ```
**/
template <class Expr, class GridView>
auto integrate(Expr&& expr, GridView const& gridView)
{
......@@ -47,7 +56,7 @@ namespace AMDiS
#if AMDIS_HAS_CXX_CONSTEXPR_IF
if constexpr(expr_has_order)
return Impl::integrateImpl(gridFct, gridView, makeQuad);
return Impl::integrateImpl(std::forward<decltype(gridFct)>(gridFct), gridView, makeQuad);
else
return 0.0;
#else
......@@ -58,6 +67,13 @@ namespace AMDiS
}
/// Integrate expression with quadrature rule provided
/**
* **Example:**
* ```
* auto quad = Dune::QuadratureRules<double,DIM>::rule(Dune::GeometryTypes::triangle, 1);
* double i4 = integrate([](auto const& x) { return x[0]; }, prob.gridView(), quad); // OK
* ```
**/
template <class Expr, class GridView,
class QuadratureRule = Dune::QuadratureRule<typename GridView::ctype, GridView::dimension>>
auto integrate(Expr&& expr, GridView const& gridView, QuadratureRule const& quad)
......@@ -68,6 +84,12 @@ namespace AMDiS
}
/// Integrate expression with quadrature rule determined by provided polynomial `degree`
/**
* **Example:**
* ```
* double i5 = integrate([](auto const& x) { return x[0]; }, prob.gridView(), 1); // OK
* ```
**/
template <class Expr, class GridView>
auto integrate(Expr&& expr, GridView const& gridView, int degree,
Dune::QuadratureType::Enum qt = Dune::QuadratureType::GaussLegendre)
......
......@@ -24,6 +24,10 @@ dune_add_test(SOURCES FiniteElementTypeTest.cpp
dune_add_test(SOURCES FilesystemTest.cpp
LINK_LIBRARIES amdis)
dune_add_test(SOURCES IntegrateTest.cpp
LINK_LIBRARIES amdis
CMD_ARGS "${CMAKE_SOURCE_DIR}/examples/init/ellipt.dat.2d")
dune_add_test(SOURCES MarkerTest.cpp
LINK_LIBRARIES amdis
CMD_ARGS "${CMAKE_SOURCE_DIR}/examples/init/marker.dat.2d")
......
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#include <iostream>
#include <dune/geometry/type.hh>
#include <amdis/AMDiS.hpp>
#include <amdis/ProblemStat.hpp>
#include <amdis/Operators.hpp>
#include <amdis/common/Literals.hpp>
#include <amdis/common/FieldMatVec.hpp>
#include <amdis/gridfunctions/Integrate.hpp>
#include "Tests.hpp"
using namespace AMDiS;
using ElliptParam = YaspGridBasis<2, 2>;
using ElliptProblem = ProblemStat<ElliptParam>;
int main(int argc, char** argv)
{
AMDiS::init(argc, argv);
ElliptProblem prob("ellipt");
prob.initialize(INIT_ALL);
auto u = prob.getSolution(0);
auto gv = u.basis().gridView();
u << 1.0;
auto i1 = integrate(1.0, gv);
auto i2 = integrate(u, gv);
auto i3 = integrate([](auto const& x) { return 1.0; }, gv, 0);
AMDIS_TEST(i1 == 1.0);
AMDIS_TEST(i2 == 1.0);
AMDIS_TEST(i3 == 1.0);
auto quad = Dune::QuadratureRules<double,2>::rule(Dune::GeometryTypes::quadrilateral, 0);
auto i4 = integrate([](auto const& x) { return 1.0; }, gv, quad);
AMDIS_TEST(i4 == 1.0);
AMDiS::finalize();
return 0;
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment