Skip to content
Snippets Groups Projects
Commit 7c7b7339 authored by oliver.sander_at_tu-dresden.de's avatar oliver.sander_at_tu-dresden.de
Browse files

Merge branch 'feature/boundary-dof-power-basis' into 'master'

Boundary dofs: extent to generic dune-functions basis

See merge request !64
parents 658617b3 ffe2c269
No related branches found
No related tags found
1 merge request!64Boundary dofs: extent to generic dune-functions basis
Pipeline #27831 failed
# Master (will become release 2.8) # Master (will become release 2.8)
- ... - constructBoundaryDofs:
- Small interface change in the template parameters: drop `blocksize` and replace it by `BitSetVector`
- The method can now handle generic `dune-functions` basis types, as long as we have consistency in the data types
## Deprecations and removals ## Deprecations and removals
......
...@@ -10,59 +10,56 @@ ...@@ -10,59 +10,56 @@
#include <dune/fufem/boundarypatch.hh> #include <dune/fufem/boundarypatch.hh>
#include <dune/functions/backends/concepts.hh>
#include <dune/functions/backends/istlvectorbackend.hh>
#include <dune/functions/functionspacebases/subentitydofs.hh>
/** \brief For a given basis and boundary patch, determine all degrees /** \brief For a given basis and boundary patch, determine all degrees
of freedom on the patch. of freedom on the patch.
\param boundaryDofs Bit is set if corresponding dofs belongs to the boundary patch A bit of boundaryDofs is set if corresponding dof belongs to the boundary patch.
\param boundaryDofs BitSetVector corresponding to the basis DOFs
\param sfinae Dummy parameter to instantiate this method only for dune-functions bases \param sfinae Dummy parameter to instantiate this method only for dune-functions bases
*/ */
template <class GridView, class Basis, int blocksize> template <class GridView, class Basis, class BitSetVector>
void constructBoundaryDofs(const BoundaryPatch<GridView>& boundaryPatch, void constructBoundaryDofs(const BoundaryPatch<GridView>& boundaryPatch,
const Basis& basis, const Basis& basis,
Dune::BitSetVector<blocksize>& boundaryDofs, BitSetVector& boundaryDofs,
typename Basis::LocalView* sfinae = nullptr) typename Basis::LocalView* sfinae = nullptr)
{ {
// Check consistency of the input // Check consistency of the input
static_assert((std::is_same<GridView, typename Basis::GridView>::value), static_assert((std::is_same<GridView, typename Basis::GridView>::value),
"BoundaryPatch and global basis must be for the same grid view!"); "BoundaryPatch and global basis must be for the same grid view!");
// Small helper function to wrap vectors using istlVectorBackend
// if they do not already satisfy the VectorBackend interface.
auto toVectorBackend = [&](auto& v) -> decltype(auto) {
if constexpr (Dune::models<Dune::Functions::Concept::VectorBackend<Basis>, decltype(v)>()) {
return v;
} else {
return Dune::Functions::istlVectorBackend(v);
}
};
////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////
// Init output bitfield // Init output bitfield
////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////
boundaryDofs.resize(basis.size({})); // ensure IstlVectorBackend interface
boundaryDofs.unsetAll(); auto boundaryDofsBackend = toVectorBackend(boundaryDofs);
auto localView = basis.localView(); boundaryDofsBackend.resize(basis);
#if DUNE_VERSION_LT(DUNE_FUNCTIONS,2,6) boundaryDofsBackend = false;
auto localIndexSet = basis.localIndexSet();
#endif
for (auto it = boundaryPatch.begin(); it != boundaryPatch.end(); ++it)
{
localView.bind(it->inside());
#if DUNE_VERSION_LT(DUNE_FUNCTIONS,2,6)
localIndexSet.bind(localView);
#endif
const auto& localCoefficients = localView.tree().finiteElement().localCoefficients();
for (size_t i=0; i<localCoefficients.size(); i++) auto localView = basis.localView();
auto seDOFs = Dune::Functions::subEntityDOFs(basis);
for(const auto& intersection : boundaryPatch)
{ {
//////////////////////////////////////////////////// localView.bind(intersection.inside());
// Test whether dof is on the boundary face for(auto localIndex: seDOFs.bind(localView, intersection))
//////////////////////////////////////////////////// boundaryDofsBackend[localView.index(localIndex)] = true;
unsigned int entity = localCoefficients.localKey(i).subEntity();
unsigned int codim = localCoefficients.localKey(i).codim();
if (it.containsInsideSubentity(entity, codim))
#if DUNE_VERSION_LT(DUNE_FUNCTIONS,2,6)
boundaryDofs[localIndexSet.index(i)] = true;
#else
boundaryDofs[localView.index(i)] = true;
#endif
}
} }
} }
......
...@@ -8,6 +8,7 @@ set(GRID_BASED_TESTS ...@@ -8,6 +8,7 @@ set(GRID_BASED_TESTS
boundarypatchprolongatortest boundarypatchprolongatortest
coarsegridfunctionwrappertest coarsegridfunctionwrappertest
composedfunctiontest composedfunctiontest
constructboundarydofstest
dunefunctionsassemblertest dunefunctionsassemblertest
dunefunctionsipdgassemblertest dunefunctionsipdgassemblertest
functionintegratortest functionintegratortest
......
#include <config.h>
#include <dune/common/parallel/mpihelper.hh>
#include <dune/common/test/testsuite.hh>
#include <dune/common/bitsetvector.hh>
#include <dune/common/fvector.hh>
#include <dune/fufem/functiontools/boundarydofs.hh>
#include <dune/fufem/boundarypatch.hh>
#include <dune/fufem/test/common.hh>
#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
#include <dune/functions/functionspacebases/lagrangebasis.hh>
#include <dune/functions/functionspacebases/powerbasis.hh>
#include <dune/istl/bvector.hh>
using namespace Dune;
template<int order, class GW>
bool checkConstructBoundaryDofs(const GW& gw)
{
// create active vertices
const auto dim = GW::dimension;
const auto size = gw.size(dim);
BitSetVector<dim> vertices(size);
// set some values
vertices[0].set();
vertices[size / 2].set();
vertices[size-1].set();
// create BoundaryPatch
BoundaryPatch bp(gw,vertices);
// create scalar and powerBasis
using namespace Functions::BasisFactory;
auto lagrangeBasis = makeBasis(gw, lagrange<order>());
auto powerBasis = makeBasis(gw, power<3>(lagrange<order>()));
// create common bitsetvector
BitSetVector<1> lagrangeDOFs;
BitSetVector<3> powerDOFs;
// ... and functions-compatible bitVectors
Dune::BlockVector<Dune::FieldVector<char,1>> lagrangeBits;
Dune::BlockVector<Dune::FieldVector<char,3>> powerBits;
auto lagrangeBackend = Dune::Functions::istlVectorBackend(lagrangeBits);
auto powerBackend = Dune::Functions::istlVectorBackend(powerBits);
// check whether this compiles
constructBoundaryDofs(bp,lagrangeBasis,lagrangeDOFs);
constructBoundaryDofs(bp,powerBasis,powerDOFs);
constructBoundaryDofs(bp,lagrangeBasis,lagrangeBackend);
constructBoundaryDofs(bp,powerBasis,powerBackend);
return true;
}
struct Suite
{
template<class GridType>
bool check(const GridType& grid)
{
auto maxLevel = grid.maxLevel();
auto gridView = grid.levelGridView(maxLevel);
bool passed = checkConstructBoundaryDofs<0>(gridView);
passed = passed and checkConstructBoundaryDofs<1>(gridView);
passed = passed and checkConstructBoundaryDofs<2>(gridView);
return passed;
}
};
int main(int argc, char** argv)
{
Dune::MPIHelper::instance(argc, argv);
Suite tests;
bool passed = checkWithStructuredGrid<2>(tests, 3);
passed = passed and checkWithStructuredGrid<3>(tests, 3);
return passed ? 0 : 1;
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment