Rodin is a lightweight and modular finite element framework which provides many of the associated functionalities that are needed when implementing shape and topology optimization algorithms. These functionalities range from refining and remeshing the underlying shape, to providing elegant mechanisms to specify and solve variational problems.
It is named after the French sculptor Auguste Rodin, considered the founder of modern sculpture.
The library is still in development. It is primarily maintained by Carlos Brito-Pacheco and was developed to generate examples for his ongoing PhD.
Any contributors are warmly encouraged and any help or comments are always appreciated!
Branch | Matrix | Tests | Code Coverage | Benchmarks | Documentation |
---|---|---|---|---|---|
master | |||||
develop |
- Building the project
- Features
- Third-Party integrations
- Requirements
- CMake options
- Building the documentation
git clone --recursive https://github.com/carlos-brito-pacheco/rodin
cd rodin
mkdir build && cd build
cmake ..
make -j4
Rodin comes with documentation that is built automatically on each merge, hence it's always up to date.
The documentation can be found here.
Rodin comes with a native C++20 form language for assembling and solving variational formulations.
For example, given a domain
has the associated weak formulation:
which can be quickly implemented via the following lines of code:
#include <Rodin/Solver.h>
#include <Rodin/Geometry.h>
#include <Rodin/Variational.h>
using namespace Rodin;
using namespace Rodin::Geometry;
using namespace Rodin::Variational;
int main(int, char**)
{
Mesh Omega;
Omega = Omega.UniformGrid(Polytope::Type::Triangle, 16, 16);
mesh.getConnectivity().compute(1, 2);
P1 Vh(Omega);
TrialFunction u(Vh);
TestFunction v(Vh);
Solver::SparseLU solver;
Problem poisson(u, v);
poisson = Integral(Grad(u), Grad(v))
- Integral(v)
+ DirichletBC(u, Zero());
poisson.solve(solver);
return 0;
}
Solution of the Poisson equation. |
The API offers full support for iteration over all polytopes of the mesh of some given dimension:
Mesh mesh;
mesh = mesh.UniformGrid(Polytope::Type::Triangle, 16, 16); // 2D Mesh
for (auto it = mesh.getCell(); it; ++it)
{
// Access information about the cell
}
for (auto it = mesh.getFace(); it; ++it)
{
// Access information about the face
}
for (auto it = mesh.getVertex(); it; ++it)
{
// Access information about the vertex
}
for (auto it = mesh.getPolytope(1); it; ++it)
{
// Access information about the face (face dimension in 2D is equal to 1)
}
Rodin is able to compute any connectivity information on the mesh. For example, the following computes the adjacency information from faces to cells:
Mesh mesh;
mesh = mesh.UniformGrid(Polytope::Type::Triangle, 16, 16); // 2D Mesh
mesh.getConnectivity().compute(1, 2);
In general, this means that given a face, we are able to obtain the incident (neighboring) cells.
However, one can also compute any connectivity information on different dimensions.
For example, for a mesh
// Compute connectivity between vertices and faces
// i.e. Given a vertex, give me the incident edges
mesh.getConnectivity().compute(0, 1);
// Compute connectivity between faces and cells
// i.e. Given a vertex, give me the incident cells
mesh.getConnectivity().compute(0, 2);
// Compute connectivity between faces
// i.e. Given a face, give me the adjacent faces
mesh.getConnectivity().compute(1, 1);
// Compute connectivity between cells
// i.e. Given a cell, give me the adjacent cells
mesh.getConnectivity().compute(2, 2);
// Compute connectivity between cells and faces
// i.e. Given a cell, give me the adjacent faces
mesh.getConnectivity().compute(2, 1);
// Etc.
- MFEM
- MEDIT
Rodin supports different kinds of quadrature.
- Grundmann-Moeller
MMG is an open source software for bidimensional and tridimensional surface and volume remeshing.
-
Loading the mesh:
MMG::Mesh Omega; Omega.load(meshFile, IO::FileFormat::MEDIT);
-
Optimizing the mesh:
MMG::Optimizer().setHMax(hmax) // maximal edge size .setHMin(hmin) // minimal edge size .setGradation(hgrad) // ratio between two edges .setHausdorff(hausd) // curvature refinement .optimize(Omega);
List of features and modules that are in the works:
- Discontinuous Galerkin methods
Rodin::Plot
module- H1
- L2
- HDiv
- HCurl
- P2
- P0
- PETSc
- METIS
Any of these should be available for quick install from your standard package manager.
Option | Description |
---|---|
RODIN_BUILD_EXAMPLES | Builds the examples in the examples/ directory. |
RODIN_BUILD_DOC | Builds the documentation using Doxygen |
RODIN_USE_MCSS | Builds the documentation using Doxygen and m.css |
RODIN_BUILD_SRC | Build the Rodin source code |
RODIN_BUILD_EXAMPLES | Build the Rodin examples |
RODIN_BUILD_DOC | Build the Rodin documentation |
RODIN_USE_MCSS | Use m.css style documentation |
RODIN_WITH_PLOT | Build the Rodin::Plot module |
RODIN_USE_MPI | Build with MPI support |
RODIN_USE_OPENMP | Build with OpenMP support |
RODIN_USE_SUITESPARSE | Build with SuiteSparse support |
RODIN_SILENCE_WARNINGS | Silence warnings outputted by Rodin |
RODIN_BUILD_PY | Build Python bindings |
See this page to see how to build the documentation.