v0.14.0
User data operators table
Note
Work in progress. Tables not complete.

Template arguments:

  • BASE_DIM is dimension of base function. Typically bases for L2 and H1 spaces are scalar bases. Thus BASE_DIM is 1. Typically H-div and H-curl bases are vectorial bases and BASE_DIM 3. For some special case could be exceptions.
  • FIELD_DIM is dimension of the field. We can have scalar field, vectorial field. For example, the scalar field for scalar base function has one coefficient for base function and FIELD_DIM = 1. Vectorial field with scalar base in 3D has three coefficients for each base function and FIELD_DIM = 3. On the other hand, a vectorial field with a vectorial base has one coefficient for each base function and FIELD_DIM = 3, and tensorial field of rank two in 3D on vectorial base has three coefficients for each base function and FIELD_DIM = 9.
  • SPACE_DIM is dismension of space. Typically, it is used to indicate how many elements have a gradient of base function or field.
  • ASSEMBLE_TYPE by this template assembly method is indicated. For example, if the assembly is to CST compressed matrix assembly type is PETSc. For some other cases, for example, integration to nested matrices, or other matrices formats could be added and implemented.
  • INTEGRATION_TYPE indicate method used to assembly elements. If standard numerical integration is used, then the integration type is GAUSS. Some other more sophisticated integration methods can be implemented to exploit some base functions for fast integration.
  • OP type of opetaor for voulme, face, edge or vertex finite element.

Field operators

Table of user data operators
Operator Usage
Scalar fields

Evaluate scalar field at integration points. MoFEM::OpCalculateScalarFieldValues

\[ u^h(\pmb\xi^g) = \sum_i^N \phi^i(\pmb\xi^g) \overline{u}_i \]

where \(u^h(\pmb\xi^g)\) is the field value at integration point \(g\), \(\phi^i(\pmb\xi^g)\) is value of base function \(i\) at integtation point coordinate \(\pmb\xi^g\), and \(\overline{u}_i\) is \(i\)-th degree of freedom.

This function has a similar variant that is used to calculate time derivative, when a time solver is used, see MoFEM::OpCalculateScalarFieldValuesDot

\[ \frac{\partial u^h}{\partial t} (\pmb\xi^g) = \sum_i^N \phi^i(\pmb\xi^g) \frac{\partial \overline{u}}{\partial t} \]

where \( \frac{\partial u^h}{\partial t} (\pmb\xi^g) \) (or \(\dot{u}^h(\pmb\xi^g)\) in short) is the rate value at integration point \(g\) and \( \frac{\partial \overline{u}}{\partial t} \) (or \(\dot{\overline{u}}_i\) in short) is the rate of \(i\)-th degree of freedom. Note that this funtion will work only with TS solver which is the numerical scheme used to calculate time derivative depends on time solver method.

Basic usage:

auto vals_U_at_integration_pts = boost::make_shared<VectorDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateScalarFieldValues("U", vals_U_at_integration_pts));

Operator for a given field calulates values of a scalar field at integration points. Values are stored in, for example, vals_U_at_integration_pts. Then value can be accessed as follows

auto nb_integration_pts = vals_U_at_integration_pts->size();
auto t_val = getFTensor0FromVec(vals_U_at_integration_pts);
for(auto gg = 0;gg!=nb_integration_pts;++gg) {
// do something here
++t_val;
}

Evaluate scalar field gradient

MoFEM::OpCalculateScalarFieldGradient

\[ {u}^h_{,j} (\pmb\xi^g) = \sum_i^N \phi^i_{,j}(\pmb\xi^g) \overline{u}_i \]

where

\[ \phi^i_{,j} = \left. \frac{\partial \phi^i}{\partial \xi_j} \right|_{\pmb\xi^g} \]

Depending on code context, derivatives are calculated in reference element or on the current physical element. For most cases, the gradient is calculated on the current physical element. Then derivative of base function is interpreted as

\[ \phi^i_{,j} = \left. \frac{\partial \phi^i}{\partial \xi_ks} \frac{\partial \xi_k}{\partial x_j} \right|_{\mathbf{x}(\pmb\xi^g)} \]

