HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
HArD::Core2D Documentation

Table of Contents

HArD::Core (sources: https://github.com/jdroniou/HArDCore) provides a suite of C++ tools to implement numerical schemes whose unknowns are polynomials in the cells, on the edges, and on the faces. The focus is on dealing on generic polytopal meshes. This document addresses the 2D version of HArD::Core. The tools are similar in 2D and 3D, and we refer to the main page of documentation of the 3D version for a more thorough introduction to the library. We only present here the specific way of loading a 2D mesh structure (which differs from the 3D version), and a brief description of the schemes available in this 2D library.

Transferring a scheme's implementation from 3D to 2D or vice-versa is very straightforward, provided that the scheme's mathematical definition does not depend on the dimension and that the generic types provided in basis.hpp are used; see README file of the HArD::Core github depository https://github.com/jdroniou/HArDCore.

  • Loading a 2D mesh – How to load a mesh.
  • Schemes – The list of schemes currently implemented in HArD::Core2D, and scripts to run them.

Loading a mesh

HArDCore2D currently reads meshes in the typ2 format designed for the FVCA5 Benchmark. A short documentation describing this format is provided in the typ2_meshes directory (see README.pdf therein). Several meshes can also be found in this directory.

A mesh file must be read using an instance of the MeshReaderTyp2 class, and then built using MeshBuilder. A working example is given below (assuming the executable will be in build/Schemes for example).

#include "mesh.hpp"
#include "import_mesh.hpp"
#include "mesh_builder.hpp"
using namespace HArDCore2D;
int main() {
// Build mesh
MeshBuilder builder = MeshBuilder(mesh_file);
std::unique_ptr<Mesh> mesh_ptr = builder.build_the_mesh();
// Get the BC and re-order the edges (useful to set BC for hybrid schemes)
std::string bc_id = vm["bc_id"].as<std::string>();
BoundaryConditions BC(bc_id, *mesh_ptr.get());
BC.reorder_edges();
// Create an HybridCore instance
HybridCore hho(mesh_ptr.get(), K+1, K, use_threads, output);
}
The BoundaryConditions class provides definition of boundary conditions.
Definition: BoundaryConditions.hpp:45
std::string bc_id
Definition: HHO_DiffAdvecReac.hpp:44
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
size_t K
Definition: HHO_DiffAdvecReac.hpp:46
int main()
Definition: compute_rates.cpp:8
Definition: ddr-klplate.hpp:27

Note: the typ2 format allows for meshes with very generic polygonal cells, including non-convex cells. However, the builder assumes that each cell is star-shaped with respect to the isobarycenter of its vertices – otherwise, the calculation of the center of mass may be incorrect. Similarly, the quadrature rules assume that each cell is star-shaped with respect to its center of mass.

Schemes

The following schemes are currently available in HArD::Core2D. The Hybrid High-Order schemes follow the implementation principles described in Appendix B of the book available at https://hal.archives-ouvertes.fr/hal-02151813.

  • HHO_diffusion: Hybrid High-Order (HHO) for \(-\mathrm{div}(K\nabla u)=f\), for Dirichlet, Neumann or mixed boundary conditions, with \(K\) a diffusion tensor that is piecewise constant on the mesh.
  • HHO_locvardiff: HHO for \(-\mathrm{div}(K\nabla u)=f\), for Dirichlet, Neumann or mixed boundary conditions, with \(K\) a diffusion tensor that can vary in each cell.
  • HHO_diffadvecreac: Hybrid High-Order (HHO) for \(-\mathrm{div}(K\nabla u+\beta u)+\mu u=f\), for Dirichlet or mixed boundary conditions, with \(K\) a diffusion tensor that is piecewise constant on the mesh.
  • LEPNC_diffusion, in module LEPNC: Locally Enriched Polytopal Non-Conforming (LEPNC) method for the pure diffusion problem \(-\mathrm{div}(K\nabla\zeta(u))=f\).
  • LEPNC_StefanPME, in module LEPNC: Locally Enriched Polytopal Non-Conforming (LEPNC) method for the stationnary Stefan/PME problem \(u-\mathrm{div}(K\nabla\zeta(u))=f\).
  • LEPNC_StefanPME_Transient, in module LEPNC: LEPNC for the transient Stefan/PME problem \(\partial_t u-\mathrm{div}(K\nabla\zeta(u))=f\) (may not be working any longer due to changes in the management of quadrature rules)
  • HMM_StefanPME_Transient, in module HMM: Hybrid Mimetic Mixed (HMM) method for the transient Stefan/PME problem \(\partial_t u-\mathrm{div}(K\nabla\zeta(u))=f\).
  • DDR_rmplate: Discrete de Rham (DDR) scheme for the Reissner-Mindlin plate bending problem.
  • DDR_klplate: Scheme for the Kirchhoff-Love plate problem based on the DDR principles (with design of a specific complex for plates).
  • HHO_fullgradientdiff: HHO scheme with full gradient; similar to HHO_locvardiff but implemented using the HHOSpace module instead of the HybridCore module.

Besides the schemes, some codes/libraries provide generic useful tools:

  • BoundaryConditions: tools to handle boundary conditions (selection of boundary mesh quantities for natural/essential BCs, tools to create maps DOFs to Unkowns (removing certains DOFs corresponding to essential boundary conditions), etc.).
  • MeshCoarsen: tools to coarsen (agglomerate) a mesh. Mostly one executable (mesh-coarsen) is useful.

The directory runs contains BASH to run series of tests on families of meshes. The files data.sh describe the parameters of the test cases (polynomial degrees, boundary conditions, mesh families, etc.). The script produces results in the output directory, including a pdf file rate.pdf describing the rates of convergence in various energy norms.

To run the scripts as they are, you will need pdflatex and (for the LEPNC and HMM schemes) a FORTRAN compiler (adjust the Makefile to your compiler).