v0.15.5
Loading...
Searching...
No Matches
src
finite_elements
BiLinearFormsIntegrators.hpp
Go to the documentation of this file.
1
/** \file BiLinearFormsIntegrators.hpp
2
* \brief Bilinear forms integrators for finite element assembly
3
* \ingroup mofem_form
4
*
5
* This file provides template aliases for bilinear form integrators used in finite
6
* element assembly. Bilinear forms represent the system matrices (left-hand side) in
7
* finite element formulations through the mathematical expression \f$(u, v)_\Omega\f$,
8
* where \f$u\f$ and \f$v\f$ are trial and test functions respectively.
9
*
10
* The integrators encompass fundamental finite element operations:
11
* - **Gradient operators**: \f$(\nabla v, \nabla u)_\Omega\f$ for diffusion and elasticity
12
* - **Mass matrices**: \f$(v, u)_\Omega\f$ for dynamic and parabolic problems
13
* - **Mixed formulations**: Coupling different field types for constrained problems
14
* - **Tensor operations**: Anisotropic material behavior with full tensor coupling
15
* - **Hybridization**: Interface constraints for discontinuous Galerkin methods
16
*
17
* All template aliases reference concrete implementations in BiLinearFormsIntegratorsImpl.hpp,
18
* providing a clean interface while maintaining computational efficiency. The operators
19
* support various coordinate systems, field dimensions, and specializations for different
20
* problem types in computational mechanics and multi-physics simulations.
21
22
\todo Some operators could be optimized. To do that, we need to write tests
23
and use Valgrind to profile code, checking cache misses. For example, some
24
operators should have iteration first over columns, then rows for optimal
25
memory access patterns. Since these operators are used in many problems,
26
implementation efficiency is critical for overall performance.
27
28
*/
29
30
#ifndef __BILINEAR_FORMS_INTEGRATORS_HPP__
31
#define __BILINEAR_FORMS_INTEGRATORS_HPP__
32
33
namespace
MoFEM
{
34
35
/**
36
* @brief Bilinear integrator form
37
* @ingroup mofem_forms
38
*
39
* @tparam EleOp
40
* @tparam A
41
* @tparam I
42
*/
43
template
<
typename
EleOp>
44
template
<AssemblyType A>
45
template
<IntegrationType I>
46
struct
FormsIntegrators
<
EleOp
>::
Assembly
<A>
::BiLinearForm
{
47
48
/**
49
* @brief Gradient-gradient bilinear form: \f$(\nabla v,\beta(\mathbf{x}) \nabla u)_\Omega\f$
50
* @ingroup mofem_forms
51
*
52
* This operator integrates the product of test function gradients with trial
53
* function gradients, scaled by a coefficient \f$\beta(\mathbf{x})\f$. This is
54
* fundamental for diffusion, heat conduction, and elasticity problems. Applications include:
55
* - Laplacian operators in diffusion equations
56
* - Stiffness matrices in structural mechanics
57
* - Gradient-based regularization terms
58
*
59
* @note \f$\beta(\mathbf{x})\f$ is a spatially varying coefficient
60
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
61
*
62
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
63
* @tparam FIELD_DIM Dimension of the field
64
* @tparam SPACE_DIM Spatial dimension of the problem
65
*/
66
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM>
67
using
OpGradGrad
=
OpGradGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase>
;
68
69
/**
70
* @brief Mass matrix bilinear form: \f$(v,\beta(\mathbf{x}) u)_\Omega\f$
71
* @ingroup mofem_forms
72
*
73
* This operator integrates the product of test and trial functions, scaled by
74
* a coefficient \f$\beta(\mathbf{x})\f$. This forms the mass matrix in dynamic
75
* problems and identity-type operators in various formulations. Applications include:
76
* - Mass matrices in dynamic finite element analysis
77
* - Identity operators in parabolic equations
78
* - Penalty terms and stabilization methods
79
* - Time-dependent term discretization
80
*
81
* @note FIELD_DIM = 4 or 9 assumes OpMass is for tensorial field 2x2 or 3x3
82
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
83
*
84
* @todo Consider another template parameter RANK to distinguish between
85
* scalar, vectorial and tensorial fields
86
*
87
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
88
* @tparam FIELD_DIM Dimension of the field
89
*/
90
template
<
int
BASE_DIM,
int
FIELD_DIM>
91
using
OpMass
=
OpMassImpl<BASE_DIM, FIELD_DIM, I, OpBase>
;
92
93
/**
94
* @brief Cached mass matrix bilinear form (same as OpMass with caching)
95
* @ingroup mofem_forms
96
*
97
* This is a cached version of OpMass that stores intermediate computations
98
* to improve performance in repeated evaluations. Useful when the same
99
* mass matrix structure is evaluated multiple times.
100
*
101
* @copydetails OpMass
102
*/
103
template
<
int
BASE_DIM,
int
FIELD_DIM>
104
using
OpMassCache
=
OpMassCacheImpl<BASE_DIM, FIELD_DIM, I, OpBase>
;
105
106
/**
107
* @brief Symmetric tensor gradient bilinear form: \f$(v_{i,j},D_{ijkl} u_{k,l})_\Omega\f$
108
* @ingroup mofem_forms
109
*
110
* This operator integrates test function gradients with trial function gradients through
111
* a symmetric fourth-order tensor \f$D_{ijkl}\f$ with minor and major symmetry.
112
* Essential for anisotropic diffusion and linear elasticity. Applications include:
113
* - Elasticity stiffness matrices with material symmetry
114
* - Anisotropic diffusion operators
115
* - Constitutive relations in continuum mechanics
116
* - Stress-strain coupling in solid mechanics
117
*
118
* @note \f$D_{ijkl}\f$ is a tensor with minor & major symmetry
119
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
120
*
121
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
122
* @tparam FIELD_DIM Dimension of the field
123
* @tparam SPACE_DIM Spatial dimension of the problem
124
* @tparam S Additional parameter for specialization (default: 0)
125
*/
126
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM,
int
S = 0>
127
using
OpGradSymTensorGrad
=
128
OpGradSymTensorGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>
;
129
130
/**
131
* @brief Second-order gradient symmetric tensor bilinear form: \f$(v_{,ij},D_{ijkl} u_{,kl})_\Omega\f$
132
* @ingroup mofem_forms
133
*
134
* This operator integrates second-order gradients (Hessians) of test and trial
135
* functions through a symmetric fourth-order tensor \f$D_{ijkl}\f$. Used in
136
* higher-order partial differential equations and plate/shell theories. Applications include:
137
* - Biharmonic operators in plate bending problems
138
* - Fourth-order diffusion equations
139
* - Higher-order regularization in optimization
140
* - Gradient damage models in mechanics
141
*
142
* @note \f$D_{ijkl}\f$ is a tensor with minor & major symmetry
143
* @note Requires C1 continuous finite elements
144
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
145
*
146
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
147
* @tparam FIELD_DIM Dimension of the field
148
* @tparam SPACE_DIM Spatial dimension of the problem
149
* @tparam S Additional parameter for specialization (default: 0)
150
*/
151
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM,
int
S = 0>
152
using
OpGradGradSymTensorGradGrad
=
153
OpGradGradSymTensorGradGradImpl
<
BASE_DIM
,
FIELD_DIM
,
SPACE_DIM
, S,
I
,
154
OpBase
>;
155
156
/**
157
* @brief Second-order gradient tensor bilinear form: \f$(v_{,ij},D_{ijkl} u_{,kl})_\Omega\f$
158
* @ingroup mofem_forms
159
*
160
* This operator integrates second-order gradients (Hessians) of test and trial
161
* functions through a fourth-order tensor \f$D_{ijkl}\f$. Used in
162
* higher-order partial differential equations and plate/shell theories. Applications include:
163
* - Biharmonic operators in plate bending problems
164
* - Fourth-order diffusion equations
165
* - Higher-order regularization in optimization
166
* - Gradient damage models in mechanics
167
*
168
* @note \f$D_{ijkl}\f$ is a tensor with no symmetries
169
* @note Requires C1 continuous finite elements
170
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
171
*
172
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
173
* @tparam FIELD_DIM Dimension of the field
174
* @tparam SPACE_DIM Spatial dimension of the problem
175
* @tparam S Additional parameter for specialization (default: 0)
176
*/
177
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM,
int
S = 0>
178
using
OpGradGradTensorGradGrad
=
179
OpGradGradTensorGradGradImpl
<
BASE_DIM
,
FIELD_DIM
,
SPACE_DIM
, S,
I
,
180
OpBase
>;
181
182
/**
183
* @brief General tensor gradient bilinear form: \f$(v_{i,j},D_{ijkl} u_{k,l})_\Omega\f$
184
* @ingroup mofem_forms
185
*
186
* This operator integrates test function gradients with trial function gradients
187
* through a general fourth-order tensor \f$D_{ijkl}\f$ with no symmetry assumptions.
188
* This is the most general form for anisotropic material behavior. Applications include:
189
* - General anisotropic elasticity without symmetry
190
* - Non-symmetric material tensors in advanced materials
191
* - Coupled multi-physics problems with complex constitutive relations
192
* - Plasticity with kinematic hardening
193
*
194
* @note \f$D_{ijkl}\f$ is a tensor with no symmetries
195
* @note More computationally expensive than symmetric versions
196
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
197
*
198
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
199
* @tparam FIELD_DIM Dimension of the field
200
* @tparam SPACE_DIM Spatial dimension of the problem
201
* @tparam S Additional parameter for specialization (default: 0)
202
*/
203
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM,
int
S = 0>
204
using
OpGradTensorGrad
=
205
OpGradTensorGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, S, I, OpBase>
;
206
207
/**
208
* @brief Mixed divergence-scalar bilinear form: \f$(\lambda_{i,i},u)_\Omega\f$
209
* @ingroup mofem_forms
210
*
211
* This operator integrates the divergence of vector test functions with scalar
212
* trial functions. Essential for mixed formulations enforcing divergence constraints
213
* such as incompressibility. Applications include:
214
* - Pressure-velocity coupling in incompressible flow
215
* - Mixed formulations for Darcy flow
216
* - Constraint enforcement in optimization problems
217
* - Divergence-free condition implementation
218
*
219
* @note \f$\lambda_{i,i}\f$ is the divergence of vector test function
220
* @note \f$u\f$ is a scalar trial function
221
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
222
*
223
* @tparam SPACE_DIM Spatial dimension of the problem
224
*/
225
template
<
int
SPACE_DIM>
226
using
OpMixDivTimesScalar
=
OpMixDivTimesScalarImpl<SPACE_DIM, I, OpBase>
;
227
228
/**
229
* @brief Mixed tensor divergence-vector bilinear form: \f$(\lambda_{ij,j},u_{i})_\Omega\f$
230
* @ingroup mofem_forms
231
*
232
* This operator integrates the divergence of tensor test functions with vector
233
* trial functions. Used in mixed formulations with stress-displacement coupling
234
* and equilibrium enforcement. Applications include:
235
* - Stress-displacement mixed formulations in elasticity
236
* - Mixed finite element methods for mechanics
237
* - Equilibrium constraint enforcement
238
* - Multi-physics coupling with tensor fields
239
*
240
* @note \f$\lambda_{ij,j}\f$ is the divergence of tensor test function
241
* @note \f$u_i\f$ is a vector trial function
242
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
243
*
244
* @tparam SPACE_DIM Spatial dimension of the problem
245
* @tparam CoordSys Coordinate system type (default: CARTESIAN)
246
*/
247
template
<
int
SPACE_DIM, CoordinateTypes CoordSys = CARTESIAN>
248
using
OpMixDivTimesVec
=
OpMixDivTimesVecImpl<SPACE_DIM, I, OpBase, CoordSys>
;
249
250
/**
251
* @brief Mixed scalar-divergence bilinear form: \f$(\lambda,u_{i,i})_\Omega\f$
252
* @ingroup mofem_forms
253
*
254
* This operator integrates scalar test functions with the divergence of vector
255
* trial functions. This is the transpose of OpMixDivTimesScalar and enforces
256
* the same constraints but with different test/trial function roles. Applications include:
257
* - Symmetric mixed formulations
258
* - Lagrange multiplier methods for divergence constraints
259
* - Saddle point problems in fluid mechanics
260
* - Mixed Poisson formulations
261
*
262
* @note \f$\lambda\f$ is a scalar test function
263
* @note \f$u_{i,i}\f$ is the divergence of vector trial function
264
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
265
*
266
* @tparam SPACE_DIM Spatial dimension of the problem
267
* @tparam COORDINATE_SYSTEM Coordinate system type (default: CARTESIAN)
268
*/
269
template
<
int
SPACE_DIM, CoordinateTypes COORDINATE_SYSTEM = CARTESIAN>
270
using
OpMixScalarTimesDiv
=
271
OpMixScalarTimesDivImpl<SPACE_DIM, I, OpBase, COORDINATE_SYSTEM>
;
272
273
/**
274
* @brief Mixed vector-gradient bilinear form: \f$(\lambda_{i},u_{,i})_\Omega\f$
275
* @ingroup mofem_forms
276
*
277
* This operator integrates vector test functions with the gradient of trial
278
* functions. Used in mixed formulations where vector fields couple with
279
* gradient terms. Applications include:
280
* - Mixed formulations with vector Lagrange multipliers
281
* - Gradient constraint enforcement
282
* - Vector-scalar coupling in multi-physics problems
283
* - Stabilization terms in finite element methods
284
*
285
* @note \f$\lambda_i\f$ is a vector test function
286
* @note \f$u_{,i}\f$ is the gradient of trial function
287
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
288
*
289
* @tparam BASE_DIM Dimension of the basis functions (1D, 2D, 3D)
290
* @tparam FIELD_DIM Dimension of the field
291
* @tparam SPACE_DIM Spatial dimension of the problem
292
*/
293
template
<
int
BASE_DIM,
int
FIELD_DIM,
int
SPACE_DIM>
294
using
OpMixVectorTimesGrad
=
295
OpMixVectorTimesGradImpl<BASE_DIM, FIELD_DIM, SPACE_DIM, I, OpBase>
;
296
297
/**
298
* @brief Mixed tensor-gradient bilinear form: \f$(\lambda_{ij},u_{i,j})_\Omega\f$
299
* @ingroup mofem_forms
300
*
301
* This operator integrates tensor test functions with the gradient of vector
302
* trial functions. Essential for mixed formulations in solid mechanics where
303
* stress tensors couple with displacement gradients. Applications include:
304
* - Stress-strain mixed formulations
305
* - Mixed finite element methods in elasticity
306
* - Constitutive relation enforcement
307
* - Coupled stress-displacement problems
308
*
309
* @note \f$\lambda_{ij}\f$ is a tensor test function
310
* @note \f$u_{i,j}\f$ is the gradient of vector trial function
311
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
312
*
313
* @tparam SPACE_DIM Spatial dimension of the problem
314
*/
315
template
<
int
SPACE_DIM>
316
using
OpMixTensorTimesGrad
=
OpMixTensorTimesGradImpl<SPACE_DIM, I, OpBase>
;
317
318
/**
319
* @brief Hybridization constraint bilinear form for broken spaces
320
* @ingroup mofem_forms
321
*
322
* This operator assembles constraint matrices during hybridization of broken
323
* finite element spaces. Hybridization introduces Lagrange multipliers on
324
* element interfaces to enforce continuity in of broken spaces fields.
325
* The operator handles the coupling between interior element solutions through
326
* Lagrange multipliers at interfaces or skelton. Applications include:
327
* - Hybridizable discontinuous Galerkin (HDG) methods
328
*
329
* @note Used specifically for broken space hybridization procedures
330
* @note Handles interface coupling in DG methods
331
* @note Implementation details are in BiLinearFormsIntegratorsImpl.hpp
332
*
333
* @tparam FIELD_DIM Dimension of the field
334
*/
335
template
<
int
FIELD_DIM>
336
using
OpBrokenSpaceConstrain
=
337
OpBrokenSpaceConstrainImpl<FIELD_DIM, I, OpBrokenBase>
;
338
};
339
340
}
// namespace MoFEM
341
342
#endif
//__BILINEAR_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
BiLinearForm
EleOp
BASE_DIM
constexpr int BASE_DIM
Definition
dg_projection.cpp:35
MoFEM
implementation of Data Operators for Forces and Sources
Definition
Common.hpp:10
I
constexpr IntegrationType I
Definition
operators_tests.cpp:31
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition
FormsIntegrators.hpp:461
MoFEM::FormsIntegrators
Integrator forms.
Definition
FormsIntegrators.hpp:450
MoFEM::OpBaseImpl
Definition
FormsIntegrators.hpp:292
MoFEM::OpBrokenSpaceConstrainImpl
Definition
FormsBrokenSpaceConstraintImpl.hpp:384
MoFEM::OpGradGradImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:23
MoFEM::OpGradGradSymTensorGradGradImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:236
MoFEM::OpGradGradTensorGradGradImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:259
MoFEM::OpGradSymTensorGradImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:188
MoFEM::OpGradTensorGradImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:214
MoFEM::OpMassCacheImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:165
MoFEM::OpMassImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:64
MoFEM::OpMixDivTimesScalarImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:281
MoFEM::OpMixDivTimesVecImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:305
MoFEM::OpMixScalarTimesDivImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:343
MoFEM::OpMixTensorTimesGradImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:413
MoFEM::OpMixVectorTimesGradImpl
Definition
BiLinearFormsIntegratorsImpl.hpp:368
Generated by
Doxygen
1.9.8 and hosted at