v0.14.0
Loading...
Searching...
No Matches
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
11namespace 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
20protected:
24
25 boost::shared_ptr<MatrixDouble> baseMassPtr;
26 boost::shared_ptr<EntitiesFieldData> dataL2;
27
28 using GetOrderFun = boost::function<int()>;
31
32};
33
34template <int BASE_DIM>
36
37template <> struct OpBaseDerivativesMass<1> : public OpBaseDerivativesBase {
38
40
41 MoFEMErrorCode doWork(int side, EntityType type,
43};
44
45template <> struct OpBaseDerivativesMass<3> : public OpBaseDerivativesMass<1> {
46 using OpBaseDerivativesMass<1>::OpBaseDerivativesMass;
47};
48
49template <int DIM> struct OpBaseDerivativesSetHOInvJacobian;
50
51template <>
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
62protected:
63 boost::shared_ptr<EntitiesFieldData> dataL2;
64};
65
66template <int BASE_DIM>
68
69/**
70 * @brief Specialisation for calculate directives for scalar base functions
71 *
72 * @tparam
73 */
74template <> 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
85protected:
88
89 template <int BASE_DIM>
90 MoFEMErrorCode doWorkImpl(int side, EntityType type,
92
93 template <int BASE_DIM, int SPACE_DIM>
96};
97
98/**
99 * @brief Specialisation for calculate directives for scalar base functions
100 *
101 * @tparam
102 */
103template <> struct OpBaseDerivativesNext<3> : public OpBaseDerivativesNext<1> {
104
105 using OpBaseDerivativesNext<1>::OpBaseDerivativesNext;
106 MoFEMErrorCode doWork(int side, EntityType type,
108};
109
110} // namespace MoFEM
111
112#endif //__BASE_DERIVATIVES_DATA_OPERATORS_HPP__
@ QUIET
Definition: definitions.h:208
FieldApproximationBase
approximation base
Definition: definitions.h:58
FieldSpace
approximation spaces
Definition: definitions.h:82
SeverityLevel
Severity levels.
Definition: LogManager.hpp:33
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
Data on single entity (This is passed as argument to DataOperator::doWork)
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)
boost::shared_ptr< EntitiesFieldData > dataL2
MoFEMErrorCode calculateBase(GetOrderFun get_order)
boost::shared_ptr< MatrixDouble > baseMassPtr
Transform local reference derivatives of shape functions to global derivatives.