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 
17 
18 namespace MoFEM {
19 
20 /**
21  * @brief Bilinear 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>::BiLinearForm {
32 
33  /**
34  * @brief Integrate \f$(v_{,i},\beta(\mathbf{x}) u_{,j}))_\Omega\f$
35  * @ingroup mofem_forms
36  *
37  * @tparam SPACE_DIM
38  */
39  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
41 
42  /**
43  * @brief Integrate \f$(v_i,\beta(\mathbf{x}) u_j)_\Omega\f$
44  * @ingroup mofem_forms
45  *
46  * @note That FIELD_DIM = 4 or 9 is assumed that OpMass is for tensorial field
47  * 2x2 or 3x3
48  *
49  * @todo It should be considered another template parameter RANK which will
50  * allow to distinguish between scalar, vectorial and tensorial fields
51  *
52  * @tparam BASE_DIM dimension of base
53  * @tparam FIELD_DIM dimension of field
54  */
55  template <int BASE_DIM, int FIELD_DIM>
57 
58  /**
59  * @brief Integrate \f$(v_k,D_{ijkl} u_{,l})_\Omega\f$
60  *
61  * \note \f$D_{ijkl}\f$ is * tensor with minor & major symmetry
62  *
63  * @ingroup mofem_forms
64  *
65  * @tparam SPACE_DIM
66  * @tparam S
67  */
68  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
69  using OpGradSymTensorGrad =
71 
72  /**
73  * @brief Integrate \f$(v_{,ij},D_{ijkl} u_{,kl})_\Omega\f$
74  *
75  * \note \f$D_{ijkl}\f$ is * tensor with no symmetries
76  *
77  * @ingroup mofem_forms
78  *
79  * @tparam SPACE_DIM
80  * @tparam S
81  */
82  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
86 
87  /**
88  * @brief Integrate \f$(v_{,j},D_{ijkl} u_{,l})_\Omega\f$
89  *
90  * \note \f$D_{ijkl}\f$ is * tensor with no symmetries
91  *
92  * @ingroup mofem_forms
93  *
94  * @tparam SPACE_DIM
95  * @tparam S
96  */
97  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM, int S = 0>
98  using OpGradTensorGrad =
100 
101  /**
102  * @brief Integrate \f$(\lambda_{i,i},u)_\Omega\f$
103  *
104  * @tparam SPACE_DIM
105  */
106  template <int SPACE_DIM>
108 
109  /**
110  * @brief Integrate \f$(\lambda_{ij,j},u_{i})_\Omega\f$
111  *
112  * @tparam SPACE_DIM
113  */
114  template <int SPACE_DIM, CoordinateTypes CoordSys = CARTESIAN>
116 
117  /**
118  * @brief Integrate \f$(\lambda,u_{i,i})_\Omega\f$
119  *
120  * @tparam SPACE_DIM
121  */
122  template <int SPACE_DIM, CoordinateTypes COORDINATE_SYSTEM = CARTESIAN>
123  using OpMixScalarTimesDiv =
125 
126  /**
127  * @brief Integrate \f$(\lambda_{i},u_{,j})_\Omega\f$
128  *
129  * @tparam SPACE_DIM
130  */
131  template <int BASE_DIM, int FIELD_DIM, int SPACE_DIM>
132  using OpMixVectorTimesGrad =
134 
135  /**
136  * @brief Integrate \f$(\lambda_{ij},u_{i,j})_\Omega\f$
137  *
138  * @tparam SPACE_DIM
139  */
140  template <int SPACE_DIM>
142 };
143 
144 } // namespace MoFEM
145 
146 #endif //__BILINEAR_FORMS_INTEGRATORS_HPP__
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:147
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
EleOp
MoFEM::OpGradGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:20
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:173
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:310
MoFEM::OpMixDivTimesVecImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:242
MoFEM::OpMixScalarTimesDivImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:280
MoFEM::OpMixVectorTimesGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:305
MoFEM::OpMassImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:61
MoFEM::OpMixTensorTimesGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:350
MoFEM::OpMixDivTimesScalarImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:218
MoFEM::OpGradGradSymTensorGradGradImpl
Definition: BiLinearFormsIntegratorsImpl.hpp:195
BiLinearFormsIntegratorsImpl.hpp
Bilinear forms integrators (implementation)
MoFEM::FormsIntegrators
Integrator forms.
Definition: FormsIntegrators.hpp:299