v0.15.0
Loading...
Searching...
No Matches
Classes | Typedefs | Enumerations | Functions
Forms Integrators

Classes and functions used to evaluate fields at integration pts, jacobians, etc.. More...

Collaboration diagram for Forms Integrators:

Classes

struct  MoFEM::EssentialBC< EleOp >
 Essential boundary conditions. More...
 
struct  MoFEM::EssentialBC< EleOp >::Assembly< A >
 Assembly methods. More...
 
struct  MoFEM::EssentialBC< EleOp >::Assembly< A >::LinearForm< I >
 
struct  MoFEM::EssentialBC< EleOp >::Assembly< A >::BiLinearForm< I >
 
struct  MoFEM::NaturalBC< EleOp >
 Natural boundary conditions. More...
 
struct  MoFEM::NaturalBC< EleOp >::Assembly< A >
 Assembly methods. More...
 
struct  MoFEM::NaturalBC< EleOp >::Assembly< A >::LinearForm< I >
 
struct  MoFEM::NaturalBC< EleOp >::Assembly< A >::BiLinearForm< I >
 
struct  MoFEM::OpSetBc
 Set indices on entities on finite element. More...
 
struct  MoFEM::OpUnSetBc
 Unset indices on entities on finite element. More...
 
struct  MoFEM::AssemblyTypeSelector< A >
 Template selector for assembly type. More...
 
struct  MoFEM::FormsIntegrators< EleOp >
 Integrator forms. More...
 
struct  MoFEM::FormsIntegrators< EleOp >::Assembly< A >
 Assembly methods. More...
 
struct  MoFEM::FormsIntegrators< EleOp >::Assembly< A >::LinearForm< I >
 Linear form. More...
 
struct  MoFEM::FormsIntegrators< EleOp >::Assembly< A >::BiLinearForm< I >
 Bi linear form. More...
 
struct  MoFEM::FormsIntegrators< EleOp >::Assembly< A >::TriLinearForm< I >
 Tri linear form. More...
 
struct  MoFEM::FormsIntegrators< EleOp >::Assembly< A >
 Bilinear integrator form. More...
 

Typedefs

template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradGrad = OpGradGradImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase >
 Gradient-gradient bilinear form: \((v_{,i},\beta(\mathbf{x}) u_{,j})_\Omega\).
 
template<int BASE_DIM, int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMass = OpMassImpl< BASE_DIM, FIELD_DIM, I, OpBase >
 Mass matrix bilinear form: \((v_i,\beta(\mathbf{x}) u_j)_\Omega\).
 
template<int BASE_DIM, int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMassCache = OpMassCacheImpl< BASE_DIM, FIELD_DIM, I, OpBase >
 Cached mass matrix bilinear form (same as OpMass with caching)
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradSymTensorGrad = OpGradSymTensorGradImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase >
 Symmetric tensor gradient bilinear form: \((v_k,D_{ijkl} u_{,l})_\Omega\).
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradGradSymTensorGradGrad = OpGradGradSymTensorGradGradImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase >
 Second-order gradient symmetric tensor bilinear form: \((v_{,ij},D_{ijkl} u_{,kl})_\Omega\).
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradTensorGrad = OpGradTensorGradImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase >
 General tensor gradient bilinear form: \((v_{,j},D_{ijkl} u_{,l})_\Omega\).
 
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixDivTimesScalar = OpMixDivTimesScalarImpl< SPACE_DIM, I, OpBase >
 Mixed divergence-scalar bilinear form: \((\lambda_{i,i},u)_\Omega\).
 
template<int SPACE_DIM, CoordinateTypes CoordSys = CARTESIAN>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixDivTimesVec = OpMixDivTimesVecImpl< SPACE_DIM, I, OpBase, CoordSys >
 Mixed tensor divergence-vector bilinear form: \((\lambda_{ij,j},u_{i})_\Omega\).
 
