The following explains how to solve the Poisson equation using dune-functions. The full example poisson-pq2.cc can be found in the examples/ subdirectory.
Local assemblers
The program first includes the required headers and then defines free functions for assembling the Poisson problem. The function getLocalMatrix() implements the assembler of the local stiffness matrix for the bilinear form
.
template <class LocalView, class MatrixType>
void getLocalMatrix(const LocalView& localView, MatrixType& elementMatrix)
{
typedef typename LocalView::Element Element;
const Element& element = localView.element();
const int dim = Element::dimension;
auto geometry = element.geometry();
const auto& localFiniteElement = localView.tree().finiteElement();
elementMatrix.setSize(localFiniteElement.localBasis().size(),localFiniteElement.localBasis().size());
elementMatrix = 0;
int order = 2*(dim*localFiniteElement.localBasis().order()-1);
const QuadratureRule<double, dim>& quad = QuadratureRules<double, dim>::rule(element.type(), order);
for (size_t pt=0; pt < quad.size(); pt++) {
const FieldVector<double,dim>& quadPos = quad[pt].position();
const auto& jacobianInverse = geometry.jacobianInverse(quadPos);
const double integrationElement = geometry.integrationElement(quadPos);
std::vector<FieldMatrix<double,1,dim> > referenceJacobians;
localFiniteElement.localBasis().evaluateJacobian(quadPos, referenceJacobians);
std::vector<FieldMatrix<double,1,dim> > jacobians(referenceJacobians.size());
for (size_t i=0; i<jacobians.size(); i++)
jacobians[i] = referenceJacobians[i] * jacobianInverse;
for (size_t i=0; i<elementMatrix.N(); i++)
for (size_t j=0; j<elementMatrix.M(); j++ )
elementMatrix[i][j] += (jacobians[i] * transpose(jacobians[j])) * quad[pt].weight() * integrationElement;
}
}
The getVolumeTerm() functions implements the local assembler for the volume right hand side term
.
template <class LocalView, class LocalVolumeTerm>
void getVolumeTerm( const LocalView& localView,
BlockVector<double>& localRhs,
LocalVolumeTerm&& localVolumeTerm)
{
typedef typename LocalView::Element Element;
const Element& element = localView.element();
const int dim = Element::dimension;
const auto& localFiniteElement = localView.tree().finiteElement();
localRhs.resize(localFiniteElement.localBasis().size());
localRhs = 0;
int order = dim*localFiniteElement.localBasis().order();
const QuadratureRule<double, dim>& quad = QuadratureRules<double, dim>::rule(element.type(), order);
for ( size_t pt=0; pt < quad.size(); pt++ ) {
const FieldVector<double,dim>& quadPos = quad[pt].position();
const double integrationElement = element.geometry().integrationElement(quadPos);
double functionValue = localVolumeTerm(quadPos);
std::vector<FieldVector<double,1> > shapeFunctionValues;
localFiniteElement.localBasis().evaluateFunction(quadPos, shapeFunctionValues);
for (size_t i=0; i<localRhs.size(); i++)
localRhs[i] += shapeFunctionValues[i] * functionValue * quad[pt].weight() * integrationElement;
}
}
Global assembler
Assemble the global matrix pattern.
template <class FEBasis>
void getOccupationPattern(const FEBasis& feBasis, MatrixIndexSet& nb)
{
auto n = feBasis.size();
nb.resize(n, n);
auto localView = feBasis.localView();
for(
const auto& e :
elements(feBasis.gridView()))
{
localView.bind(e);
for (size_t i=0; i<localView.tree().size(); i++) {
for (size_t j=0; j<localView.tree().size(); j++) {
auto iIdx = localView.index(i);
auto jIdx = localView.index(j);
nb.add(iIdx, jIdx);
}
}
}
}
Assembly of matrix and right-hand-side.
template <class FEBasis, class VolumeTerm>
void assembleLaplaceMatrix(const FEBasis& feBasis,
BCRSMatrix<double>& matrix,
BlockVector<double>& rhs,
VolumeTerm&& volumeTerm)
{
typedef typename FEBasis::GridView GridView;
GridView gridView = feBasis.gridView();
MatrixIndexSet occupationPattern;
getOccupationPattern(feBasis, occupationPattern);
occupationPattern.exportIdx(matrix);
rhs.resize(feBasis.size());
matrix = 0;
rhs = 0;
auto localView = feBasis.localView();
{
localView.bind(e);
Matrix<double> elementMatrix;
getLocalMatrix(localView, elementMatrix);
for (size_t i=0; i<elementMatrix.N(); i++) {
auto row = localView.index(i);
for (size_t j=0; j<elementMatrix.M(); j++ ) {
auto col = localView.index(j);
matrix[row][col] += elementMatrix[i][j];
}
}
BlockVector<double> localRhs;
localVolumeTerm.bind(e);
getVolumeTerm(localView, localRhs, localVolumeTerm);
for (size_t i=0; i<localRhs.size(); i++) {
auto row = localView.index(i);
rhs[row] += localRhs[i];
}
}
}
Helper functions
Treatment of boundary condition.
template <class FEBasis>
void boundaryTreatment (const FEBasis& feBasis, std::vector<char>& dirichletNodes )
{
dirichletNodes.clear();
dirichletNodes.resize(feBasis.size(), false);
dirichletNodes[index] = true;
});
}
Create a mixed grid containing triangles and quadrilaterals.
auto createMixedGrid()
{
using Grid = Dune::UGGrid<2>;
Dune::GridFactory<Grid> factory;
for(unsigned int k : Dune::range(9))
factory.insertVertex({0.5*(k%3), 0.5*(k/3)});
factory.insertElement(Dune::GeometryTypes::cube(2), {0, 1, 3, 4});
factory.insertElement(Dune::GeometryTypes::cube(2), {1, 2, 4, 5});
factory.insertElement(Dune::GeometryTypes::simplex(2), {3, 4, 6});
factory.insertElement(Dune::GeometryTypes::simplex(2), {4, 7, 6});
factory.insertElement(Dune::GeometryTypes::simplex(2), {4, 5, 7});
factory.insertElement(Dune::GeometryTypes::simplex(2), {5, 8, 7});
return factory.createGrid();
}
Setup and initialization
Initialize MPI
int main (int argc, char *argv[]) try
{
MPIHelper::instance(argc, argv);
Create grid and finite element basis
Create a mixed grid and obtain a leaf grid view.
- Note
- A GridView is a view to a subset of the grid's elements, vertices, ... that should be stored by value and can be copied cheaply. Grids in dune are in general hierarchical and composed by elements on several levels. The discretization usually lives on the set of most refined elements that is denote the leaf GridView in dune.
auto grid = createMixedGrid();
grid->globalRefine(4);
auto gridView = grid->leafGridView();
using GridView = decltype(gridView);
As a next step the program creates a global finite element basis on the GridView.
FEBasis feBasis(gridView);
And here is the rest of the file:
using VectorType = BlockVector<double>;
using MatrixType = BCRSMatrix<double>;
VectorType rhs;
MatrixType stiffnessMatrix;
std::cout << "Number of DOFs is " << feBasis.dimension() << std::endl;
auto rightHandSide = [] (const auto& x) { return 10;};
Dune::Timer timer;
assembleLaplaceMatrix(feBasis, stiffnessMatrix, rhs, rightHandSide);
std::cout << "Assembling the problem took " << timer.elapsed() << "s" << std::endl;
VectorType x(feBasis.size());
x = 0;
std::vector<char> dirichletNodes;
boundaryTreatment(feBasis, dirichletNodes);
auto pi = std::acos(-1.0);
auto dirichletValueFunction = [pi](const auto& x){ return std::sin(2*pi*x[0]); };
interpolate(feBasis, x, dirichletValueFunction, dirichletNodes);
stiffnessMatrix.mmv(x, rhs);
for (size_t i=0; i<stiffnessMatrix.N(); i++) {
if (dirichletNodes[i]) {
rhs[i] = x[i];
for (auto&& [entry, idx] : sparseRange(stiffnessMatrix[i]))
entry = (i==idx) ? 1.0 : 0.0;
}
else {
for (auto&& [entry, idx] : sparseRange(stiffnessMatrix[i]))
if (dirichletNodes[idx])
entry = 0.0;
}
}
MatrixAdapter<MatrixType,VectorType,VectorType> op(stiffnessMatrix);
SeqILDL<MatrixType,VectorType,VectorType> ildl(stiffnessMatrix,1.0);
CGSolver<VectorType> cg(op,
ildl,
1e-4,
100,
2);
InverseOperatorResult statistics;
cg.apply(x, rhs, statistics);
SubsamplingVTKWriter<GridView> vtkWriter(gridView, Dune::refinementLevels(2));
vtkWriter.addVertexData(xFunction, VTK::FieldInfo("x", VTK::FieldInfo::Type::scalar, 1));
vtkWriter.write("poisson-pq2");
}
catch (Exception& e) {
std::cout << e.what() << std::endl;
}