v0.14.0
BaseDerivativesDataOperators.hpp
Go to the documentation of this file.
1 /** \file BaseDerivativesDataOperators.hpp
2  * \brief Base derivatives data operators
3 
4 */
5 
6 
7 
8 #ifndef __BASE_DERIVATIVES_DATA_OPERATORS_HPP__
9 #define __BASE_DERIVATIVES_DATA_OPERATORS_HPP__
10 
11 namespace MoFEM {
12 
14 
15  OpBaseDerivativesBase(boost::shared_ptr<MatrixDouble> base_mass_ptr,
16  boost::shared_ptr<EntitiesFieldData> data_l2,
17  const FieldApproximationBase b, const FieldSpace s,
18  int verb = QUIET, Sev sev = Sev::verbose);
19 
20 protected:
22  int verbosity;
24 
25  boost::shared_ptr<MatrixDouble> baseMassPtr;
26  boost::shared_ptr<EntitiesFieldData> dataL2;
27 
28  using GetOrderFun = boost::function<int()>;
31 
32 };
33 
34 template <int BASE_DIM>
36 
37 template <> struct OpBaseDerivativesMass<1> : public OpBaseDerivativesBase {
38 
40 
41  MoFEMErrorCode doWork(int side, EntityType type,
43 };
44 
45 template <> struct OpBaseDerivativesMass<3> : public OpBaseDerivativesMass<1> {
47 };
48 
49 template <int DIM> struct OpBaseDerivativesSetHOInvJacobian;
50 
51 template <>
53  : public OpSetInvJacSpaceForFaceImpl<2, 1> {
54 
56  boost::shared_ptr<EntitiesFieldData> data_l2,
57  boost::shared_ptr<MatrixDouble> inv_jac_ptr);
58 
59  MoFEMErrorCode doWork(int side, EntityType type,
61 
62 protected:
63  boost::shared_ptr<EntitiesFieldData> dataL2;
64 };
65 
66 template <int BASE_DIM>
68 
69 /**
70  * @brief Specialisation for calculate directives for scalar base functions
71  *
72  * @tparam
73  */
74 template <> struct OpBaseDerivativesNext<1> : public OpBaseDerivativesBase {
75 
76  OpBaseDerivativesNext(int derivative,
77  boost::shared_ptr<MatrixDouble> base_mass_ptr,
78  boost::shared_ptr<EntitiesFieldData> data_l2,
79  const FieldApproximationBase b, const FieldSpace s,
80  int verb = QUIET, Sev sev = Sev::verbose);
81 
82  MoFEMErrorCode doWork(int side, EntityType type,
84 
85 protected:
88 
89  template <int BASE_DIM>
90  MoFEMErrorCode doWorkImpl(int side, EntityType type,
92 
93  template <int BASE_DIM, int SPACE_DIM>
95  EntitiesFieldData::EntData &ent_data);
96 };
97 
98 /**
99  * @brief Specialisation for calculate directives for scalar base functions
100  *
101  * @tparam
102  */
103 template <> struct OpBaseDerivativesNext<3> : public OpBaseDerivativesNext<1> {
104 
106  MoFEMErrorCode doWork(int side, EntityType type,
108 };
109 
110 } // namespace MoFEM
111 
112 #endif //__BASE_DERIVATIVES_DATA_OPERATORS_HPP__
MoFEM::OpBaseDerivativesBase::calculateBase
MoFEMErrorCode calculateBase(GetOrderFun get_order)
Definition: BaseDerivativesDataOperators.cpp:20
MoFEM::OpBaseDerivativesBase::base
FieldApproximationBase base
Definition: BaseDerivativesDataOperators.hpp:21
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
MoFEM::OpBaseDerivativesBase::OpBaseDerivativesBase
OpBaseDerivativesBase(boost::shared_ptr< MatrixDouble > base_mass_ptr, boost::shared_ptr< EntitiesFieldData > data_l2, const FieldApproximationBase b, const FieldSpace s, int verb=QUIET, Sev sev=Sev::verbose)
Definition: BaseDerivativesDataOperators.cpp:12
MoFEM::OpBaseDerivativesBase::baseMassPtr
boost::shared_ptr< MatrixDouble > baseMassPtr
Definition: BaseDerivativesDataOperators.hpp:25
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::OpBaseDerivativesNext
Definition: BaseDerivativesDataOperators.hpp:67
MoFEM::OpBaseDerivativesBase::dataL2
boost::shared_ptr< EntitiesFieldData > dataL2
Definition: BaseDerivativesDataOperators.hpp:26
MoFEM::OpBaseDerivativesSetHOInvJacobian< 2 >::dataL2
boost::shared_ptr< EntitiesFieldData > dataL2
Definition: BaseDerivativesDataOperators.hpp:63
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
MoFEM::OpBaseDerivativesBase::calculateMass
MoFEMErrorCode calculateMass()
Definition: BaseDerivativesDataOperators.cpp:54
MoFEM::OpBaseDerivativesMass
Definition: BaseDerivativesDataOperators.hpp:35
MoFEM::OpBaseDerivativesNext< 1 >::nF
MatrixDouble nF
Definition: BaseDerivativesDataOperators.hpp:87
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
convert.type
type
Definition: convert.py:64
MoFEM::OpBaseDerivativesSetHOInvJacobian
Definition: BaseDerivativesDataOperators.hpp:49
MoFEM::LogManager::SeverityLevel
SeverityLevel
Severity levels.
Definition: LogManager.hpp:33
MoFEM::OpBaseDerivativesBase
Definition: BaseDerivativesDataOperators.hpp:13
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::OpBaseDerivativesNext< 1 >::calcBaseDerivative
int calcBaseDerivative
Definition: BaseDerivativesDataOperators.hpp:86
MoFEM::OpBaseDerivativesBase::severityLevel
Sev severityLevel
Definition: BaseDerivativesDataOperators.hpp:23
QUIET
@ QUIET
Definition: definitions.h:208
convert.int
int
Definition: convert.py:64
MoFEM::OpSetInvJacSpaceForFaceImpl
Transform local reference derivatives of shape functions to global derivatives.
Definition: UserDataOperators.hpp:2837
MoFEM::OpBaseDerivativesBase::GetOrderFun
boost::function< int()> GetOrderFun
Definition: BaseDerivativesDataOperators.hpp:28
MoFEM::OpBaseDerivativesBase::verbosity
int verbosity
Definition: BaseDerivativesDataOperators.hpp:22