Skip to content
Snippets Groups Projects
Praetorius, Simon's avatar
Praetorius, Simon authored
removed several warnings, by including Eigen as a system library

See merge request spraetor/dune-dec!13
927c225b
History

Dune-DEC

A discrete exterior calculus (DEC) framework for the solution of partial differential equations. The toolbox provides essential tools to handle the corresponding meshes and geometric quantities, as well as the assembling of various discrete operators, like the discrete exterior derivative, codifferential, Hodge star, discrete wedge product, discrete flat and sharp operators and thus also higher order operators, like the Laplace-Beltrami operator or the Laplace-deRham operator.

The implementation allows well-centered triangulations and provides a weighted triangulation (see, F. de Goes, P. Memari, P. Mullen, M. Desbrun: "Weighted Triangulations for Geometry Processing").

As references for discrete exterior calculus, see, e.g., Anil N. Hirani: "Discrete Exterior Calculus", PhD Thesis (2003) Link.

Dune Grid

The framework is implemented as dune-grid adapter, i.e., it reads an arbitrary (simplex) grid and builds its own data-structures on top of it.

Example: Create a FoamGrid using an AlbertaReader and from this, create a DecGrid:

using GridBase = Dune::FoamGrid<DEC_DIM, DEC_DOW>;
    
Dune::GridFactory<GridBase> gridFactory;
Dune::AlbertaReader<GridBase>().readGrid(filename, gridFactory);
std::unique_ptr<GridBase> gridBase( gridFactory.createGrid() );

// create a grid adapter
DecGrid<GridBase> grid(*gridBase);

This grid adapter provides iterators for incidence relations, like edges of a vertex or faces of an edge.

Example: Iterate over all edges of the grid (grid-view) and store the coefficients of a Grad-Div laplacian in a sparse matrix:

auto gv = grid.leafGridView();
DOFMatrix<double> A(gv.size(1), gv.size(1)); // sparse matrix with num_rows=#edges

auto const& I = gv.indexSet();

{ auto ins = A.inserter();

// iterator over all edges in the grid
for (auto const& e : edges(gv)) {

  // iterate over all vertices incident to edge e
  for (auto const& v : vertices(e)) {
    auto factor1 = v.sign(e) * gv.dual_volume(v);
    
    // iterate over all edges incident to vertex v
    for (auto const& e_ : edges(v)) {
      auto factor2 = v.sign(e_) * gv.dual_volume(e_) / gv.volume(e_);
      
      ins(I.index(e), I.index(e_)) << factor1 * factor2;
    }
  }
}

} // finish insertion, i.e. compress sparse matrix

Installation

We provide a cmake-based configuration and use the dunecontrol build system. Simply run

dunecontrol --current all
cmake --build build-cmake --target examples

to build all the example problems. The dunecontrol script searches for the required (and suggested) dune modules this library depends on. These include:

(See the file dune.module for an up-to-date list of dependencies). The dune modules can be obtained from https://gitlab.dune-project.org and need to be found in a subdirectory of DUNE_CONTROL_PATH. See also https://dune-project.org/doc/installation for details about the installation of dune modules.

Additionally we require the following libraries to be found:

And a compiler that supports the C++14 standard, e.g. g++ >= 4.9 and clang >= 3.6.

See also the Dokerfiles in iwr/docker-images/dune-latest for a docker container based installation.

Documentation

Currently only a doxygen-based documentation of the source files is available. Generate a html version by

cmake --build build-cmake --target documentation