Code Overview
MFEM consists of the following closely interconnected modules:
- General - general data structures and algorithms.
- Linear Algebra - linear algebra, linear and nonlinear solvers, time steppers.
- Mesh - mesh class, mesh readers, mesh manipulation.
- FEM - finite elements, spaces, linear and bilinear forms, etc.
General
Array
- Class
Array<T>. - Similar to
std::vector<T>in many respects. - Can allocate and manage data or hold external data.
Hash table
- Class
HashTable<T>. - Used in the
NCMeshclass.
Table
- Class
Table. - Maps each row,
{0,1,..,n-1}, to a list of integers (columns). - This is CSR-like data structure without data, only the '
I' and 'J' arrays. - Represents relations like
vertex-to-element,element-to-dof, etc.
Dynamic symmetric table
- Class
DSTableis a dynamic symmetric version of theTableclass. - Used for building relations like
vertex-to-vertex- mesh edges.
Dynamic symmetric 3D table
- Class
STable3D. - Used for building
vertex-to-vertex-to-vertexrelations - mesh faces in 3D.
Communication
- Classes
GroupTopology,GroupCommunicator. - Groups are sets of MPI ranks that need to exchange common data at processor boundaries.
- Reduce (gather) and broadcast (scatter) MPI communications within groups.
Socket stream
- Class
socketstream - Two-way TCP sockets as c++ streams.
- Can be compiled with GnuTLS for security.
- Primarily used for sending data to GLVis.
Timers
- Class
StopWatch. - Has various "backends" -
std::clock, POSIX clocks, Windows'QueryPerformanceCounteretc.
Options parser
- Class
OptionsParser. - Makes it easy to define and parse command line parameters.
- Used in all examples and miniapps.
Linear Algebra
Vector
- Class
Vector- a vector ofdoubles. - Can allocate and manage data or warp external data.
- Defines a number of vector operations on the data.
Operator
- Class
Operator - An abstract base class for all linear and non-linear operators.
- Virtual method
Mult(const Vector &, Vector &). - Optional virtual method
Operator &GetGradient(const Vector &x).
Dense matrix
- Class
DenseMatrix- a matrix ofdoubles. - Can allocate and manage data or wrap external data.
- Uses column-major storage.
- Defines a number of matrix operations, matrix-vector, etc.
- Inherits from
Operator.
Dense tensor
- Class
DenseTensor - Can be viewed as an array of
DenseMatrix(of the same size). - Can be used in batched matrix operations.
Sparse matrix
- Class
SparseMatrix-intindices,doubledata. - Compressed sparse row (CSR) or linked list (LIL) storage.
- Various operations: assembly, matrix-vector, smoothers, etc.
- Inherits from
Operator.
Parallel hypre vector
- Class
HypreParVector- wraps hypre's data structure & operations.
Parallel hypre matrix
- Class
HypreParMatrix- wraps hypre's data structure & operations.
Solvers
- Abstract base class
Solver- extends (inherits) classOperator. - Adds virtual method
void SetOperator(const Operator &).
Direct dense solver
- Class
DenseMatrixInverse- inheritsSolver. - Inverts dense matrices, class
DenseMatrix. - Uses standard LU factorization with pivoting.
Iterative solvers
- Krylov methods, Newton method.
Direct sparse solvers
- Classes
UMFPackSolver,KLUSolver- wraps UMFPACK and KLU from SuiteSparse; can be used withSparseMatrix(serial). - Class
SuperLUSolver- wraps SuperLU_DIST; the parallel matrix needs to be converted toSuperLURowLocMatrix. - Class
STRUMPACKSolver- wraps STRUMPACK; the parallel matrix needs to be converted toSTRUMPACKLURowLocMatrix.
Hypre preconditioners and solvers
- Classes
HypreBoomerAMG,HypreAMS, ect. - Classes
HyprePCG,HypreGMRES.
Time dependent operator
- Class
TimeDependentOperator, inherits fromOperator. - Implements basic virtual methods
double GetTime()andvoid SetTime(const double). - Optional virtual method
void ImplicitSolve(const double dt, const Vector &x, Vector &k)- solve backward Euler system; required for implicit time steppers.
ODE solvers
- Abstract base class
ODESolver. - Has virtual method
void Init(TimeDependentOperator &). - Has pure virtual method
void Step(Vector &x, double &t, double &dt). - Derived classes for explicit Runge-Kutta and implicit (SDIRK) methods.
Symplectic Integrators for Hamiltonian Systems
- Abstract base class 'SIASolver'.
- Has virtual method
void Init(Operator &, TimeDependentOperator &). - Has pure virtual method
void Step(Vector &q, Vector &p, double &t, double &dt). - Derived classes for explicit first and second order integrators.
- Derived class supporting integration orders from 1 to 4.
Constraint operator
- Class
ConstrainedOperator, inherits fromOperator. - Impose essential boundary conditions on any
Operator(matrix-free). - Used by
Operator::FormLinearSystem().
Block vector
- Class
BlockVector, inherits fromVector. - Holds a set of multiple contiguously allocated vectors.
- Useful for systems of equations with multiple components using different finite element spaces.
Block operator
- Class
BlockOperator, inherits fromOperator. - Each block is itself an
Operator.
Block matrix
- Class
BlockMatrix, inheritsAbstractSparseMatrix. - Each block is a
SparseMatrix. - Supports more operations than
BlockOperator.
Block diagonal preconditioner
- Class
BlockDiagonalPreconditioner, inheritsSolver. - Similar to
BlockOperatorbut with diagonal block structure and square diagonal blocks.
Mesh
Mesh
- Class
Mesh. - The mesh topology/connectivity is given by
element-to-vertexrelation. - Elements have type (triangle, quadrilateral, tetrahedron, hexahedron, etc) and
attribute (
int). - Boundary elements can be included allowing tagging of boundary subsets, e.g. for boundary conditions, by the boundary element attribute.
- Edges, faces, and other connectivity are derived automatically based on the element type.
- A high-order mesh uses a vector FE function, i.e. a vector
GridFunction, to represent its high-order nodes. - Hanging/slave vertices are regular vertices - the mesh is "cut" along non-conforming edges and faces.
- Conforming constraints and interpolation are added at the level of the
FiniteElementSpacebased on additional data from theMeshandNCMeshobjects. - Local conforming refinement for triangles and tets.
- Local non-conforming refinement for triangles, quads, and hexes.
- De-refinement and parallel rebalancing for non-conforming meshes.
- Supports curve and surface meshes.
- Periodic meshes: periodic topology with a DG
GridFunctionas nodes, cut along the periodic edges/faces.
Non-conforming mesh
- Class
NCMesh. - Used indirectly through the class
Mesh. - Supports triangles, quads, and hexes including anisotropic refinement for quads and hexes.
- Arbitrary level of hanging nodes and full refinement hierarchy.
- Generates the "cut"
Meshfrom the leaf elements. - See HowTo Article on NC AMR
NURBS mesh
- Class
NURBSExtension. - Used through class
Mesh. - The NURBS patch connectivity is itself a quad/hex
Mesh. - Supports knot insertion, degree elevation, (serial, uniform) h-refinement.
- Generates a quad/hex
Mesh. - Easy to convert to a polynomial high-order mesh.
Mesh readers and writers
- Own formats, read and write, for:
MFEM mesh v1.0,MFEM mesh v1.1(extension for non-conforming meshes),MFEM NURBS mesh v1.0, andMFEM INLINE mesh v1.0(boxes). - Readers for (some) Netgen, TrueGrid, VTK, Gmsh, and CUBIT mesh files.
- Write linear and quadratic VTK meshes.
Parallel mesh
- Classes
ParMesh,ParNCMesh,ParNURBSExtension. - Inherit and extend classes
Mesh,NCMesh,NURBSExtension. ParMeshis constructed from a serialMeshavailable on all tasks.- Built-in mesh partitioning is based on METIS.
- Parallel conforming and non-conforming refinement.
Mesh operators
- Classes
ThresholdRefiner,ThresholdDerefiner,Rebalancer.
Finite Elements
Quadrature formulas
- Class
IntegrationPoint- coordinates plus weights. - Class
IntegrationRule- an array ofIntegrationPoints. - Class
IntegrationRules- container forIntegrationRules.
Element transformation
- Class
ElementTransformation,IsoparametricTransformation. - Defined through a
FiniteElement - Transforms reference
IntegrationPoints into physicalVectors. - On demand, computes and stores Jacobian matrix, and weight.
Finite elements
- Abstract base class
FiniteElement. - Arbitrary order
H1_*,L2_*,RT_*, andND_*elements on segment, triangles, quads, tets, and hexes. - Abstract method
CalcShape,CalcDShapefor scalar FE;CalcVShape,CalcDivShape,CalcCurlShapefor vector H(div)/H(curl) FE. - Other interpolation and projection methods.
Finite element collections
- Base class
FiniteElementCollection. - Associates
FiniteElements with elements, faces, edges, vertices. - Degrees of freedom on faces/edges/vertices are shared between adjacent elements.
- Derived classes for arbitrary order
H1_*,L2_*,RT_*, andND_*collections.
Finite element space
- Class
FiniteElementSpace. - Constructed from a
Meshand aFiniteElementCollection. - Defines the mappings
elements-to-dofs,faces-to-dofs, etc. - Defines, when necessary, a prolongation,
P, and a restriction,R, matrices:R.P = Ito constrain "slave" dofs. - On non-conforming meshes, the space is "cut" or "partially conforming"
(before applying
P). - The domain of
Pdefines the "true" or "conforming" dofs. - The range of
Pis a sub-set of the "partially conforming" dofs.
Grid function
- Class
GridFunction, extends classVector. - A container
Vectoron the "partially conforming" dofs. - Defines a number of useful operations like computing values, gradient, etc
at quadrature points (
IntegrationPointorIntegrationRule). - Methods for projecting (interpolating)
Coefficient,VectorCoefficient. - Methods for computing error norms with respect to a
Coefficient.
Linear form
- Class
LinearForm, extends classVector. - Assembles r.h.s. vector.
- Uses a sum of local
LinearFormIntegrators.
Bilinear form
- Class
BilinearForm. - Assembles linear system matrix.
- Uses a sum of local
BilinearFormIntegrators. - Method
FormLinearSystemapplies necessary transformations, e.g. $P^T A P$.
Mixed bilinear form
- Class
MixedBilinearForm.
Coefficients
- Abstract base classes:
Coefficient,VectorCoefficient, andMatrixCoefficient. - Derived classes include:
ConstantCoefficient,FunctionCoefficient,GridFunctionCoefficient;VectorFunctionCoefficient,VectorGridFunctionCoefficient, etc. - Easy to derive new coefficient classes.
Parallel versions
ParFiniteElementSpaceParGridFunctionParLinearFormParBilinearForm- etc
Error estimators
- Classes
ZienkiewiczZhuEstimator,L2ZienkiewiczZhuEstimator.