template<int SPACE_DIM, CoordinateTypes COORDINATE_SYSTEM = CARTESIAN>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixScalarTimesDiv = OpMixScalarTimesDivImpl< SPACE_DIM, I, OpBase, COORDINATE_SYSTEM >
 Mixed scalar-divergence bilinear form: \((\lambda,u_{i,i})_\Omega\).
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixVectorTimesGrad = OpMixVectorTimesGradImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase >
 Mixed vector-gradient bilinear form: \((\lambda_{i},u_{,j})_\Omega\).
 
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixTensorTimesGrad = OpMixTensorTimesGradImpl< SPACE_DIM, I, OpBase >
 Mixed tensor-gradient bilinear form: \((\lambda_{ij},u_{i,j})_\Omega\).
 
template<int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBrokenSpaceConstrain = OpBrokenSpaceConstrainImpl< FIELD_DIM, I, OpBrokenBase >
 Hybridization constraint bilinear form for broken spaces.
 
using MoFEM::ScalarFun = boost::function< double(const double, const double, const double)>
 Scalar function type.
 
template<int DIM>
using MoFEM::VectorFun = boost::function< FTensor::Tensor1< double, DIM >(const double, const double, const double)>
 Vector function type.
 
template<int BASE_DIM, int FIELD_DIM, typename S = SourceFunctionSpecialization>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpSource = OpSourceImpl< BASE_DIM, FIELD_DIM, I, typename S::template S< OpBase > >
 Integrate \((v,f(\mathbf{x}))_\Omega\) where f is a scalar function.
 
template<int BASE_DIM, int S = 1>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBaseTimesScalar = OpBaseTimesScalarImpl< BASE_DIM, S, I, OpBase >
 Scalar field integrator \((v_i,f)_\Omega\) where f is a scalar.
 
template<int BASE_DIM, int FIELD_DIM, int S>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBaseTimesVector = OpBaseTimesVectorImpl< BASE_DIM, FIELD_DIM, S, I, OpBase >
 Vector field integrator \((v,f_i)_\Omega\) where f is a vector.
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 1>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradTimesTensor = OpGradTimesTensorImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase >
 [Source operator]
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 1>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradTimesSymTensor = OpGradTimesSymTensorImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase >
 Integrate \((v_{,i},f_{ij})_\Omega\) where f is a symmetric tensor.
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, CoordinateTypes CoordSys = CARTESIAN>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixDivTimesU = OpMixDivTimesUImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase, CoordSys >
 Mixed formulation integrator: \((\lambda_{ij,j},u_{i})_\Omega\).
 
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixTensorTimesGradU = OpMixTensorTimesGradUImpl< SPACE_DIM, I, OpBase >
 Mixed formulation integrator: \((\lambda_{ij},u_{i,j})_\Omega\).
 
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixVecTimesDivLambda = OpMixVecTimesDivLambdaImpl< SPACE_DIM, I, OpBase >
 Mixed formulation integrator: \((u_{i},\lambda_{ij,j})_\Omega\).
 
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpNormalMixVecTimesScalar = OpNormalMixVecTimesScalarImpl< SPACE_DIM, I, OpBase >
 Boundary integrator: normal vector times scalar function.
 
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpNormalMixVecTimesVectorField = OpNormalMixVecTimesVectorFieldImpl< SPACE_DIM, I, OpBase >
 Boundary integrator: normal vector times vector field.
 
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpConvectiveTermRhs = OpConvectiveTermRhsImpl< BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase >
 Convective term integrator: \((v, u_i \mathbf{y}_{,i})\).
 
template<int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBrokenSpaceConstrainDHybrid = OpBrokenSpaceConstrainDHybridImpl< FIELD_DIM, I, OpBase >
 Hybridization constraint integrator for broken spaces.
 

Enumerations

enum  MoFEM::AssemblyType {
  MoFEM::PETSC , MoFEM::SCHUR , MoFEM::BLOCK_MAT , MoFEM::BLOCK_SCHUR ,
  MoFEM::BLOCK_PRECONDITIONER_SCHUR , MoFEM::USER_ASSEMBLE , MoFEM::LAST_ASSEMBLE
}
 [Storage and set boundary conditions] More...
 
enum  MoFEM::IntegrationType { MoFEM::GAUSS , MoFEM::USER_INTEGRATION , MoFEM::LAST_INTEGRATION }
 Form integrator integration types. More...
 