Basic usage:

constexpr int SPACE_DIM = 3;
auto grad_U_at_integration_pts = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateVectorFieldValues<SPACE_DIM>("U", grad_U_at_integration_pts));

Operator for a given field calculate values of scalar field gradient at integration points. Values are stored in, for example, grad_U_at_integration_pts. Then value can be accessed as follows

constexpr int FIELD_DIM = 3;
auto nb_integration_pts = grad_U_at_integration_pts->size2();
auto t_grad = getFTensor1FromMat<FIELD_DIM>((grad_U_at_integration_pts);
for(auto gg = 0;gg!=nb_integration_pts;++gg) {
auto dot = t_grad(i)*t_grad(i);
++t_grad;
}

Vector fields

Evaluate vector field at integration points MoFEM::OpCalculateVectorFieldValues

\[ u_j^h(\pmb\xi^g) = \sum_i^N \phi^i(\pmb\xi^g) \overline{u}_i^j \]

where \(j\) is vector element index.

This function has a similar variant to calculate rates when a time solver is used; see MoFEM::OpCalculateVectorFieldValuesDot, which can be used with a time solver (TS). Also is a variant to calculate accelerations (second derivatives with respect to time), MoFEM::OpCalculateVectorFieldValuesDotDot, which can be used with a second-order time solver (TS2). PETSc has currently implemented the alpha method for second-order time differential equations.

Basic usage:

constexpr int FIELD_DIM = 3;
auto vals_U_at_integration_pts = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateVectorFieldValues<FIELD_DIM>("U", vals_U_at_integration_pts));

Operator for a given field calulates values of vector field at integration points. Values are stored in, for example, vals_U_at_integration_pts. Then value can be accessed as follows

constexpr int FIELD_DIM = 3;
auto nb_integration_pts = vals_U_at_integration_pts->size2();
auto t_val = getFTensor1FromMat<FIELD_DIM>((vals_U_at_integration_pts);
for(auto gg = 0;gg!=nb_integration_pts;++gg) {
auto dot = t_val(i)*t_val(i);
++t_val;
}

Evaluate vector field gradient at integration points MoFEM::OpCalculateVectorFieldGradient

\[ u_{jk}^h(\pmb\xi^g) = \sum_i^N \phi^i_{,k}(\pmb\xi^g) \overline{u}_i^{j} \]

where \(j\) is vector element index, \(k\) is derivative and \(i\) is index of base function.

This function has a similar variant to calculate rates, when a time solver is used, see MoFEM::OpCalculateVectorFieldGradientDot, which can be used with time solver (TS).

Basic usage:

constexpr int FIELD_DIM = 3;
constexpr int SPACE_DIM = 3;
auto grad_U_at_integration_pts = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateVectorFieldGradient<FIELD_DIM, SPACE_DIM>("U", grad_U_at_integration_pts));

Operator for a given field calulates gradient of a vector field at integration points. Values are stored in for example in grad_U_at_integration_pts. Then value can be accessed as follows

constexpr int FIELD_DIM = 3;
constexpr int SPACE_DIM = 3;
auto nb_integration_pts = vals_U_at_integration_pts->size2();
auto t_grad = getFTensor2FromMat<FIELD_DIM, SPACE_DIM>((grad_U_at_integration_pts);
for(auto gg = 0;gg!=nb_integration_pts;++gg) {
auto trace = t_grad(i,j)*t_kd(i,j);
++t_grad;
}

Tensor fields

Evaluate tensor field at integration points MoFEM::OpCalculateTensor2FieldValues

\[ u_{jk}^h(\pmb\xi^g) = \sum_i^N \phi^i(\pmb\xi^g) \overline{u}_i^{jk} \]

where \(j\) and \(k\) are tensor element (coefficient) indices.

This function has a similar variant to calculate rates, when a time solver is used, see MoFEM::OpCalculateTensor2FieldValuesDot, which can be used with time solver (TS).

Basic usage:

constexpr int FIELD_DIM = 3;
auto vals_U_at_integration_pts = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateTensor2FieldValues<FIELD_DIM, FIELD_DIM>("U", vals_U_at_integration_pts));

