v0.14.0
BiLinearFormsIntegrators.hpp
Go to the documentation of this file.
1 /** \file BiLinearFormsIntegrators.hpp
2  * \brief Bilinear forms integrators
3  * \ingroup mofem_form
4 
5  \todo SSome operators could be optimised. To do that, we need to write tests
6  and use Valgrind to profile code, shaking cache misses. For example, some
7  operators should have iteration first over columns, then rows. ome operators.
8  Since those operators are used in many problems, an implementation must be
9  efficient.
10 
11 */
12 
13 #ifndef __BILINEAR_FORMS_INTEGRATORS_HPP__
14 #define __BILINEAR_FORMS_INTEGRATORS_HPP__
15 
16 namespace MoFEM {
17 
18 /**
19  * @brief Bilinear integrator form
20  * @ingroup mofem_forms
21  *
22  * @tparam EleOp
23  * @tparam A
24  * @tparam I
25  */
26 template <typename EleOp>
27 template <AssemblyType A>
28 template <IntegrationType I>
29 struct FormsIntegrators<EleOp>::Assembly<A>::BiLinearForm {
30 
31  /**
32  * @brief Integrate \f$(v_{,i},\beta(\mathbf{x}) u_{,j}))_\Omega\f$
33  * @ingroup mofem_forms
34  *
35  * @tparam SPACE_DIM
36  */
37  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
39 
40  /**
41  * @brief Integrate \f$(v_i,\beta(\mathbf{x}) u_j)_\Omega\f$
42  * @ingroup mofem_forms
43  *
44  * @note That FIELD_DIM = 4 or 9 is assumed that OpMass is for tensorial field
45  * 2x2 or 3x3
46  *
47  * @todo It should be considered another template parameter RANK which will
48  * allow to distinguish between scalar, vectorial and tensorial fields
49  *
50  * @tparam BASE_DIM dimension of base
51  * @tparam FIELD_DIM dimension of field
52  */
53  template <int BASE_DIM, int FIELD_DIM>
55 
56  /** @copydoc OpMass */
57  template <int BASE_DIM, int FIELD_DIM>
59 
60  /**
61  * @brief Integrate \f$(v_k,D_{ijkl} u_{,l})_\Omega\f$
62  *
63  * \note \f$D_{ijkl}\f$ is * tensor with minor & major symmetry
64  *
65  * @ingroup mofem_forms
66  *
67  * @tparam SPACE_DIM
68  * @tparam S
69  */
70  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
71  using OpGradSymTensorGrad =
73 
74  /**
75  * @brief Integrate \f$(v_{,ij},D_{ijkl} u_{,kl})_\Omega\f$
76  *
77  * \note \f$D_{ijkl}\f$ is * tensor with no symmetries
78  *
79  * @ingroup mofem_forms
80  *
81  * @tparam SPACE_DIM
82  * @tparam S
83  */
84  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
88 
89  /**
90  * @brief Integrate \f$(v_{,j},D_{ijkl} u_{,l})_\Omega\f$
91  *
92  * \note \f$D_{ijkl}\f$ is * tensor with no symmetries
93  *
94  * @ingroup mofem_forms
95  *
96  * @tparam SPACE_DIM
97  * @tparam S
98  */
99  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
100  using OpGradTensorGrad =
102 
103  /**
104  * @brief Integrate \f$(\lambda_{i,i},u)_\Omega\f$
105  *
106  * @tparam SPACE_DIM
107  */
108  template <int SPACE_DIM>
110 
111  /**
112  * @brief Integrate \f$(\lambda_{ij,j},u_{i})_\Omega\f$
113  *
114  * @tparam SPACE_DIM
115  */
116  template <int SPACE_DIM, CoordinateTypes CoordSys = CARTESIAN>
118 
119  /**
120  * @brief Integrate \f$(\lambda,u_{i,i})_\Omega\f$
121  *
122  * @tparam SPACE_DIM
123  */
124  template <int SPACE_DIM, CoordinateTypes COORDINATE_SYSTEM = CARTESIAN>
125  using OpMixScalarTimesDiv =
127 
128  /**
129  * @brief Integrate \f$(\lambda_{i},u_{,j})_\Omega\f$
130  *
131  * @tparam SPACE_DIM
132  */
133  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
134  using OpMixVectorTimesGrad =
136 
137  /**
138  * @brief Integrate \f$(\lambda_{ij},u_{i,j})_\Omega\f$
139  *
140  * @tparam SPACE_DIM
141  */
142  template <int SPACE_DIM>
144 
145  /**
146  * @brief Assemble constraint matrix while hybridization
147  *
148  * @tparam FIELD_DIM
149  */
150  template <int FIELD_DIM>
151  using OpBrokenSpaceConstrain =
153 };
154 
155 } // namespace MoFEM
156 
157 #endif //__BILINEAR_FORMS_INTEGRATORS_HPP__
MoFEM::OpBrokenSpaceConstrainImpl
Definition: FormsBrokenSpaceConstraintImpl.hpp:334
BASE_DIM
constexpr int BASE_DIM
Definition: dg_projection.cpp:15
FIELD_DIM
constexpr int FIELD_DIM
Definition: child_and_parent.cpp:15
MoFEM::OpGradSymTensorGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:183
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
EleOp
MoFEM::OpGradGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:23
MoFEM::OpBaseImpl
Definition: FormsIntegrators.hpp:178
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
MoFEM::OpGradTensorGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:209
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
MoFEM::OpMixDivTimesVecImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:278
MoFEM::OpMixScalarTimesDivImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:316
MoFEM::OpMixVectorTimesGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:341
MoFEM::OpMassImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:64
MoFEM::OpMixTensorTimesGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:386
MoFEM::OpMassCacheImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:160
MoFEM::OpMixDivTimesScalarImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:254
MoFEM::OpGradGradSymTensorGradGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:231
MoFEM::FormsIntegrators
Integrator forms.
Definition: FormsIntegrators.hpp:306