Functions

template<>
MoFEMErrorCode MoFEM::VecSetValues< EssentialBcStorage > (Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
 Set values to vector in operator.
 
double MoFEM::scalar_fun_one (const double, const double, const double)
 Default scalar function returning 1.
 

Detailed Description

Classes and functions used to evaluate fields at integration pts, jacobians, etc..

Typedef Documentation

◆ OpBaseTimesScalar

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int S = 1>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBaseTimesScalar = OpBaseTimesScalarImpl<BASE_DIM, S, I, OpBase>

Scalar field integrator \((v_i,f)_\Omega\) where f is a scalar.

This operator integrates the product of vector test functions with a scalar field value. Each component of the test function is multiplied by the same scalar value. Useful for:

  • Uniform pressure loads in structural mechanics
  • Constant temperature fields in thermal problems
  • Scalar material properties affecting vector fields
Note
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
SStride or scaling factor (default: 1)

Definition at line 74 of file LinearFormsIntegrators.hpp.

◆ OpBaseTimesVector

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int S>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBaseTimesVector = OpBaseTimesVectorImpl<BASE_DIM, FIELD_DIM, S, I, OpBase>

Vector field integrator \((v,f_i)_\Omega\) where f is a vector.

This operator integrates the product of test functions with vector field values. Each test function component is multiplied by the corresponding component of the vector field. Applications include:

  • Vector body forces (wind loads, electromagnetic forces)
  • Distributed loads with directional components
  • Multi-component source terms in coupled problems
Note
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the vector field
SStride or scaling factor

Definition at line 104 of file LinearFormsIntegrators.hpp.

◆ OpBrokenSpaceConstrain

template<typename EleOp >
template<typename EleOp >
template<int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBrokenSpaceConstrain = OpBrokenSpaceConstrainImpl<FIELD_DIM, I, OpBrokenBase>

Hybridization constraint bilinear form for broken spaces.

This operator assembles constraint matrices during hybridization of broken finite element spaces. Hybridization introduces Lagrange multipliers on element interfaces to enforce continuity in discontinuous Galerkin methods. The operator handles the coupling between interior element solutions and interface unknowns. Applications include:

  • Discontinuous Galerkin (DG) methods
  • Hybridizable discontinuous Galerkin (HDG) methods
  • Mixed finite element methods with broken spaces
  • Static condensation procedures
Note
Used specifically for broken space hybridization procedures
Handles interface coupling in DG methods
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
FIELD_DIMDimension of the field

Definition at line 313 of file BiLinearFormsIntegrators.hpp.

◆ OpBrokenSpaceConstrainDHybrid

template<typename EleOp >
template<typename EleOp >
template<int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpBrokenSpaceConstrainDHybrid = OpBrokenSpaceConstrainDHybridImpl<FIELD_DIM, I, OpBase>

Hybridization constraint integrator for broken spaces.

This operator assembles the right-hand side for constraint matrix C during hybridization of broken finite element spaces. Hybridization introduces Lagrange multipliers to enforce continuity across element boundaries in discontinuous Galerkin methods. Applications include:

  • Discontinuous Galerkin methods
  • Mixed finite element methods with broken spaces
  • Hybridizable discontinuous Galerkin (HDG) methods
Note
Used specifically for broken space hybridization
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
FIELD_DIMDimension of the field

Definition at line 303 of file LinearFormsIntegrators.hpp.

◆ OpConvectiveTermRhs

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpConvectiveTermRhs = OpConvectiveTermRhsImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase>

Convective term integrator: \((v, u_i \mathbf{y}_{,i})\).

This operator integrates convective terms where a velocity field \(u_i\) transports a scalar or vector field \(\mathbf{y}\). The convective operator \(u_i \mathbf{y}_{,i}\) represents advection. Applications include:

  • Advection in fluid mechanics
  • Convective transport in heat transfer
  • Material transport in multi-physics problems
Note
\(u_i\) is the velocity field at integration points
\(\mathbf{y}_{,i}\) is the gradient of the transported field
\(\mathbf{y}\) can be scalar or vector field
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem

Definition at line 282 of file LinearFormsIntegrators.hpp.

◆ OpGradGrad

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradGrad = OpGradGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase>

Gradient-gradient bilinear form: \((v_{,i},\beta(\mathbf{x}) u_{,j})_\Omega\).

This operator integrates the product of test function gradients with trial function gradients, scaled by a coefficient \(\beta(\mathbf{x})\). This is fundamental for diffusion, heat conduction, and elasticity problems. Applications include:

  • Laplacian operators in diffusion equations
  • Stiffness matrices in structural mechanics
  • Gradient-based regularization terms
Note
\(\beta(\mathbf{x})\) is a spatially varying coefficient
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem

Definition at line 67 of file BiLinearFormsIntegrators.hpp.

◆ OpGradGradSymTensorGradGrad

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradGradSymTensorGradGrad = OpGradGradSymTensorGradGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>

Second-order gradient symmetric tensor bilinear form: \((v_{,ij},D_{ijkl} u_{,kl})_\Omega\).

This operator integrates second-order gradients (Hessians) of test and trial functions through a symmetric fourth-order tensor \(D_{ijkl}\). Used in higher-order partial differential equations and plate/shell theories. Applications include:

  • Biharmonic operators in plate bending problems
  • Fourth-order diffusion equations
  • Higher-order regularization in optimization
  • Gradient damage models in mechanics
Note
\(D_{ijkl}\) is a tensor with no symmetries
Requires C1 continuous finite elements
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem
SAdditional parameter for specialization (default: 0)

Definition at line 152 of file BiLinearFormsIntegrators.hpp.

◆ OpGradSymTensorGrad

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradSymTensorGrad = OpGradSymTensorGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>

Symmetric tensor gradient bilinear form: \((v_k,D_{ijkl} u_{,l})_\Omega\).

This operator integrates test functions with trial function gradients through a symmetric fourth-order tensor \(D_{ijkl}\) with minor and major symmetry. Essential for anisotropic diffusion and linear elasticity. Applications include:

  • Elasticity stiffness matrices with material symmetry
  • Anisotropic diffusion operators
  • Constitutive relations in continuum mechanics
  • Stress-strain coupling in solid mechanics
Note
\(D_{ijkl}\) is a tensor with minor & major symmetry
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem
SAdditional parameter for specialization (default: 0)

Definition at line 127 of file BiLinearFormsIntegrators.hpp.

◆ OpGradTensorGrad

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradTensorGrad = OpGradTensorGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>

General tensor gradient bilinear form: \((v_{,j},D_{ijkl} u_{,l})_\Omega\).

This operator integrates test function gradients with trial function gradients through a general fourth-order tensor \(D_{ijkl}\) with no symmetry assumptions. This is the most general form for anisotropic material behavior. Applications include:

  • General anisotropic elasticity without symmetry
  • Non-symmetric material tensors in advanced materials
  • Coupled multi-physics problems with complex constitutive relations
  • Plasticity with kinematic hardening
Note
\(D_{ijkl}\) is a tensor with no symmetries
More computationally expensive than symmetric versions
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem
SAdditional parameter for specialization (default: 0)

Definition at line 178 of file BiLinearFormsIntegrators.hpp.

◆ OpGradTimesSymTensor

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 1>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradTimesSymTensor = OpGradTimesSymTensorImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>

Integrate \((v_{,i},f_{ij})_\Omega\) where f is a symmetric tensor.

This operator is specialized for symmetric tensors, providing computational efficiency when the tensor field has symmetry properties. Applications include:

  • Symmetric stress tensors in solid mechanics
  • Symmetric diffusion tensors in transport problems
  • Symmetric material property tensors
Note
\(f_{ij}\) is a symmetric tensor field at integration points
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem
SStride or scaling factor (default: 1)

Definition at line 152 of file LinearFormsIntegrators.hpp.

◆ OpGradTimesTensor

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 1>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpGradTimesTensor = OpGradTimesTensorImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>

[Source operator]

[Grad times tensor]

Integrate \((v_{,i},f_{ij})_\Omega\) where f is a tensor

This operator integrates the product of test function gradients with tensor field values. The gradient of test functions contracts with the tensor field. Common applications include:

  • Stress-based formulations in mechanics
  • Flux calculations in transport problems
  • Gradient-dependent source terms
Note
\(f_{ij}\) is a tensor field at integration points
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem
SStride or scaling factor (default: 1)

Definition at line 130 of file LinearFormsIntegrators.hpp.

◆ OpMass

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMass = OpMassImpl<BASE_DIM, FIELD_DIM, I, OpBase>

Mass matrix bilinear form: \((v_i,\beta(\mathbf{x}) u_j)_\Omega\).

This operator integrates the product of test and trial functions, scaled by a coefficient \(\beta(\mathbf{x})\). This forms the mass matrix in dynamic problems and identity-type operators in various formulations. Applications include:

  • Mass matrices in dynamic finite element analysis
  • Identity operators in parabolic equations
  • Penalty terms and stabilization methods
  • Time-dependent term discretization
Note
FIELD_DIM = 4 or 9 assumes OpMass is for tensorial field 2x2 or 3x3
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Todo:
Consider another template parameter RANK to distinguish between scalar, vectorial and tensorial fields
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field

Definition at line 91 of file BiLinearFormsIntegrators.hpp.

◆ OpMassCache

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMassCache = OpMassCacheImpl<BASE_DIM, FIELD_DIM, I, OpBase>

Cached mass matrix bilinear form (same as OpMass with caching)

This is a cached version of OpMass that stores intermediate computations to improve performance in repeated evaluations. Useful when the same mass matrix structure is evaluated multiple times.

This operator integrates the product of test and trial functions, scaled by a coefficient \(\beta(\mathbf{x})\). This forms the mass matrix in dynamic problems and identity-type operators in various formulations. Applications include:

  • Mass matrices in dynamic finite element analysis
  • Identity operators in parabolic equations
  • Penalty terms and stabilization methods
  • Time-dependent term discretization
Note
FIELD_DIM = 4 or 9 assumes OpMass is for tensorial field 2x2 or 3x3
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Todo:
Consider another template parameter RANK to distinguish between scalar, vectorial and tensorial fields
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field

Definition at line 104 of file BiLinearFormsIntegrators.hpp.

◆ OpMixDivTimesScalar

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixDivTimesScalar = OpMixDivTimesScalarImpl<SPACE_DIM, I, OpBase>

Mixed divergence-scalar bilinear form: \((\lambda_{i,i},u)_\Omega\).

This operator integrates the divergence of vector test functions with scalar trial functions. Essential for mixed formulations enforcing divergence constraints such as incompressibility. Applications include:

  • Pressure-velocity coupling in incompressible flow
  • Mixed formulations for Darcy flow
  • Constraint enforcement in optimization problems
  • Divergence-free condition implementation
Note
\(\lambda_{i,i}\) is the divergence of vector test function
\(u\) is a scalar trial function
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem

Definition at line 200 of file BiLinearFormsIntegrators.hpp.

◆ OpMixDivTimesU

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM, CoordinateTypes CoordSys = CARTESIAN>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixDivTimesU = OpMixDivTimesUImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase, CoordSys>

Mixed formulation integrator: \((\lambda_{ij,j},u_{i})_\Omega\).

This operator integrates the divergence of a tensor field \(\lambda_{ij}\) with vector test functions \(u_i\). Essential for mixed formulations where stress equilibrium or vector field divergence is enforced. Applications include:

  • Stress-displacement coupling in solid mechanics
  • Mixed formulations with tensor constraints
  • Multi-physics problems with vector-tensor coupling
Note
\(\lambda_{ij,j}\) is the divergence of tensor field at integration points
\(u_i\) is the vector test function
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem
CoordSysCoordinate system type (default: CARTESIAN)

Definition at line 177 of file LinearFormsIntegrators.hpp.

◆ OpMixDivTimesVec

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM, CoordinateTypes CoordSys = CARTESIAN>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixDivTimesVec = OpMixDivTimesVecImpl<SPACE_DIM, I, OpBase, CoordSys>

Mixed tensor divergence-vector bilinear form: \((\lambda_{ij,j},u_{i})_\Omega\).

This operator integrates the divergence of tensor test functions with vector trial functions. Used in mixed formulations with stress-displacement coupling and equilibrium enforcement. Applications include:

  • Stress-displacement mixed formulations in elasticity
  • Mixed finite element methods for mechanics
  • Equilibrium constraint enforcement
  • Multi-physics coupling with tensor fields
Note
\(\lambda_{ij,j}\) is the divergence of tensor test function
\(u_i\) is a vector trial function
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem
CoordSysCoordinate system type (default: CARTESIAN)

Definition at line 222 of file BiLinearFormsIntegrators.hpp.

◆ OpMixScalarTimesDiv

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM, CoordinateTypes COORDINATE_SYSTEM = CARTESIAN>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixScalarTimesDiv = OpMixScalarTimesDivImpl<SPACE_DIM, I, OpBase, COORDINATE_SYSTEM>

Mixed scalar-divergence bilinear form: \((\lambda,u_{i,i})_\Omega\).

This operator integrates scalar test functions with the divergence of vector trial functions. This is the transpose of OpMixDivTimesScalar and enforces the same constraints but with different test/trial function roles. Applications include:

  • Symmetric mixed formulations
  • Lagrange multiplier methods for divergence constraints
  • Saddle point problems in fluid mechanics
  • Mixed Poisson formulations
Note
\(\lambda\) is a scalar test function
\(u_{i,i}\) is the divergence of vector trial function
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem
COORDINATE_SYSTEMCoordinate system type (default: CARTESIAN)

Definition at line 244 of file BiLinearFormsIntegrators.hpp.

◆ OpMixTensorTimesGrad

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixTensorTimesGrad = OpMixTensorTimesGradImpl<SPACE_DIM, I, OpBase>

Mixed tensor-gradient bilinear form: \((\lambda_{ij},u_{i,j})_\Omega\).

This operator integrates tensor test functions with the gradient of vector trial functions. Essential for mixed formulations in solid mechanics where stress tensors couple with displacement gradients. Applications include:

  • Stress-strain mixed formulations
  • Mixed finite element methods in elasticity
  • Constitutive relation enforcement
  • Coupled stress-displacement problems
Note
\(\lambda_{ij}\) is a tensor test function
\(u_{i,j}\) is the gradient of vector trial function
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem

Definition at line 290 of file BiLinearFormsIntegrators.hpp.

◆ OpMixTensorTimesGradU

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixTensorTimesGradU = OpMixTensorTimesGradUImpl<SPACE_DIM, I, OpBase>

Mixed formulation integrator: \((\lambda_{ij},u_{i,j})_\Omega\).

This operator integrates tensor field values with the gradient of vector test functions. The tensor field contracts with the gradient of test functions. Applications include:

  • Stress-strain coupling in mechanics
  • Mixed formulations with gradient constraints
  • Tensor-gradient field interactions
Note
\(\lambda_{ij}\) is a tensor field at integration points
\(u_{i,j}\) is the gradient of vector test function
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem

Definition at line 198 of file LinearFormsIntegrators.hpp.

◆ OpMixVecTimesDivLambda

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixVecTimesDivLambda = OpMixVecTimesDivLambdaImpl<SPACE_DIM, I, OpBase>

Mixed formulation integrator: \((u_{i},\lambda_{ij,j})_\Omega\).

This operator integrates vector test functions with the divergence of a tensor field. This is the adjoint form of OpMixDivTimesU, commonly used in mixed formulations for equilibrium enforcement. Applications include:

  • Adjoint stress equilibrium in mechanics
  • Dual formulations with tensor divergence
  • Vector-tensor coupling in mixed problems
Note
\(u_i\) is a vector test function
\(\lambda_{ij,j}\) is the divergence of tensor field at integration points
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem

Definition at line 218 of file LinearFormsIntegrators.hpp.

◆ OpMixVectorTimesGrad

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpMixVectorTimesGrad = OpMixVectorTimesGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase>

Mixed vector-gradient bilinear form: \((\lambda_{i},u_{,j})_\Omega\).

This operator integrates vector test functions with the gradient of trial functions. Used in mixed formulations where vector fields couple with gradient terms. Applications include:

  • Mixed formulations with vector Lagrange multipliers
  • Gradient constraint enforcement
  • Vector-scalar coupling in multi-physics problems
  • Stabilization terms in finite element methods
Note
\(\lambda_i\) is a vector test function
\(u_{,j}\) is the gradient of trial function
Implementation details are in BiLinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field
SPACE_DIMSpatial dimension of the problem

Definition at line 268 of file BiLinearFormsIntegrators.hpp.

◆ OpNormalMixVecTimesScalar

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpNormalMixVecTimesScalar = OpNormalMixVecTimesScalarImpl<SPACE_DIM, I, OpBase>

Boundary integrator: normal vector times scalar function.

This operator multiplies vector test functions by the normal vector on a face and integrates with a scalar field. Essential for natural boundary conditions in mixed formulations. Applications include:

  • Natural boundary conditions in fluid mechanics
  • Surface flux conditions in transport problems
  • Mixed boundary value problems
Note
This operator is used on domain boundaries (faces)
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem

Definition at line 238 of file LinearFormsIntegrators.hpp.

◆ OpNormalMixVecTimesVectorField

template<typename EleOp >
template<typename EleOp >
template<int SPACE_DIM>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpNormalMixVecTimesVectorField = OpNormalMixVecTimesVectorFieldImpl<SPACE_DIM, I, OpBase>

Boundary integrator: normal vector times vector field.

This operator multiplies vector test functions by the normal vector on a face and integrates with a vector field. Used for natural boundary conditions involving vector quantities. Applications include:

  • Vector-valued boundary conditions in mechanics
  • Surface traction conditions in solid mechanics
  • Vector flux conditions in multi-physics problems
Note
This operator is used on domain boundaries (faces)
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
SPACE_DIMSpatial dimension of the problem

Definition at line 258 of file LinearFormsIntegrators.hpp.

◆ OpSource

template<typename EleOp >
template<typename EleOp >
template<int BASE_DIM, int FIELD_DIM, typename S = SourceFunctionSpecialization>
using MoFEM::FormsIntegrators< EleOp >::Assembly< A >::OpSource = OpSourceImpl<BASE_DIM, FIELD_DIM, I, typename S::template S<OpBase> >

Integrate \((v,f(\mathbf{x}))_\Omega\) where f is a scalar function.

This operator integrates the product of test functions with a scalar source function defined at integration points. Commonly used for:

  • Body forces in mechanics (gravitational, thermal expansion)
  • Source terms in heat transfer and diffusion problems
  • Prescribed loads and boundary conditions
Note
\(f(\mathbf{x})\) is a scalar function of coordinates
Implementation details are in LinearFormsIntegratorsImpl.hpp
Template Parameters
BASE_DIMDimension of the basis functions (1D, 2D, 3D)
FIELD_DIMDimension of the field being integrated
SSource function specialization type (default: SourceFunctionSpecialization)

Definition at line 54 of file LinearFormsIntegrators.hpp.

◆ ScalarFun

using MoFEM::ScalarFun = typedef boost::function<double(const double, const double, const double)>

#include <src/finite_elements/FormsIntegrators.hpp>

Scalar function type.

Definition at line 248 of file FormsIntegrators.hpp.

◆ VectorFun

template<int DIM>
using MoFEM::VectorFun = typedef boost::function<FTensor::Tensor1<double, DIM>( const double, const double, const double)>

#include <src/finite_elements/FormsIntegrators.hpp>

Vector function type.

Template Parameters
DIMdimension of the return

Definition at line 289 of file FormsIntegrators.hpp.

Enumeration Type Documentation

◆ AssemblyType

#include <src/finite_elements/FormsIntegrators.hpp>

[Storage and set boundary conditions]

Form integrator assembly types

Enumerator
PETSC 

Standard PETSc assembly.

SCHUR 

Schur complement assembly.

BLOCK_MAT 

Block matrix assembly.

BLOCK_SCHUR 

Block Schur assembly.

BLOCK_PRECONDITIONER_SCHUR 

Block preconditioner Schur assembly.

USER_ASSEMBLE 

User-defined assembly.

LAST_ASSEMBLE 

Sentinel value.

Definition at line 200 of file FormsIntegrators.hpp.

200 {
201 PETSC, ///< Standard PETSc assembly
202 SCHUR, ///< Schur complement assembly
203 BLOCK_MAT, ///< Block matrix assembly
204 BLOCK_SCHUR, ///< Block Schur assembly
205 BLOCK_PRECONDITIONER_SCHUR, ///< Block preconditioner Schur assembly
206 USER_ASSEMBLE, ///< User-defined assembly
207 LAST_ASSEMBLE ///< Sentinel value
208};
@ PETSC
Standard PETSc assembly.
@ BLOCK_PRECONDITIONER_SCHUR
Block preconditioner Schur assembly.
@ LAST_ASSEMBLE
Sentinel value.
@ BLOCK_MAT
Block matrix assembly.
@ USER_ASSEMBLE
User-defined assembly.
@ SCHUR
Schur complement assembly.
@ BLOCK_SCHUR
Block Schur assembly.

◆ IntegrationType

#include <src/finite_elements/FormsIntegrators.hpp>

Form integrator integration types.

Enumerator
GAUSS 

Gaussian quadrature integration.

USER_INTEGRATION 

User-defined integration rules.

LAST_INTEGRATION 

Sentinel value.

Definition at line 237 of file FormsIntegrators.hpp.

237 {
238 GAUSS, ///< Gaussian quadrature integration
239 USER_INTEGRATION, ///< User-defined integration rules
240 LAST_INTEGRATION ///< Sentinel value
241};
@ USER_INTEGRATION
User-defined integration rules.
@ LAST_INTEGRATION
Sentinel value.
@ GAUSS
Gaussian quadrature integration.

Function Documentation

◆ scalar_fun_one()

double MoFEM::scalar_fun_one ( const double  ,
const double  ,
const double   
)
inline

#include <src/finite_elements/FormsIntegrators.hpp>

Default scalar function returning 1.

Parameters
xcoordinate
ycoordinate
zcoordinate
Returns
1.0

Definition at line 260 of file FormsIntegrators.hpp.

260 {
261 return 1;
262}

◆ VecSetValues< EssentialBcStorage >()

template<>
MoFEMErrorCode MoFEM::VecSetValues< EssentialBcStorage > ( Vec  V,
const EntitiesFieldData::EntData data,
const double ptr,
InsertMode  iora 
)

#include <src/finite_elements/FormsIntegrators.hpp>

Set values to vector in operator.

MoFEM::FieldEntity provides MoFEM::FieldEntity::getWeakStoragePtr() storage function which allows to transfer data between FEs or operators processing the same entities.

When MoFEM::OpSetBc is pushed in weak storage indices taking in account indices which are skip to take boundary conditions are stored. Those entities are used by VecSetValues.

Parameters
V
data
ptr
iora
Returns
MoFEMErrorCode
Parameters
V
data
ptr
iora
Returns
MoFEMErrorCode
Examples
mofem/atom_tests/test_cache_on_entities.cpp.

Definition at line 75 of file FormsIntegrators.cpp.

77 {
78
79#ifndef NDEBUG
80 if (!V)
81 CHK_THROW_MESSAGE(MOFEM_INVALID_DATA, "Pointer to PETSc vector is null in "
82 "VecSetValues<EssentialBcStorage>");
83#endif
84
85 if (!data.getFieldEntities().empty()) {
86 if (auto e_ptr = data.getFieldEntities()[0]) {
87 if (auto stored_data_ptr =
88 e_ptr->getSharedStoragePtr<EssentialBcStorage>()) {
90 VecSetOption(V, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE),
91 "In VecSetValues<EssentialBcStorage> failed to set option "
92 "VEC_IGNORE_NEGATIVE_INDICES");
93 return ::VecSetValues(V, stored_data_ptr->entityIndices.size(),
94 &*stored_data_ptr->entityIndices.begin(), ptr,
95 iora);
96 }
97 }
98 }
99 // If no storage is found, use the indices from data
100 return ::VecSetValues(V, data.getIndices().size(),
101 &*data.getIndices().begin(), ptr, iora);
102}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ MOFEM_INVALID_DATA
Definition definitions.h:36