Variational namespace
Module which provides the necessary tools for constructing and solving variational problems.
The Variational module is the core of the Rodin finite element framework. It provides classes for defining finite element spaces, trial and test functions, variational forms, boundary conditions, and problems.
Finite Element Spaces
Functions
- TrialFunction — Represents the unknown solution in a variational formulation
- TestFunction — Represents the weighting function
- GridFunction — A concrete function storing DOF values in a finite element space
Variational Forms
The module supports building variational problems of the form:
which are translated into linear algebra problems .
Key expression classes:
Integral(A, B)— Domain integrals (bilinear or linear forms)BoundaryIntegral(A, B)— Boundary integralsGrad(u)— Gradient operatorDiv(u)— Divergence operatorJacobian(u)— Jacobian operatorDirichletBC(u, g)— Dirichlet boundary conditions
Problem
The Problem class assembles all forms and boundary conditions into a solvable linear system.
P1 Vh(mesh); TrialFunction u(Vh); TestFunction v(Vh); Problem poisson(u, v); poisson = Integral(Grad(u), Grad(v)) - Integral(f, v) + DirichletBC(u, Zero()); Solver::CG(poisson).solve();
Namespaces
- namespace F
- Contains built-in coordinate functions.
Classes
-
template<class Operand>class Abs
- Represents the absolute value of a value.
-
template<class NestedDerived>class Abs<FunctionBase<NestedDerived>>
- Absolute value operation for functions.
-
template<class LHS, class RHS>class AND
- Represents the logical AND expression.
-
template<class LHSDerived, class RHSDerived>class AND<BooleanFunctionBase<LHSDerived>, BooleanFunctionBase<RHSDerived>>
- Logical AND operator between two boolean functions.
-
template<class FunctionDerived>class Average<FunctionBase<FunctionDerived>>
- Average operator for computing interface averages.
-
template<class Solution, class TrialFES, class TestFES, class OperatorType>class BilinearForm
- Represents a bilinear form on a trial and test space.
-
template<class Solution, class TrialFES, class TestFES>class BilinearForm<Solution, TrialFES, TestFES, ::Mat>
- Bilinear form specialization that assembles into a PETSc matrix.
-
template<class Solution, class TrialFES, class TestFES, class Scalar>class BilinearForm<Solution, TrialFES, TestFES, Math::SparseMatrix<Scalar>>
- Speciallization of BilinearForm for a matrix type.
-
template<class Operator>class BilinearFormBase
- Base class for bilinear form representations.
-
template<class Number, class Derived>class BilinearFormIntegratorBase
- Abstract base class for bilinear form integrators.
-
template<class T>class BooleanFunction
-
template<class Derived>class BooleanFunctionBase
- Base class for boolean-valued functions.
- struct BoundaryHit
- Information passed to a boundary policy when a trajectory reaches the domain boundary.
-
template<class Integrand>class BoundaryIntegral
- Represents expressions of the integral operator on the boundary of a domain.
-
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>class BoundaryIntegral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>>
- Integration of the dot product of a trial and test operators.
-
template<class NestedDerived, class FES>class BoundaryIntegral<ShapeFunctionBase<NestedDerived, FES, TestSpace>>
- Integration of a test operator.
- class BoundaryNormal
- Outward unit normal vector on domain boundaries.
-
template<class ... Values>class ComplexFunction
-
template<>class ComplexFunction<Complex>
- Represents a constant scalar function with type Complex.
-
template<class F>class ComplexFunction<F>
- Represents a scalar function given by an arbitrary scalar function.
-
template<class Derived>class ComplexFunctionBase
- Base class for complex-valued scalar functions.
-
template<class OperandDerived>class Component<FunctionBase<OperandDerived>, size_t>
- Extracts a single component from a vector-valued function.
-
template<class OperandDerived>class Component<FunctionBase<OperandDerived>, size_t, size_t>
- Extracts a single entry from a matrix-valued function.
-
template<class FES, class Data>class Component<GridFunction<FES, Data>>
- Extracts a component from a vector-valued GridFunction.
-
template<class OperandDerived, class FES, ShapeFunctionSpaceType Space>class Component<ShapeFunctionBase<OperandDerived, FES, Space>>
- Extracts a component from a vector-valued ShapeFunction.
-
template<class LHS, class RHS>class Composition
-
template<class NestedDerived>class Conjugate<FunctionBase<NestedDerived>>
- Complex conjugate of a function.
-
template<class NestedDerived, class FES, ShapeFunctionSpaceType Space>class Conjugate<ShapeFunctionBase<NestedDerived, FES, Space>>
- Specialization for shape functions.
-
template<class Operand>class Cos
- Represents the cosine function.
-
template<class NestedDerived>class Cos<FunctionBase<NestedDerived>>
- Cosine function operator for real-valued scalar functions.
-
template<class NestedDerived>class Cosh<FunctionBase<NestedDerived>>
- Hyperbolic cosine function operator for real-valued scalar functions.
- class DefaultBoundaryPolicy
- Default boundary policy for flow maps.
-
template<class LinearSystem, class U, class V>class DenseProblem<LinearSystem, U, V>
- General class to assemble linear systems with
Math::andMatrix Math::types in a serial context.Vector -
template<class Range, class Data, class Mesh>class Derivative<GridFunction<P1<Range, Mesh>, Data>>
- Derivative of a P1 GridFunction.
-
template<class Operand, class Derived>class DerivativeBase
- Base class for directional derivative operators.
-
template<class Operand, class Value>class DirichletBC
-
template<class Solution, class FES, class ValueDerived>class DirichletBC<TrialFunction<Solution, FES>, FunctionBase<ValueDerived>>
- Represents a Dirichlet boundary condition on a ShapeFunction object.
-
template<class Scalar>class DirichletBCBase
- Abstract base class for a Dirichlet boundary condition.
-
template<class Operand>class Div
- Represents the divergence of a vector valued function.
-
template<size_t K, class Scalar, class Data, class Mesh>class Div<GridFunction<H1<K, Math::Vector<Scalar>, Mesh>, Data>>
- Divergence of an H1 vector GridFunction.
-
template<class Scalar, class Data, class Mesh>class Div<GridFunction<P0g<Math::Vector<Scalar>, Mesh>, Data>>
- Divergence of a P0g vector GridFunction (identically zero).
-
template<class Scalar, class Data, class Mesh>class Div<GridFunction<P1<Math::Vector<Scalar>, Mesh>, Data>>
- Divergence of a P1 vector GridFunction.
-
template<size_t K, class NestedDerived, class Number, class Mesh, ShapeFunctionSpaceType Space>class Div<ShapeFunction<NestedDerived, H1<K, Math::Vector<Number>, Mesh>, Space>>
- Divergence of an H1 vector ShapeFunction.
-
template<class NestedDerived, class Scalar, class Mesh, ShapeFunctionSpaceType Space>class Div<ShapeFunction<NestedDerived, P0g<Math::Vector<Scalar>, Mesh>, Space>>
- Divergence of a P0g vector ShapeFunction (identically zero).
-
template<class NestedDerived, class Scalar, class Mesh, ShapeFunctionSpaceType Space>class Div<ShapeFunction<NestedDerived, P1<Math::Vector<Scalar>, Mesh>, Space>>
- Divergence of a P1 vector ShapeFunction.
-
template<class Operand, class Derived>class DivBase
- Base class for divergence operator implementations.
-
template<class FES, class Data, class Derived>class DivBase<GridFunction<FES, Data>, Derived>
- Divergence of a P1 GridFunction.
-
template<class LHS, class RHS>class Division
- Represents the division operation.
-
template<class LHSDerived, class RHSDerived>class Division<FunctionBase<LHSDerived>, FunctionBase<RHSDerived>>
- Division of a function by another function.
-
template<class LHS, class RHS>class Dot
- Represents the dot product between two objects.
-
template<class LHSDerived, class RHSDerived>class Dot<FunctionBase<LHSDerived>, FunctionBase<RHSDerived>>
- Dot product between two functions.
-
template<class LHSDerived, class RHSDerived, class FES, ShapeFunctionSpaceType Space>class Dot<FunctionBase<LHSDerived>, ShapeFunctionBase<RHSDerived, FES, Space>>
- Dot product between a function and a shape function.
-
template<class KernelType, class LHSDerived, class TrialFES, class RHSDerived, class TestFES>class Dot<Potential<KernelType, ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>
- Dot product of a potential with a test shape function.
-
template<class LHSDerived, class RHSDerived, class FES, ShapeFunctionSpaceType Space>class Dot<ShapeFunctionBase<LHSDerived, FES, Space>, FunctionBase<RHSDerived>>
- Dot product between a shape function and a function.
-
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>class Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>
- Dot product of trial and test shape functions for bilinear forms.
-
template<size_t K>class DubinerTetrahedron
- Dubiner orthogonal modal basis on the reference tetrahedron.
-
template<size_t K>class DubinerTriangle
- Dubiner orthogonal modal basis on the reference triangle.
-
template<class LHS, class RHS>class EQ
-
template<class LHSDerived, class RHSDerived>class EQ<FunctionBase<LHSDerived>, FunctionBase<RHSDerived>>
- Logical EQ operator between two instances of FunctionBase.
-
template<class NestedDerived>class Exp<FunctionBase<NestedDerived>>
- Exponential function operation.
-
template<class Integrand>class FaceIntegral
- Represents expressions of the integral operator on the faces of a mesh.
-
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>class FaceIntegral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>>
- Integration of the dot product of a trial and test operators.
-
template<class NestedDerived, class FES>class FaceIntegral<ShapeFunctionBase<NestedDerived, FES, TestSpace>>
- Integration of a test operator.
- class FaceNormal
- Unit normal vector on mesh faces.
-
template<size_t K>class FeketeTetrahedron
- Cached Fekete-type nodes on the reference tetrahedron for degree K.
-
template<size_t K>class FeketeTriangle
- Cached Fekete-type nodes on the reference triangle for degree K.
-
template<class Derived>class FiniteElementBase
- Base class for finite elements.
-
template<class Mesh, class Derived>class FiniteElementSpace
- Represernts a finite element space.
-
template<class Derived>class FiniteElementSpace<Geometry::Mesh<Context::MPI>, Derived>
- Represents a finite element space defined on a distributed mesh.
- class FiniteElementSpaceBase
- Base class for finite element spaces.
-
template<class Derived>class FiniteElementSpacePullbackBase
- Base class for mappings taking functions defined on physical elements to reference elements.
-
template<class Derived>class FiniteElementSpacePushforwardBase
- Base class for inverse mappings taking functions defined on reference elements to physical elements.
-
template<class Derived, class VectorField, class Step, class BoundaryPolicy>class Flow<FunctionBase<Derived>, VectorField, Step, BoundaryPolicy>
- Flow map operator for functions.
-
template<class Operand>class Frobenius
- Represents the Frobenius norm.
-
template<class NestedDerived>class Frobenius<FunctionBase<NestedDerived>>
- Frobenius norm of a function.
-
template<class Derived>class FunctionBase
- Base class for function objects which can be evaluated over a mesh.
-
template<class StrictType>class FunctionBaseCopy
-
template<class LHS, class RHS>class GEQ
-
template<size_t K>class GLL
- Compile-time Gauss–Lobatto–Legendre (GLL) nodes on [-1, 1].
-
template<size_t K>class GLL01
- Gauss–Lobatto–Legendre nodes mapped from [-1, 1] to [0, 1].
-
template<class Operand>class Grad
- Represents the gradient of a scalar function .
-
template<size_t K, class Scalar, class Mesh, class Data>class Grad<GridFunction<H1<K, Scalar, Mesh>, Data>>
- Gradient of a GridFunction on H1<K> space.
-
template<class Range, class Data, class Mesh>class Grad<GridFunction<P0<Range, Mesh>, Data>>
- Gradient of a P0 GridFunction.
-
template<size_t K, class NestedDerived, class Scalar, class Mesh, ShapeFunctionSpaceType SpaceType>class Grad<ShapeFunction<NestedDerived, H1<K, Scalar, Mesh>, SpaceType>>
- Gradient of a ShapeFunction on H1<K> space.
-
template<class NestedDerived, class Scalar, class Mesh, ShapeFunctionSpaceType SpaceType>class Grad<ShapeFunction<NestedDerived, P0<Scalar, Mesh>, SpaceType>>
- Gradient of a P0 ShapeFunction.
-
template<class NestedDerived, class Scalar, class Mesh, ShapeFunctionSpaceType SpaceType>class Grad<ShapeFunction<NestedDerived, P1<Scalar, Mesh>, SpaceType>>
- Gradient of a ShapeFunction.
-
template<class Operand, class Derived>class GradBase
- Base class for gradient operator implementations.
-
template<class FES, class Data, class Derived>class GradBase<GridFunction<FES, Data>, Derived>
- Gradient operator for grid functions.
-
template<class FES, class Data>class GridFunction
- Represents a grid function belonging to some finite element space.
-
template<class FES>class GridFunction<FES, ::Vec>
- Grid function specialization storing DOF coefficients in a PETSc vector (
Vec). -
template<class Derived, class FES = typename FormLanguage::class GridFunctionBase
Traits<Derived>::FESType, class Data = typename FormLanguage:: Traits<Derived>::DataType> - Base class for grid function objects.
-
template<class Derived>class GridFunctionBaseReference
- Base class for discrete finite element functions.
-
template<class LHS, class RHS>class GT
-
template<size_t K, class Range, class Mesh>class H1
- Degree K H1-conforming Lagrange finite element space.
-
template<size_t K, class Scalar>class H1<K, Math::Vector<Scalar>, Geometry::Mesh<Context::Local>>
- Vector-valued continuous piecewise polynomial (degree K) Lagrange finite element space.
-
template<size_t K, class Scalar>class H1<K, Scalar, Geometry::Mesh<Context::Local>>
- Degree K H1-conforming Lagrange finite element space.
-
template<size_t K, class Scalar>class H1Element
- Degree k H1-conforming Lagrange element.
-
template<size_t K, class Scalar>class H1Element<K, Math::Vector<Scalar>>
- Continuous H1-conforming piecewise polynomial (degree k) vector Lagrange element.
- class IdentityMatrix
- Represents the identity matrix function .
-
template<class Operand>class Im
- Extracts the imaginary part of a complex-valued function.
-
template<class NestedDerived>class Im<FunctionBase<NestedDerived>>
- Specialization for FunctionBase operands.
-
template<class T>class IntegerFunction
-
template<class Derived>class IntegerFunctionBase
- Base class for objects representing integer functions.
-
template<class Integrand>class Integral
- Represents mathematical expressions of the integral operator on a domain.
-
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>class Integral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>>
- Integration of the dot product of a trial and test operators.
-
template<class FES, class Data>class Integral<GridFunction<FES, Data>>
- Integration of a GridFunction object.
-
template<class NestedDerived, class FES>class Integral<ShapeFunctionBase<NestedDerived, FES, TestSpace>>
- Integration of a test operator.
- class Integrator
- Abstract base class for integrators in variational formulations.
-
template<class Integrand>class InterfaceIntegral
- Represents expressions of the integral operator on the interface of a domain.
-
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>class InterfaceIntegral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>>
- Integration of the dot product of a trial and test operators.
-
template<class NestedDerived, class FES>class InterfaceIntegral<ShapeFunctionBase<NestedDerived, FES, TestSpace>>
- Integration of a test operator.
-
template<class T>struct IsTestFunction
- Type trait to check if a type is a TestFunction.
-
template<class FES>struct IsTestFunction<PETSc::Variational::TestFunction<FES>>
- Marks PETSc test functions as valid test functions.
-
template<class FES>struct IsTestFunction<TestFunction<FES>>
- Specialization for TestFunction types.
-
template<class T>struct IsTrialFunction
- Type trait to check if a type is a TrialFunction.
-
template<class Solution, class FES>struct IsTrialFunction<PETSc::Variational::TrialFunction<Solution, FES>>
- Marks PETSc trial functions as valid trial functions.
-
template<class Solution, class FES>struct IsTrialFunction<TrialFunction<Solution, FES>>
- Specialization for TrialFunction types.
-
template<class Operand>class Jacobian
- Represents the Jacobian matrix of a type.
-
template<size_t K, class Scalar, class Data, class Mesh>class Jacobian<GridFunction<H1<K, Math::Vector<Scalar>, Mesh>, Data>>
- Jacobian of an H1 vector GridFunction.
-
template<class Data, class Mesh, class Scalar>class Jacobian<GridFunction<P0g<Math::Vector<Scalar>, Mesh>, Data>>
- Jacobian of a P0g vector GridFunction (identically zero).
-
template<class Data, class Mesh, class Scalar>class Jacobian<GridFunction<P1<Math::Vector<Scalar>, Mesh>, Data>>
- Jacobian of a P1 vector GridFunction.
-
template<size_t K, class NestedDerived, class Scalar, class Mesh, ShapeFunctionSpaceType Space>class Jacobian<ShapeFunction<NestedDerived, H1<K, Math::Vector<Scalar>, Mesh>, Space>>
- Jacobian of an H1 ShapeFunction object.
-
template<class ShapeFunctionDerived, class Scalar, class Mesh, ShapeFunctionSpaceType Space>class Jacobian<ShapeFunction<ShapeFunctionDerived, P0g<Math::Vector<Scalar>, Mesh>, Space>>
- Jacobian of a P0g vector ShapeFunction (identically zero).
-
template<class ShapeFunctionDerived, class Range, class Mesh, ShapeFunctionSpaceType Space>class Jacobian<ShapeFunction<ShapeFunctionDerived, P1<Range, Mesh>, Space>>
- Jacobian of an P1 ShapeFunction object.
-
template<class Operand, class Derived>class JacobianBase
- Base class for Jacobian matrix operator implementations.
-
template<class FES, class Data, class Derived>class JacobianBase<GridFunction<FES, Data>, Derived>
- Jacobian of a P1 GridFunction.
-
template<size_t K>class JacobiPolynomial
- Jacobi polynomial evaluator.
-
template<class Range, class Context>class L2
- Arbitrary order broken Sobolev space.
-
template<size_t K>class LagrangeBasis1D
- Generic 1D Lagrange basis on arbitrary nodes.
-
template<size_t K>class LagrangeBasisPoint
- Lagrange basis on the reference point (0D element).
-
template<size_t K>class LagrangeBasisQuadrilateral
- Lagrange basis on the reference quadrilateral [0,1]².
-
template<size_t K>class LagrangeBasisSegment
- Lagrange basis on the reference segment [0,1] using GLL nodes.
-
template<size_t K>class LagrangeBasisTetrahedron
- Lagrange basis on the reference tetrahedron.
-
template<size_t K>class LagrangeBasisTriangle
- Lagrange basis on the reference triangle.
-
template<size_t K>class LagrangeBasisWedge
- Lagrange basis on the reference wedge (triangular prism).
-
template<size_t K>class LegendrePolynomial
- Compile-time degree Legendre polynomial evaluator.
-
template<class LHS, class RHS>class LEQ
-
template<class Solution, class FES>class LinearElasticityIntegral
- Helper class to construct linear elasticity integrators.
-
template<class Solution, class FES, class LambdaDerived, class MuDerived>class LinearElasticityIntegrator
- Bilinear form integrator for linear elasticity.
-
template<class Solution, class MuDerived, class LambdaDerived, class Range, class Mesh>class LinearElasticityIntegrator<Solution, P1<Range, Mesh>, MuDerived, LambdaDerived>
- P1 linear elasticity bilinear form integrator.
-
template<class FES, class VectorType>class LinearForm
- Represents a linear form on some finite element space.
-
template<class FES>class LinearForm<FES, ::Vec>
- Linear form specialization that assembles into a PETSc vector.
-
template<class FES>class LinearForm<FES, Math::Vector<typename FormLanguage::Traits<FES>::ScalarType>>
- Represents a linear form defined over some finite element space.
-
template<class Vector>class LinearFormBase
- Base class for linear form objects.
-
template<class Number>class LinearFormIntegratorBase
- Base class for linear form integrators.
-
template<class LHS, class RHS>class LT
-
template<class T>class MatrixFunction
-
template<class Scalar, class Derived>class MatrixFunctionBase
- Base class for matrix-valued functions defined on a mesh.
-
template<class LHSDerived, class RHSDerived>class Max<FunctionBase<LHSDerived>, FunctionBase<RHSDerived>>
- Represents the maximum operation between two functions.
-
template<class LHSDerived, class RHSDerived>class Min<FunctionBase<LHSDerived>, FunctionBase<RHSDerived>>
- Represents the minimum operation between two functions.
-
template<class LHS, class RHS>class Mult
- Represents the multiplication operation.
-
template<class LHSDerived, class RHSDerived>class Mult<FunctionBase<LHSDerived>, FunctionBase<RHSDerived>>
- Multiplication of two FunctionBase instances.
-
template<class LHSDerived, class RHSDerived, class FES, ShapeFunctionSpaceType Space>class Mult<FunctionBase<LHSDerived>, ShapeFunctionBase<RHSDerived, FES, Space>>
- Left Multiplication of a ShapeFunctionBase by a FunctionBase.
-
template<class LHSDerived, class RHSDerived, class FES, ShapeFunctionSpaceType Space>class Mult<ShapeFunctionBase<LHSDerived, FES, Space>, FunctionBase<RHSDerived>>
- Right multiplication of a ShapeFunctionBase by a FunctionBase.
-
template<class LHS, class RHS>class NEQ
-
template<class LHS, class RHS>class OR
-
template<class LHSDerived, class RHSDerived>class OR<BooleanFunctionBase<LHSDerived>, BooleanFunctionBase<RHSDerived>>
- Logical OR operator between two instances of BooleanFunctionBase.
-
template<class Range, class Mesh>class P0
- Degree 0 Lagrange finite element space.
-
template<>class P0<Real, Geometry::Mesh<Context::Local>>
- Scalar-valued piecewise constant (P0) Lagrange finite element space.
-
template<class Scalar>class P0Element
- Degree 0 Lagrange element.
-
template<class Scalar>class P0Element<Math::Vector<Scalar>>
- Piecewise constant (degree 0) vector Lagrange element.
-
template<class Range, class Mesh>class P0g
- Global constant (P0g) finite element space.
-
template<class Range, class Mesh>class P1
- Degree 1 Lagrange finite element space.
-
template<class Scalar>class P1<Math::Vector<Scalar>, Geometry::Mesh<Context::Local>>
- Vector-valued continuous piecewise linear Lagrange finite element space.
-
template<class Range>class P1<Range, Geometry::Mesh<Context::MPI>>
- Distributed P1 finite element space for MPI meshes.
-
template<class Scalar>class P1<Scalar, Geometry::Mesh<Context::Local>>
- Real valued Lagrange finite element space.
-
template<class Scalar>class P1Element
- Degree 1 Lagrange element.
-
template<class Scalar>class P1Element<Math::Vector<Scalar>>
- Continuous piecewise linear (degree 1) vector Lagrange element.
-
template<class Solution, class FES>class PeriodicBC<TrialFunction<Solution, FES>, IndexMap<IndexSet>>
- Represents a Peridodic boundary condition on a ShapeFunction object.
-
template<class Scalar>class PeriodicBCBase
- Abstract base class for a periodic boundary condition.
-
template<class LHS, class RHSDerived, class FES, ShapeFunctionSpaceType SpaceType>class Potential<LHS, ShapeFunctionBase<ShapeFunction<RHSDerived, FES, SpaceType>, FES, SpaceType>>
-
template<class Base, class Exponent>class Pow
- Represents the power function.
-
template<class BaseDerived, class Number>class Pow<FunctionBase<BaseDerived>, Number>
- Represents the power function.
-
template<class ... Parameters>class Problem
- Represents a variational problem.
-
template<class U, class V>class Problem<PETSc::Math::LinearSystem, U, V>
- PETSc variational problem for a single trial / test function pair.
-
template<class U1, class U2, class U3, class... Us>class Problem<PETSc::Math::LinearSystem, U1, U2, U3, Us...>
- PETSc variational problem for multiple coupled trial / test functions.
-
template<class LinearSystem>class ProblemBase
- Base class for variational problem objects.
-
template<class Scalar>class ProblemBodyBase
- Base class representing the body of a variational problem.
-
template<class LinearSystem, class U, class V>class ProblemUVBase
- General class to assemble linear systems with
OperatorandVectorgeneric types in a sequential context. -
template<class Integrand>class QuadratureRule
-
template<class Kernel, class Range, class Mesh, class LHSDerived, class RHSDerived>class QuadratureRule<Dot<Potential<Kernel, ShapeFunctionBase<LHSDerived, P1<Range, Mesh>, TrialSpace>>, ShapeFunctionBase<RHSDerived, P1<Range, Mesh>, TestSpace>>>
- Specialization for in the case of P1 trial/test shape functions.
-
template<size_t KTrial, size_t KTest, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<Div<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class LHSDerived, class RHSDerived, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<Div<ShapeFunction<LHSDerived, P1<Math::Vector<Real>, LHSMesh>, TrialSpace>>, P1<Math::Vector<Real>, LHSMesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<Real, RHSMesh>, TestSpace>, P1<Real, RHSMesh>, TestSpace>>>
- Specialization for in the case of P1 shape functions.
-
template<size_t KTrial, size_t KTest, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<Grad<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<Grad<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<Grad<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<Grad<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>>, P1<RHSRange, RHSMesh>, TestSpace>>>
- Integration of the isotropic Dot product of two instances of the P1 Grad of ShapeFunction.
-
template<size_t KTrial, size_t KTest, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<Jacobian<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<Jacobian<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<Jacobian<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<Jacobian<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>>, P1<RHSRange, RHSMesh>, TestSpace>>>
- Integration of the isotropic Frobenius inner product two instances of the P1 Jacobian of ShapeFunction.
-
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>class QuadratureRule<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>>
- Approximation of the integral of the the dot product between a trial shape function and a test shape function.
-
template<size_t KTrial, size_t KTest, class CoefficientDerived, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<FunctionBase<CoefficientDerived>, ShapeFunctionBase<Grad<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<Grad<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class CoefficientDerived, class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<FunctionBase<CoefficientDerived>, ShapeFunctionBase<Grad<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<Grad<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>>, P1<RHSRange, RHSMesh>, TestSpace>>>
- Integration of the anisotropic Dot product of two instances of the P1 Grad of ShapeFunction.
-
template<size_t KTrial, size_t KTest, class CoefficientDerived, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<FunctionBase<CoefficientDerived>, ShapeFunctionBase<Jacobian<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<Jacobian<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class CoefficientDerived, class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<FunctionBase<CoefficientDerived>, ShapeFunctionBase<Jacobian<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<Jacobian<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>>, P1<RHSRange, RHSMesh>, TestSpace>>>
- Integration of the anisotropic Frobenius inner product two instances of the P1 Jacobian of ShapeFunction.
-
template<size_t KTrial, size_t KTest, class CoefficientDerived, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<FunctionBase<CoefficientDerived>, ShapeFunctionBase<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class CoefficientDerived, class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<FunctionBase<CoefficientDerived>, ShapeFunctionBase<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>, P1<RHSRange, RHSMesh>, TestSpace>>>
- Integration of the anisotropic Dot product of two instances of the P1 ShapeFunction.
-
template<size_t KTrial, size_t KTest, class CoefficientDerived, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<ShapeFunctionBase<Jacobian<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, FunctionBase<CoefficientDerived>>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class CoefficientDerived, class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<Mult<ShapeFunctionBase<Jacobian<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>>, P1<LHSRange, LHSMesh>, TrialSpace>, FunctionBase<CoefficientDerived>>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>, P1<RHSRange, RHSMesh>, TestSpace>>>
- Specialization for ∫ (∇u · f)·v in the case of P1 shape functions.
-
template<size_t KTrial, size_t KTest, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<Div<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<size_t KTrial, size_t KTest, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Dot<ShapeFunctionBase<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>, H1<KTest, Scalar, Mesh>, TestSpace>>>
- Specialization for with H1 trial and test shape functions.
-
template<class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>, P1<RHSRange, RHSMesh>, TestSpace>>>
- Integration of the isotropic Dot product of two instances of the P1 ShapeFunction.
-
template<class LHSDerived, class RHSDerived, class LHSMesh, class RHSMesh>class QuadratureRule<Dot<ShapeFunctionBase<ShapeFunction<LHSDerived, P1<Real, LHSMesh>, TrialSpace>, P1<Real, LHSMesh>, TrialSpace>, ShapeFunctionBase<Div<ShapeFunction<RHSDerived, P1<Math::Vector<Real>, RHSMesh>, TestSpace>>, P1<Math::Vector<Real>, RHSMesh>, TestSpace>>>
- Specialization for in the case of P1 shape functions.
-
template<class FunctionDerived>class QuadratureRule<FunctionBase<FunctionDerived>>
- Quadrature rule for integrating functions on mesh polytopes.
-
template<class FES, class Data>class QuadratureRule<GridFunction<FES, Data>>
- Integration of a GridFunction object.
-
template<size_t KTrial, size_t KTest, class CoefficientDerived, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<Mult<FunctionBase<CoefficientDerived>, Dot<ShapeFunctionBase<ShapeFunction<LHSDerived, H1<KTrial, Scalar, Mesh>, TrialSpace>, H1<KTrial, Scalar, Mesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, H1<KTest, Scalar, Mesh>, TestSpace>, H1<KTest, Scalar, Mesh>, TestSpace>>>>
- Specialization for with H1 trial and test shape functions.
-
template<class CoefficientDerived, class LHSDerived, class RHSDerived, class LHSRange, class RHSRange, class LHSMesh, class RHSMesh>class QuadratureRule<Mult<FunctionBase<CoefficientDerived>, Dot<ShapeFunctionBase<ShapeFunction<LHSDerived, P1<LHSRange, LHSMesh>, TrialSpace>, P1<LHSRange, LHSMesh>, TrialSpace>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<RHSRange, RHSMesh>, TestSpace>, P1<RHSRange, RHSMesh>, TestSpace>>>>
- Specialization for in the case of P1 shape functions.
-
template<size_t K, class LHSDerived, class RHSDerived, class Scalar, class Mesh>class QuadratureRule<ShapeFunctionBase<Dot<FunctionBase<LHSDerived>, ShapeFunctionBase<ShapeFunction<RHSDerived, H1<K, Scalar, Mesh>, TestSpace>, H1<K, Scalar, Mesh>, TestSpace>>, H1<K, Scalar, Mesh>, TestSpace>>
- Specialization for with an H1 test shape function.
-
template<class LHSDerived, class RHSDerived, class Range, class Mesh>class QuadratureRule<ShapeFunctionBase<Dot<FunctionBase<LHSDerived>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<Range, Mesh>, TestSpace>, P1<Range, Mesh>, TestSpace>>, P1<Range, Mesh>, TestSpace>>
- Integration of the Dot product of some coefficient function and a P1 ShapeFunction.
-
template<class NestedDerived, class FES>class QuadratureRule<ShapeFunctionBase<NestedDerived, FES, TestSpace>>
- Approximation of the integral of a test shape function.
-
template<size_t K, class NestedDerived, class Scalar, class Mesh>class QuadratureRule<ShapeFunctionBase<ShapeFunction<NestedDerived, H1<K, Scalar, Mesh>, TestSpace>, H1<K, Scalar, Mesh>, TestSpace>>
- Specialization for with an H1 test shape function.
-
template<class NestedDerived, class Range, class Mesh>class QuadratureRule<ShapeFunctionBase<ShapeFunction<NestedDerived, P1<Range, Mesh>, TestSpace>, P1<Range, Mesh>, TestSpace>>
- Integration of a P1 ShapeFunction.
-
template<class Operand>class Re
- Extracts the real part of a complex-valued function.
-
template<class NestedDerived>class Re<FunctionBase<NestedDerived>>
- Specialization for FunctionBase operands.
-
template<class ... Values>class RealFunction
-
template<class F>class RealFunction<F>
- Represents a scalar function given by an arbitrary scalar function.
-
template<>class RealFunction<Real>
- Represents a constant scalar function with type Real.
-
template<class Derived>class RealFunctionBase
- Base class for real-valued functions defined on a mesh.
- class RelativeError
- Utility class for computing relative errors.
-
template<class Scalar, class Derived>class ScalarFunctionBase
- Base class for scalar-valued functions with templated scalar type.
-
template<class Derived, class FES, ShapeFunctionSpaceType Space>class ShapeFunction
-
template<class Derived, class FES = typename FormLanguage::class ShapeFunctionBase
Traits<Derived>:: FESType, ShapeFunctionSpaceType SpaceType = FormLanguage:: Traits<Derived>::SpaceType> - Base class for shape function objects.
-
template<class NestedDerived>class Sin<FunctionBase<NestedDerived>>
- Sine function operator for real-valued scalar functions.
-
template<class NestedDerived>class Sinh<FunctionBase<NestedDerived>>
- Hyperbolic sine function operator for real-valued scalar functions.
-
template<class TrialFES, class TestFES>class SparseProblem<TrialFES, TestFES, Math::Matrix<typename FormLanguage::Mult<typename FormLanguage::Traits<TrialFES>::ScalarType, typename FormLanguage::Traits<TestFES>::ScalarType>::Type>, Math::Vector<typename FormLanguage::Traits<TestFES>::ScalarType>>
- General class to assemble linear systems with
Math::andSparseMatrix Math::types in a serial context.Vector -
template<class NestedDerived>class Sqrt<FunctionBase<NestedDerived>>
- Square root operation for functions.
-
template<class LHS, class RHS>class Sum
- Represents the sum operation.
-
template<class LHSDerived, class RHSDerived>class Sum<FunctionBase<LHSDerived>, FunctionBase<RHSDerived>>
- Addition of two functions.
-
template<class LHSDerived, class RHSDerived, class FES, ShapeFunctionSpaceType Space>class Sum<ShapeFunctionBase<LHSDerived, FES, Space>, ShapeFunctionBase<RHSDerived, FES, Space>>
- Addition of two shape functions.
-
template<class Operand>class Tan
- Represents the tangent function.
-
template<class NestedDerived>class Tan<FunctionBase<NestedDerived>>
- Tangent function operator for real-valued scalar functions.
-
template<class FES>class TestFunction
- Represents a function which belongs to a test space.
-
template<size_t K>class TetrahedronBlend
- Optimized blending parameter α for tetrahedral warp-blend.
-
template<class Operand>class Trace
- Represents the trace of a matrix function.
-
template<class NestedDerived>class Trace<FunctionBase<NestedDerived>>
- Trace of a FunctionBase instance.
-
template<>class TraceOperator<FunctionBase>
- Trace (boundary restriction) operator for functions.
-
template<class Operand>class Transpose
- Represents the transpose matrix of some matrix .
-
template<class NestedDerived>class Transpose<FunctionBase<NestedDerived>>
- Transpose of a matrix-valued function.
-
template<class NestedDerived, class FES, ShapeFunctionSpaceType Space>class Transpose<ShapeFunctionBase<NestedDerived, FES, Space>>
- Transpose of a matrix-valued ShapeFunction.
-
template<class Solution, class FES>class TrialFunction
- Represents a function which belongs to a trial space.
-
template<class Solution, class FES>class TrialFunctionReference
- Reference wrapper for trial functions.
-
template<size_t K>class TriangleBlend
- Optimized blending parameter α for triangular warp-blend.
-
template<class Operand>class UnaryMinus
- Represent the negation of an operand.
-
template<class Number>class UnaryMinus<FormLanguage::List<LinearFormIntegratorBase<Number>>>
- Negation of a list of linear form integrators.
-
template<class Number>class UnaryMinus<FormLanguage::List<LocalBilinearFormIntegratorBase<Number>>>
- Negation of a list of bilinear form integrators.
-
template<class NestedDerived>class UnaryMinus<FunctionBase<NestedDerived>>
- Negation of a function.
-
template<class Number>class UnaryMinus<LinearFormIntegratorBase<Number>>
- Negation of a linear form integrator.
-
template<class Number>class UnaryMinus<LocalBilinearFormIntegratorBase<Number>>
- Negation of a bilinear form integrator.
-
template<class NestedDerived, class FES, ShapeFunctionSpaceType Space>class UnaryMinus<ShapeFunctionBase<NestedDerived, FES, Space>>
- Negation of a shape function.
-
template<size_t K>class VandermondeTriangle
- Vandermonde matrix for nodal-to-modal transformation on triangles.
-
template<class ... Values>class VectorFunction
-
template<class V, class ... Values>class VectorFunction<V, Values...>
- Represents a vector function which may be constructed from values which can be converted to objects of type RealFunction.
-
template<class Scalar, class Derived>class VectorFunctionBase
- Base class for vector-valued functions defined on a mesh.
-
template<size_t K>class WarpBlendTetrahedron
- Applies warp-blend algorithm to move tetrahedron nodes toward Fekete positions.
-
template<size_t K>class WarpBlendTriangle
- Applies warp-blend algorithm to move triangle nodes toward Fekete positions.
-
template<size_t K>class WarpFactor1D
- Computes the 1D warp factor for moving equispaced nodes to GLL positions.
-
template<size_t K>class WarpShiftFace2D
- Computes 2D warp-blend shift for a triangular face.
-
template<size_t K>class WarpShiftFace3D
- Computes 3D warp-blend shift for a tetrahedral face.
-
template<class Scalar>class Zero<Math::Vector<Scalar>>
- Vector zero function.
-
template<class Scalar>class Zero<Scalar>
- Scalar zero function.
Typedefs
-
template<class Scalar>using EssentialBoundary = FormLanguage::
List<DirichletBCBase<Scalar>> - Alias for a list of Dirichlet boundary conditions.
-
template<size_t K>using RealH1Element = H1Element<K, Real>
- Convenience alias for real-valued H1 element.
-
template<size_t K>using ComplexH1Element = H1Element<K, Complex>
- Convenience alias for complex-valued H1 element.
-
template<size_t K>using VectorH1Element = H1Element<K, Math::
Vector<Real>> - Convenience alias for real vector-valued H1 element.
-
template<size_t K>using ComplexVectorH1Element = H1Element<K, Math::
Vector<Complex>> - Convenience alias for complex vector-valued H1 element.
-
template<size_t K, class Mesh>using RealH1 = H1<K, Real, Mesh>
- Alias for a scalar real-valued H1 finite element space.
-
template<size_t K, class Mesh>using ComplexH1 = H1<K, Complex, Mesh>
- Alias for a scalar complex-valued H1 finite element space.
-
template<size_t K, class Mesh>using VectorH1 = H1<K, Math::
Vector<Real>, Mesh> - Alias for a vector-valued real H1 finite element space.
- using RealP0Element = P0Element<Real>
- Alias for P0Element<Real>
- using ComplexP0Element = P0Element<Complex>
- Alias for P0Element<Complex>
-
template<class Scalar>using VectorP0Element = P0Element<Math::
Vector<Scalar>> - Alias for P0Element<Math::Vector<Real>>
- using RealVectorP0Element = VectorP0Element<Real>
- Alias for real vector P0Element.
- using ComplexVectorP0Element = VectorP0Element<Complex>
- Alias for complex vector P0Element.
-
template<class Mesh>using RealP0 = P0<Real, Mesh>
- Alias for a scalar real-valued P0 finite element space.
-
template<class Mesh>using ComplexP0 = P0<Complex, Mesh>
- Alias for a scalar complex-valued P0 finite element space.
- using RealP1Element = P1Element<Real>
- Alias for P1Element<Real>
- using ComplexP1Element = P1Element<Complex>
- Alias for P1Element<Complex>
-
template<class Scalar>using VectorP1Element = P1Element<Math::
Vector<Scalar>> - Alias for P1Element<Math::
Vector<Scalar>> -
using RealVectorP1Element = P1Element<Math::
Vector<Real>> - Alias for real vector P1Element.
-
template<class Mesh>using RealP1 = P1<Real, Mesh>
- Alias for a scalar real-valued P1 finite element space.
-
template<class Mesh>using ComplexP1 = P1<Complex, Mesh>
- Alias for a scalar complex-valued P1 finite element space.
-
template<class Mesh>using VectorP1 = P1<Math::
Vector<Real>, Mesh> - Alias for a vector-valued real P1 finite element space.
-
template<class Scalar>using PeriodicBoundary = FormLanguage::
List<PeriodicBCBase<Scalar>> - Alias for a list of peridodic boundary conditions.
-
template<class Scalar>using VectorZero = Zero<Math::
Vector<Scalar>> - Convenience typedef for vector zero function.
Variables
-
template<class LHSDerived, class RHSDerived>AND< BooleanFunctionBase< LHSDerived >, BooleanFunctionBase< RHSDerived > >
- Deduction guide for AND.
- static const ShapeFunctionSpaceType TrialSpace constexpr
- Shorthand variable for ShapeFunctionSpaceType::
Trial. - static const ShapeFunctionSpaceType TestSpace constexpr
- Shorthand variable for ShapeFunctionSpaceType::
Test. -
template<class LHSDerived, class RHSDerived>GEQ< FunctionBase< LHSDerived >, FunctionBase< RHSDerived > >
- CTAD for GEQ.
-
template<class LHSDerived, class RHSDerived>GT< FunctionBase< LHSDerived >, FunctionBase< RHSDerived > >
- CTAD for GT.
-
template<class LHSDerived, class RHSDerived>LT< FunctionBase< LHSDerived >, FunctionBase< RHSDerived > >
- CTAD for LT.
Typedef documentation
#include <Rodin/Variational/DirichletBC.h>
template<class Scalar>
using Rodin:: Variational:: EssentialBoundary = FormLanguage:: List<DirichletBCBase<Scalar>>
Alias for a list of Dirichlet boundary conditions.
template<size_t K>
using Rodin:: Variational:: RealH1Element = H1Element<K, Real>
Convenience alias for real-valued H1 element.
| Template parameters | |
|---|---|
| K | Polynomial degree |
Equivalent to H1Element<K, Real>. Used for scalar-valued problems such as the Poisson equation or heat conduction.
template<size_t K>
using Rodin:: Variational:: ComplexH1Element = H1Element<K, Complex>
Convenience alias for complex-valued H1 element.
| Template parameters | |
|---|---|
| K | Polynomial degree |
Equivalent to H1Element<K, Complex>. Used for complex-valued problems such as wave equations in frequency domain.
template<size_t K>
using Rodin:: Variational:: VectorH1Element = H1Element<K, Math:: Vector<Real>>
Convenience alias for real vector-valued H1 element.
| Template parameters | |
|---|---|
| K | Polynomial degree |
Equivalent to H1Element<K, Math::. Used for vector-valued problems such as linear elasticity or fluid mechanics.
template<size_t K>
using Rodin:: Variational:: ComplexVectorH1Element = H1Element<K, Math:: Vector<Complex>>
Convenience alias for complex vector-valued H1 element.
| Template parameters | |
|---|---|
| K | Polynomial degree |
Equivalent to H1Element<K, Math::. Used for complex vector-valued problems such as electromagnetic wave propagation.
template<class Scalar>
using Rodin:: Variational:: VectorP1Element = P1Element<Math:: Vector<Scalar>>
Alias for P1Element<Math::
| Template parameters | |
|---|---|
| Scalar | Scalar type for vector components |
#include <Rodin/Variational/PeriodicBC.h>
template<class Scalar>
using Rodin:: Variational:: PeriodicBoundary = FormLanguage:: List<PeriodicBCBase<Scalar>>
Alias for a list of peridodic boundary conditions.
#include <Rodin/Variational/Zero.h>
template<class Scalar>
using Rodin:: Variational:: VectorZero = Zero<Math:: Vector<Scalar>>
Convenience typedef for vector zero function.
| Template parameters | |
|---|---|
| Scalar | Scalar type for vector components |
Variable documentation
#include <Rodin/Variational/AND.h>
template<class LHSDerived, class RHSDerived>
Rodin:: Variational:: AND< BooleanFunctionBase< LHSDerived >, BooleanFunctionBase< RHSDerived > >
Deduction guide for AND.
static const ShapeFunctionSpaceType Rodin:: Variational:: TrialSpace constexpr
#include <Rodin/Variational/ForwardDecls.h>
Shorthand variable for ShapeFunctionSpaceType::
static const ShapeFunctionSpaceType Rodin:: Variational:: TestSpace constexpr
#include <Rodin/Variational/ForwardDecls.h>
Shorthand variable for ShapeFunctionSpaceType::
#include <Rodin/Variational/GEQ.h>
template<class LHSDerived, class RHSDerived>
Rodin:: Variational:: GEQ< FunctionBase< LHSDerived >, FunctionBase< RHSDerived > >
CTAD for GEQ.
#include <Rodin/Variational/GT.h>
template<class LHSDerived, class RHSDerived>
Rodin:: Variational:: GT< FunctionBase< LHSDerived >, FunctionBase< RHSDerived > >
CTAD for GT.
#include <Rodin/Variational/LT.h>
template<class LHSDerived, class RHSDerived>
Rodin:: Variational:: LT< FunctionBase< LHSDerived >, FunctionBase< RHSDerived > >
CTAD for LT.