Operator for a given field calulates values of tensor field at integration points. Values are stored in, for example, vals_U_at_integration_pts. Then value can be accessed as follows

constexpr int FIELD_DIM = 3;
auto nb_integration_pts = vals_U_at_integration_pts->size2();
auto t_val = getFTensor2FromMat<FIELD_DIM, FIELD_DIM>((vals_U_at_integration_pts);
for(auto gg = 0;gg!=nb_integration_pts;++gg) {
auto norm = sqrt(t_val(i,j)*t_val(i,j));
++t_val;
}

Evaluate symmetric tensor field at integration points MoFEM::OpCalculateTensor2SymmetricFieldValues

\[ u_{jk}^h(\pmb\xi^g) = \sum_i^N \phi^i(\pmb\xi^g) \overline{u}_i^{jk} \]

where \(j\) and \(k\) are tensor element (coefficient) indices.

Note that field degree of freedom represents only symmetric coefficients, so for each base function, we have a number of DOFs given by formula (FIELD_DIM * (FIELD_DIM + 1)) / 2. For example, for 3D case, six DOFs for each base function.

This function has a similar variant to calculate rates, when a time solver is used, see MoFEM::OpCalculateTensor2SymmetricFieldValuesDot, which can be used with a time solver (TS).

Basic usage:

constexpr int FIELD_DIM = 3;
auto vals_U_at_integration_pts = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateTensor2SymmetricFieldValues<FIELD_DIM>("U", vals_U_at_integration_pts));

Operator for a given field calulates values of tensor field at integration points. Values are stored in, for example, vals_U_at_integration_pts. Then value can be accessed as follows

constexpr int FIELD_DIM = 3;
auto nb_integration_pts = vals_U_at_integration_pts->size2();
auto t_val = getFTensor2SymmetricFromMat<FIELD_DIM>((vals_U_at_integration_pts);
for(auto gg = 0;gg!=nb_integration_pts;++gg) {
auto norm = sqrt(t_val(i,j)*t_val(i,j));
++t_val;
}

Evaluate tesorial field values for vectorial bases (Hdiv and Hcurl) MoFEM::OpCalculateHVecVectorField

\[ u_{j}^h(\pmb\xi^g) = \sum_i^N \phi^i_j(\pmb\xi^g) \overline{u}_i^{jk} \]

where \(j\) and \(k\) tensor element (coefficient) indices.

Example

constexpr int BASE_DIM = 3;
auto vals_integration_pts = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateHVecTensorField<3, 3>(
field_name, vals_integration_pts));

User data operators are defined in file UserDataOperators.hpp

Other operators

Table of other operators
Operator Usage
Face operators

Calculate Jacobian and inverse Jacobian on face. See MoFEM::OpCalculateJacForFace and MoFEM::OpInvertMatrix

For 2D element on plane XY, Jacobian is given by,

\[ J_{ij} = \frac{\partial X_i}{\partial \xi_j},\quad i,j=0,1 \]

Note
For 3D volume element, the operator is not implemented explicitly. In current implementation, Jacobian is calculated implicitly, and derivatives of shape functions are pushed to the current configuration.
For face element, it is uncertain what developer do, whether it calculates problem in 2D, or do integration on the adjacent face of 3D element, or face element is embedded in 3D, i.e. solve 2 and 1/2 dimension problem. See also MoFEM::OpCalculateInvJacForFaceEmbeddedIn3DSpace
Example usage
auto det_ptr = boost::make_shared<VectorDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateJacForFace(jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));

Push face base function to physical element configuration. MoFEM::OpSetHOInvJacToScalarBases<2> This operator assumes that the face element is on the XY plane.

\[ \frac{\partial \phi}{\partial X_i} = \frac{\partial \phi}{\partial \xi_j} J^{-1}_{ij} \]

Note
For face element, it is uncertain what developer do, whether it calculates problem in 2D, or do integration on the adjacent face of 3D element, or face element is embedded in 3D, i.e. solve 2 and 1/2 dimension problem. See also MoFEM::OpSetInvJacH1ForFaceEmbeddedIn3DSpace

Example usage

auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpSetHOInvJacToScalarBases<2>(inv_jac_ptr));

Apply contravariant (Piola) transfer to Hdiv space on face MoFEM::OpSetContravariantPiolaTransformOnFace2D and MoFEM::OpSetContravariantPiolaTransformOnFace2DEmbeddedIn3DSpace

Contravariant Piola transformation

\[ \psi_i = \frac{1}{\textrm{det}(J)}J_{ij}\hat{\psi}_j\\ \frac{\partial \psi_i}{\partial \xi_j} = \frac{1}{\textrm{det}(J)}J_{ik}\frac{\partial \hat{\psi}_k}{\partial \xi_j} \]

Example for 2D problem for face in XY space

auto jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateJacForFace(jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpMakeHdivFromHcurl());
pipeline_mng->getOpDomainRhsPipeline().push_back(

Example for face emerged in 3D space

auto det_ptr = boost::make_shared<VectorDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateJacForFaceEmbeddedIn3DSpace(jac_ptr));
pipeline.push_back(new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpMakeHdivFromHcurl());
pipeline_mng->getOpDomainRhsPipeline().push_back(
new

Note, in this case problem is 2D, however surface on which PDE is solved is emerged in 3D space.

Transform local reference derivatives of shape function to global derivatives for face MoFEM::OpSetInvJacHcurlFace and MoFEM::OpSetInvJacHcurlFaceEmbeddedIn3DSpace

\[ \frac{\partial \psi_i}{\partial X_j}= \frac{\partial \psi_i}{\partial \xi_k} \frac{\partial \xi_k}{\partial X_j} \]

Example for 2D problem for face in XY space

auto jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateJacForFace(jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpMakeHdivFromHcurl());
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(

Make H-div space from H-curl space in 2D MoFEM::OpMakeHdivFromHcurl

\[ \pmb \psi^\textrm{Hdiv} = \left\{ \psi^\textrm{Hcurl}_{1}, -\psi^\textrm{Hcurl}_{0}, 0 \right\} \]

where \(\psi^\textrm{Hdiv}\) are base functions for H-div space.

Note that applying function twice will transform base from H-curl space and then back to H-div space.

Example for 2D problem for face in XY space

auto jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateJacForFace(jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpMakeHdivFromHcurl());
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(

Modify integration weights on face to take in account higher-order geometry MoFEM::OpSetHOWeightsOnFace

Suppose Jacobian is not constant on the finite element, then this operator sets integration weights to take into account the geometry of the finite element, that is, a case for a non-parallelepiped quad. Also, that is a case for any face element given by higher-order geometry approximation. Note that geometry in MoFEM can be represented by an arbitrary order approximation field.

\[ w^g = \left(\frac{\|J(\pmb\xi^g)\|}{A}\right) W^g \]

where \(w^g\) is weight on update configuration, \(W^g\) integration weight at intergration point \(g\), at reference point \(\pmb\xi^g\). \(\|J(\pmb\xi^g)\|\) is determinant of the Jacobian, and \(A\) is area of face element.

Example

pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpSetHOWeightsOnFace());

Generic operators

Scale base (default L2) functions by inverses of measure of element MoFEM::OpScaleBaseBySpaceInverseOfMeasure

Scale base

\[ \overline{\phi}^m = \frac{1}{\mu(\Omega^e)} \phi^m \]

where \(\mu(\Omega^e)\) is length, area or volume, depending on dimension of the finite elemen entitiy. Note that function scales also base function derivatives.

Call that functions at front of the operators, before you evaluate field on given space, i.e. L2 space.

Example:

pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpScaleL2());

User data operators are defined in file UserDataOperators.hpp

Linear forms

Table of linear from operators
Operator Usage

Integrate source term for scalar and vectorial fields MoFEM::OpSourceImpl

\[ (v_i,f_i(\mathbf{x}))_\Omega = \int_\Omega v_i f_i (x,y,z) \textrm{d}\Omega \]

Template arguments:

using OpSource =
FormsIntegrators<OP>::
Assembly<ASSEMBLE_TYPE>::
LinearForm<INTEGRATION_TYPE>::
OpSource<BASE_DIM, FIELD_DIM>;

Example usage for scalar base and scalar field:

using OpDomainSource = FormsIntegrators<DomainEleOp>::Assembly<
PETSC>::LinearForm<GAUSS>::OpSource<1, 1>;
auto fun = [](double x, double y, double z) {
return sin(x/L)*sin(y/L)*sin(z/L);
};
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpDomainSource("U", fun));

Example usage for scalar base and vector field:

using OpDomainSource = FormsIntegrators<DomainEleOp>::Assembly<
PETSC>::LinearForm<GAUSS>::OpSource<1, 3>;
auto fun = [](double x, double y, double z) {
return FTensor::Tensor1<double>{ sin(x/L), sin(y/L), sin(z/L) };
};
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpDomainSource("U", fun));

Example usage for vector base and vector field:

using OpDomainSource = FormsIntegrators<DomainEleOp>::Assembly<
PETSC>::LinearForm<GAUSS>::OpSource<3, 3>;
auto fun = [](double x, double y, double z) {
return FTensor::Tensor1<double>{ sin(x/L), sin(y/L), sin(z/L) };
};
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpDomainSource("U", fun));

Integrate term for discrete scalar fields MoFEM::OpBaseTimesScalarImpl

using OpBaseTimesV = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm<
GAUSS>::OpBaseTimesScalar<BASE_DIM>;

Example

using OpBaseTimesV = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm<
GAUSS>::OpBaseTimesScalar<1>;
auto vec_at_gauss_pts = boost::make_shared<VectorDouble>();
// Push operator which sets vec_at_gauss_pts
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpBaseTimesV("U", vec_ptr_at_gauss_pts));

Integrate term for discrete vector fields MoFEM::OpBaseTimesScalarImpl

using OpBaseTimesM = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm<
GAUSS>::OpBaseTimesVector<BASE_DIM, FIELD_DIM, S>;

Example for scalar base and field vector field (e.g. body forces)

using OpBaseFlux = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm<
GAUSS>::OpBaseTimesVector<1, 3, 0>;
auto body_force = boost::make_shared<MatrixDouble>();
auto t_body_force = getFTensor1FromMat<3, 0>(*body_force);
t_body_force(0) = 0;
t_body_force(1) = 0;
t_body_force(2) = -9.81;
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpBaseFlux("U", body_force));

Example for base in Hdiv space and field in Hdiv space

using OpBaseFlux = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm<
GAUSS>::OpBaseTimesVector<3, 3, 1>;
auto flux_at_gauss_pts = boost::make_shared<MatrixDouble>();
// Push operator which sets flux at integration points
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpBaseFlux("U", flux_at_gauss_pts));

Integrate natural boundary condidition

using OpMixNaturalBC = FormsIntegrators<BoundaryEleOp>::Assembly<PETSC>::LinearForm<
GAUSS>::OpBaseTimesVector<3>;

Example

using OpTemperatureBC = FormsIntegrators<BoundaryEleOp>::Assembly<
PETSC>::LinearForm<GAUSS>::OpNormalMixVecTimesScalar<SPACE_DIM>;
pipeline.push_back(
new OpTemperatureBC("FLUX", [](double,double.double) { return 1;},
boost::make_shared<Range>(skin_edges)));

Example usage thermo_elastic.cpp

Note
Linear froms are defined in file LinearFormsIntegrators.hpp

Bilinear forms

Table of bilinear from operators
Operator Usage

Integrate Grad Grad operator MoFEM::OpGradGradImpl

FormsIntegrators<OP>::Assembly<ASSEMBLY_TYPE>::
BiLinearForm<INTEGRATION_TYPE>::
OpGradGrad<BASE_DIM, FIELD_DIM, SPACE_DIM>;

Example for scalar base, scalar field and 3D problem

auto beta = [](double x, double y , double z) {
return 1;
};
using OpDomainGradGrad = FormsIntegrators<DomainEleOp>::Assembly<
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpDomainGradGrad("U", "U", beta));

Example usage: helmholtz.cpp

Integrate mass MoFEM::OpMassImpl

using OpMass =
FormsIntegrators<OP>::Assembly<ASSEMBLY_TYPE>::
BiLinearForm<INTEGRATION_TYPE>::
OpMass<BASE_DIM, FIELD_DIM>;

Example usage: helmholtz.cpp

Example for scalar base and scalar field

auto volume_specific_heat_capacity = [](double x, double y , double z) {
return 1;
};
using OpCapacity =
FormsIntegrators<OP>::Assembly< PETSC>::BiLinearForm<GAUSS>::OpMass<1, 1>;
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpCapacity("T", "T", volume_specific_heat_capacity));

Example usage: helmholtz.cpp

Integrate base divergence times base MoFEM::OpMixDivTimesScalarImpl

using OpHcurlHcurl =
FormsIntegrators<OP>::Assembly<
ASSEMBLY_TYPE>::
BiLinearForm<INTEGRATION_TYPE>::
OpMixDivTimesScalarImpl<SPACE_DIM>;

Example usage: mixed_poisson.cpp

Example

using OpHdivU = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm<
GAUSS>::OpMixDivTimesScalar<2>;
auto unity = []() { return 1; };
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpHdivU("FLUX", "U", unity, true));

Example usage: mixed_poisson.cpp

Bilinear froms are defined in file BiLinearFormsIntegrators.hpp

MoFEM::OpNormalMixVecTimesScalarImpl
Multiply vector times normal on the face times scalar function.
Definition: LinearFormsIntegratorsImpl.hpp:378
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::OpSetContravariantPiolaTransformOnFace2D
OpSetContravariantPiolaTransformOnFace2DImpl< 2 > OpSetContravariantPiolaTransformOnFace2D
Definition: UserDataOperators.hpp:3076
MoFEM::OpFluxRhsImpl
Definition: Natural.hpp:39
MoFEM::PipelineManager::getOpDomainRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
Definition: PipelineManager.hpp:773
FTensor::Kronecker_Delta
Kronecker Delta class.
Definition: Kronecker_Delta.hpp:15
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
OpTemperatureBC
BoundaryNaturalBC::OpFlux< NaturalTemperatureMeshsets, 3, SPACE_DIM > OpTemperatureBC
Definition: seepage.cpp:148
MoFEM::OpScaleBaseBySpaceInverseOfMeasure
Scale base functions by inverses of measure of element.
Definition: HODataOperators.hpp:390
BASE_DIM
constexpr int BASE_DIM
Definition: dg_projection.cpp:15
FIELD_DIM
constexpr int FIELD_DIM
Definition: child_and_parent.cpp:15
MoFEM::OpSetContravariantPiolaTransformOnFace2DEmbeddedIn3DSpace
OpSetContravariantPiolaTransformOnFace2DImpl< 3 > OpSetContravariantPiolaTransformOnFace2DEmbeddedIn3DSpace
Definition: UserDataOperators.hpp:3078
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
OpMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition: seepage.cpp:57
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:136
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
MoFEM::PipelineManager::getOpDomainLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
Definition: PipelineManager.hpp:749
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
BiLinearForm
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', FIELD_DIM >
OpDomainGradGrad
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradGrad< 1, 1, SPACE_DIM > OpDomainGradGrad
Definition: helmholtz.cpp:25
OpHdivU
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMixDivTimesScalar< 2 > OpHdivU
Definition: mixed_poisson.cpp:25
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
OpCapacity
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, 1 > OpCapacity
Integrate Lhs base of temperature times (heat capacity) times base of temperature (T x T)
Definition: thermo_elastic.cpp:65
PlasticOps::trace
double trace(FTensor::Tensor2_symmetric< T, 2 > &t_stress)
Definition: PlasticOpsGeneric.hpp:80
OpDomainSource
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
Definition: child_and_parent.cpp:55
fun
auto fun
Function to approximate.
Definition: dg_projection.cpp:36
OpCalculateVectorFieldGradient