v0.15.5
Loading...
Searching...
No Matches
src
finite_elements
LinearFormsIntegrators.hpp
Go to the documentation of this file.
1
/** \file LinearFormsIntegrators.hpp
2
* \brief Linear forms integrators
3
* \ingroup mofem_form
4
*
5
* This file provides template aliases for linear form integrators used in finite
6
* element assembly. Linear forms represent the right-hand side vectors in finite
7
* element formulations \f$(v, f)_\Omega\f$ where \f$v\f$ is a test function and
8
* \f$f\f$ is a source term or field. The integrators handle various mathematical
9
* operations including scalar/vector/tensor field integration, mixed formulations,
10
* boundary conditions, and convective terms. All aliases reference implementations
11
* in LinearFormsIntegratorsImpl.hpp for actual computational details.
12
*
13
*/
14
15
#ifndef __LINEAR_FORMS_INTEGRATORS_HPP__
16
#define __LINEAR_FORMS_INTEGRATORS_HPP__
17
18
namespace
MoFEM
{
19
20
/**
21
* @brief Linear integrator form
22
* @ingroup mofem_forms
23
*
24
* @tparam EleOp
25
* @tparam A
26
* @tparam I
27
*/
28
template
<
typename
EleOp>
29
template
<AssemblyType A>
30
template
<IntegrationType I>
31
struct
FormsIntegrators<
EleOp
>::Assembly<A>
::LinearForm
{
32
33
/**
34
* @brief Integrate \f$(v,f(\mathbf{x}))_\Omega\f$ where f is a scalar function
35
* @ingroup mofem_forms
36
*
37
* This operator integrates the product of test functions with a scalar source
38
* function defined at integration points. Commonly used for:
39
* - Body forces in mechanics (gravitational, thermal expansion)
40
* - Source terms in heat transfer and diffusion problems
41
* - Prescribed loads and boundary conditions
42
*
43
* @note \f$f(\mathbf{x})\f$ is a scalar function of coordinates
44
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
45
*
46
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
47
* @tparam FIELD_DIM Dimension of the field being integrated
48
* @tparam S Source function specialization type (default: SourceFunctionSpecialization)
49
*/
50
template
<
int
BASE_DIM
,
int
FIELD_DIM
,
51
typename
S = SourceFunctionSpecialization>
52
using
OpSource
=
53
OpSourceImpl<BASE_DIM, FIELD_DIM, I, typename S::template S<OpBase>
>;
54
55
/**
56
* @brief Scalar field integrator \f$(v_i,f)_\Omega\f$ where f is a scalar
57
* @ingroup mofem_forms
58
*
59
* This operator integrates the product of vector test functions with a scalar
60
* field value. Each component of the test function is multiplied by the same
61
* scalar value. Useful for:
62
* - Uniform pressure loads in structural mechanics
63
* - Constant temperature fields in thermal problems
64
* - Scalar material properties affecting vector fields
65
*
66
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
67
*
68
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
69
* @tparam S Stride or scaling factor (default: 1)
70
*/
71
template
<
int
BASE_DIM,
int
S = 1>
72
using
OpBaseTimesScalar
=
OpBaseTimesScalarImpl<BASE_DIM, S, I, OpBase>
;
73
74
/**
75
* @brief Deprecated alias for OpBaseTimesScalar
76
* @deprecated Use OpBaseTimesScalar instead for consistency
77
*
78
* @tparam BASE_DIM Dimension of the basis functions
79
* @tparam S Stride or scaling factor (default: 1)
80
*/
81
template
<
int
BASE_DIM,
int
S = 1>
82
using
OpBaseTimesScalarField
=
OpBaseTimesScalar<BASE_DIM, S>
;
83
84
/**
85
* @brief Vector field integrator \f$(v,f_i)_\Omega\f$ where f is a vector
86
* @ingroup mofem_forms
87
*
88
* This operator integrates the product of test functions with vector field
89
* values. Each test function component is multiplied by the corresponding
90
* component of the vector field. Applications include:
91
* - Vector body forces (wind loads, electromagnetic forces)
92
* - Distributed loads with directional components
93
* - Multi-component source terms in coupled problems
94
*
95
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
96
*
97
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
98
* @tparam FIELD_DIM Dimension of the vector field
99
* @tparam S Stride or scaling factor
100
*/
101
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
S>
102
using
OpBaseTimesVector
=
103
OpBaseTimesVectorImpl<BASE_DIM, FIELD_DIM, S, I, OpBase>
;
104
//! [Source operator]
105
106
//! [Grad times tensor]
107
108
/**
109
* @brief Integrate \f$(v_{,i},f_{ij})_\Omega\f$ where f is a tensor
110
* @ingroup mofem_forms
111
*
112
* This operator integrates the product of test function gradients with tensor
113
* field values. The gradient of test functions contracts with the tensor field.
114
* Common applications include:
115
* - Stress-based formulations in mechanics
116
* - Flux calculations in transport problems
117
* - Gradient-dependent source terms
118
*
119
* @note \f$f_{ij}\f$ is a tensor field at integration points
120
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
121
*
122
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
123
* @tparam FIELD_DIM Dimension of the field
124
* @tparam SPACE_DIM Spatial dimension of the problem
125
* @tparam S Stride or scaling factor (default: 1)
126
*/
127
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM,
int
S = 1>
128
using
OpGradTimesTensor
=
129
OpGradTimesTensorImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>
;
130
131
/**
132
* @brief Integrate \f$(v_{,i},f_{ij})_\Omega\f$ where f is a symmetric tensor
133
* @ingroup mofem_forms
134
*
135
* This operator is specialized for symmetric tensors, providing computational
136
* efficiency when the tensor field has symmetry properties. Applications include:
137
* - Symmetric stress tensors in solid mechanics
138
* - Symmetric diffusion tensors in transport problems
139
* - Symmetric material property tensors
140
*
141
* @note \f$f_{ij}\f$ is a symmetric tensor field at integration points
142
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
143
*
144
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
145
* @tparam FIELD_DIM Dimension of the field
146
* @tparam SPACE_DIM Spatial dimension of the problem
147
* @tparam S Stride or scaling factor (default: 1)
148
*/
149
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM,
int
S = 1>
150
using
OpGradTimesSymTensor
=
151
OpGradTimesSymTensorImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>
;
152
153
/**
154
* @brief Mixed formulation integrator: \f$(\lambda_{ij,j},u_{i})_\Omega\f$
155
* @ingroup mofem_forms
156
*
157
* This operator integrates the divergence of a tensor field \f$\lambda_{ij}\f$
158
* with vector test functions \f$u_i\f$. Essential for mixed formulations where
159
* stress equilibrium or vector field divergence is enforced. Applications include:
160
* - Stress-displacement coupling in solid mechanics
161
* - Mixed formulations with tensor constraints
162
* - Multi-physics problems with vector-tensor coupling
163
*
164
* @note \f$\lambda_{ij,j}\f$ is the divergence of tensor field at integration points
165
* @note \f$u_i\f$ is the vector test function
166
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
167
*
168
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
169
* @tparam FIELD_DIM Dimension of the field
170
* @tparam SPACE_DIM Spatial dimension of the problem
171
* @tparam CoordSys Coordinate system type (default: CARTESIAN)
172
*/
173
template
<
int
BASE_DIM
,
int
FIELD_DIM
,
int
SPACE_DIM
,
174
CoordinateTypes
CoordSys =
CARTESIAN
>
175
using
OpMixDivTimesU
=
176
OpMixDivTimesUImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase, CoordSys>
;
177
178
/**
179
* @brief Mixed formulation integrator: \f$(\lambda_{ij},u_{i,j})_\Omega\f$
180
* @ingroup mofem_forms
181
*
182
* This operator integrates tensor field values with the gradient of vector
183
* test functions. The tensor field contracts with the gradient of test functions.
184
* Applications include:
185
* - Stress-strain coupling in mechanics
186
* - Mixed formulations with gradient constraints
187
* - Tensor-gradient field interactions
188
*
189
* @note \f$\lambda_{ij}\f$ is a tensor field at integration points
190
* @note \f$u_{i,j}\f$ is the gradient of vector test function
191
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
192
*
193
* @tparam SPACE_DIM Spatial dimension of the problem
194
*/
195
template
<
int
SPACE_DIM>
196
using
OpMixTensorTimesGradU
=
OpMixTensorTimesGradUImpl<SPACE_DIM, I, OpBase>
;
197
198
/**
199
* @brief Mixed formulation integrator: \f$(u_{i},\lambda_{ij,j})_\Omega\f$
200
* @ingroup mofem_forms
201
*
202
* This operator integrates vector test functions with the divergence of a
203
* tensor field. This is the adjoint form of OpMixDivTimesU, commonly used
204
* in mixed formulations for equilibrium enforcement. Applications include:
205
* - Adjoint stress equilibrium in mechanics
206
* - Dual formulations with tensor divergence
207
* - Vector-tensor coupling in mixed problems
208
*
209
* @note \f$u_i\f$ is a vector test function
210
* @note \f$\lambda_{ij,j}\f$ is the divergence of tensor field at integration points
211
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
212
*
213
* @tparam SPACE_DIM Spatial dimension of the problem
214
*/
215
template
<
int
SPACE_DIM>
216
using
OpMixVecTimesDivLambda
=
217
OpMixVecTimesDivLambdaImpl<SPACE_DIM, I, OpBase>
;
218
219
/**
220
* @brief Boundary integrator: normal vector times scalar function
221
* @ingroup mofem_forms
222
*
223
* This operator multiplies vector test functions by the normal vector on
224
* a face and integrates with a scalar field. Essential for natural boundary
225
* conditions in mixed formulations. Applications include:
226
* - Natural boundary conditions in fluid mechanics
227
* - Surface flux conditions in transport problems
228
* - Mixed boundary value problems
229
*
230
* @note This operator is used on domain boundaries (faces)
231
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
232
*
233
* @tparam SPACE_DIM Spatial dimension of the problem
234
*/
235
template
<
int
SPACE_DIM>
236
using
OpNormalMixVecTimesScalar
=
237
OpNormalMixVecTimesScalarImpl<SPACE_DIM, I, OpBase>
;
238
239
/**
240
* @brief Boundary integrator: normal vector times vector field
241
* @ingroup mofem_forms
242
*
243
* This operator multiplies vector test functions by the normal vector on
244
* a face and integrates with a vector field. Used for natural boundary
245
* conditions involving vector quantities. Applications include:
246
* - Vector-valued boundary conditions in mechanics
247
* - Surface traction conditions in solid mechanics
248
* - Vector flux conditions in multi-physics problems
249
*
250
* @note This operator is used on domain boundaries (faces)
251
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
252
*
253
* @tparam SPACE_DIM Spatial dimension of the problem
254
*/
255
template
<
int
SPACE_DIM>
256
using
OpNormalMixVecTimesVectorField
=
257
OpNormalMixVecTimesVectorFieldImpl<SPACE_DIM, I, OpBase>
;
258
259
/**
260
* @brief Convective term integrator: \f$(v, u_i \mathbf{y}_{,i})\f$
261
* @ingroup mofem_forms
262
*
263
* This operator integrates convective terms where a velocity field \f$u_i\f$
264
* transports a scalar or vector field \f$\mathbf{y}\f$. The convective operator
265
* \f$u_i \mathbf{y}_{,i}\f$ represents advection. Applications include:
266
* - Advection in fluid mechanics
267
* - Convective transport in heat transfer
268
* - Material transport in multi-physics problems
269
*
270
* @note \f$u_i\f$ is the velocity field at integration points
271
* @note \f$\mathbf{y}_{,i}\f$ is the gradient of the transported field
272
* @note \f$\mathbf{y}\f$ can be scalar or vector field
273
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
274
*
275
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
276
* @tparam FIELD_DIM Dimension of the field
277
* @tparam SPACE_DIM Spatial dimension of the problem
278
*/
279
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM>
280
using
OpConvectiveTermRhs
=
281
OpConvectiveTermRhsImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase>
;
282
283
/**
284
* @brief Hybridization constraint integrator for broken spaces
285
* @ingroup mofem_forms
286
*
287
* This operator assembles the right-hand side for constraint matrix C
288
* during hybridization of broken finite element spaces. Hybridization
289
* introduces Lagrange multipliers to enforce continuity across element
290
* boundaries in discontinuous Galerkin methods. Applications include:
291
* - Discontinuous Galerkin methods
292
* - Mixed finite element methods with broken spaces
293
* - Hybridizable discontinuous Galerkin (HDG) methods
294
*
295
* @note Used specifically for broken space hybridization
296
* @note Implementation details are in LinearFormsIntegratorsImpl.hpp
297
*
298
* @tparam FIELD_DIM Dimension of the field
299
*/
300
template
<
int
FIELD_DIM>
301
using
OpBrokenSpaceConstrainDHybrid
=
302
OpBrokenSpaceConstrainDHybridImpl<FIELD_DIM, I, OpBase>
;
303
304
/**
305
* @brief Assemble Rhs for contraint matrix C^T whole hybridisation
306
*
307
* @tparam FIELD_DIM
308
*/
309
template
<
int
FIELD_DIM>
310
using
OpBrokenSpaceConstrainDFlux
=
311
OpBrokenSpaceConstrainDFluxImpl<FIELD_DIM, I, OpBrokenBase>
;
312
313
/**
314
* @brief Assemble Rhs for contraint matrix C^T whole hybridisation with topological derivative
315
*
316
* @tparam FIELD_DIM
317
*/
318
template
<
int
FIELD_DIM>
319
using
OpTopoDerivativeBrokenSpaceConstrainDHybrid
=
320
OpTopoDerivativeBrokenSpaceConstrainDHybridImpl<FIELD_DIM, I, OpBase>
;
321
322
/**
323
* @brief Assemble Rhs for contraint matrix C^T whole hybridisation with topological derivative
324
*
325
* @tparam FIELD_DIM
326
*/
327
template
<
int
FIELD_DIM>
328
using
OpTopoDerivativeBrokenSpaceConstrainDFlux
=
329
OpTopoDerivativeBrokenSpaceConstrainDFluxImpl<FIELD_DIM, I, OpBrokenBase>
;
330
};
331
332
}
// namespace MoFEM
333
334
#endif
// __LINEAR_FORMS_INTEGRATORS_HPP__
SPACE_DIM
constexpr int SPACE_DIM
Definition
child_and_parent.cpp:17
FIELD_DIM
constexpr int FIELD_DIM
Definition
child_and_parent.cpp:16
EleOp
LinearForm
CoordinateTypes
CoordinateTypes
Coodinate system.
Definition
definitions.h:127
CARTESIAN
@ CARTESIAN
Definition
definitions.h:128
BASE_DIM
constexpr int BASE_DIM
Definition
dg_projection.cpp:35
MoFEM
implementation of Data Operators for Forces and Sources
Definition
Common.hpp:10
MoFEM::OpBaseTimesScalarImpl
Definition
LinearFormsIntegratorsImpl.hpp:121
MoFEM::OpBaseTimesVectorImpl
Definition
LinearFormsIntegratorsImpl.hpp:141
MoFEM::OpBrokenSpaceConstrainDFluxImpl
Definition
FormsBrokenSpaceConstraintImpl.hpp:449
MoFEM::OpBrokenSpaceConstrainDHybridImpl
Definition
FormsBrokenSpaceConstraintImpl.hpp:446
MoFEM::OpConvectiveTermRhsImpl
Definition
LinearFormsIntegratorsImpl.hpp:462
MoFEM::OpGradTimesSymTensorImpl
Definition
LinearFormsIntegratorsImpl.hpp:228
MoFEM::OpGradTimesTensorImpl
Definition
LinearFormsIntegratorsImpl.hpp:185
MoFEM::OpMixDivTimesUImpl
Definition
LinearFormsIntegratorsImpl.hpp:250
MoFEM::OpMixTensorTimesGradUImpl
Tensor field time gradient of vector field.
Definition
LinearFormsIntegratorsImpl.hpp:351
MoFEM::OpMixVecTimesDivLambdaImpl
Vector field time divergence of tensor.
Definition
LinearFormsIntegratorsImpl.hpp:317
MoFEM::OpNormalMixVecTimesScalarImpl
Multiply vector times normal on the face times scalar function.
Definition
LinearFormsIntegratorsImpl.hpp:385
MoFEM::OpNormalMixVecTimesVectorFieldImpl
Multiply vector times normal on the face times vector field.
Definition
LinearFormsIntegratorsImpl.hpp:442
MoFEM::OpSourceImpl
Definition
LinearFormsIntegratorsImpl.hpp:23
MoFEM::OpTopoDerivativeBrokenSpaceConstrainDFluxImpl
Definition
FormsBrokenSpaceConstraintImpl.hpp:591
MoFEM::OpTopoDerivativeBrokenSpaceConstrainDHybridImpl
Definition
FormsBrokenSpaceConstraintImpl.hpp:588
Generated by
Doxygen
1.9.8 and hosted at