Skip to content
Snippets Groups Projects
Commit 3a8242e7 authored by graeser's avatar graeser
Browse files

Update to changes in dune-functions

parent e4253f1f
Branches
No related tags found
No related merge requests found
......@@ -39,7 +39,7 @@
#include <dune/localfunctions/lagrange/pqkfactory.hh>
// included dune-functions headers
#include <dune/functions/functionspacebases/pqknodalbasis.hh>
#include <dune/functions/functionspacebases/lagrangebasis.hh>
#include <dune/functions/functionspacebases/interpolate.hh>
#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
......@@ -74,12 +74,10 @@ void assemblePoissonProblem(
ElementRhs elementRhs;
auto localView = basis.localView();
auto localIndexSet = basis.localIndexSet();
for(const auto& e: Dune::elements(gridView))
{
localView.bind(e);
localIndexSet.bind(localView);
const auto& localFiniteElement = localView.tree().finiteElement();
......@@ -90,10 +88,10 @@ void assemblePoissonProblem(
for(std::size_t i=0; i<localSize; ++i)
{
std::size_t globalI = localIndexSet.index(i);
std::size_t globalI = localView.index(i);
for(std::size_t j=0; j<localSize; ++j)
{
std::size_t globalJ = localIndexSet.index(j);
std::size_t globalJ = localView.index(j);
matrixBuilder[globalI][globalJ] += elementMatrix[i][j];
}
......@@ -110,12 +108,10 @@ void computeBoundaryVertices(const GridView& gridView, BitVector& isBoundary, co
using namespace Dune::FuTutorial;
auto localView = basis.localView();
auto localIndexSet = basis.localIndexSet();
for(const auto& element: elements(gridView))
{
localView.bind(element);
localIndexSet.bind(localView);
const auto& localFiniteElement = localView.tree().finiteElement();
std::size_t localSize = localFiniteElement.localBasis().size();
......@@ -130,7 +126,7 @@ void computeBoundaryVertices(const GridView& gridView, BitVector& isBoundary, co
if (localKey.codim() > 0)
for(const auto& subEntity: subEntities(referenceElement(element), insideFacet(intersection), codim(localKey.codim())))
if (subEntity == localKey.subEntity())
isBoundary[localIndexSet.index(i)] = true;
isBoundary[localView.index(i)] = true;
}
}
}
......@@ -191,7 +187,7 @@ int main(int argc, char** argv)
};
using Basis = typename Dune::Functions::PQkNodalBasis<GridView,4>;
using Basis = typename Dune::Functions::LagrangeBasis<GridView,4>;
Basis basis(gridView);
assemblePoissonProblem(gridView, A, rhs, rhsFunction, basis);
......@@ -201,7 +197,6 @@ int main(int argc, char** argv)
auto dirichletFunction = [](auto x) {
return sin(x[0] * 2.0*M_PI)*.1;
};
// Dune::Functions::interpolate(basis, Dune::TypeTree::hybridTreePath(), dirichletNodes, dirichletFunction);
Dune::Functions::interpolate(basis, x, dirichletFunction);
std::vector<bool> isBoundary;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment