|
| v0.14.0
|
Face finite element.
More...
#include <src/finite_elements/FaceElementForcesAndSourcesCore.hpp>
|
virtual MoFEMErrorCode | calculateAreaAndNormalAtIntegrationPts () |
| Calculate element area and normal of the face at integration points. More...
|
|
virtual MoFEMErrorCode | calculateAreaAndNormal () |
| Calculate element area and normal of the face. More...
|
|
virtual MoFEMErrorCode | setIntegrationPts () |
| Set integration points. More...
|
|
virtual MoFEMErrorCode | getSpaceBaseAndOrderOnElement () |
| Determine approximation space and order of base functions. More...
|
|
virtual MoFEMErrorCode | calculateCoordinatesAtGaussPts () |
| Calculate coordinate at integration points. More...
|
|
Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore |
MoFEMErrorCode | getEntitySense (const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const |
| get sense (orientation) of entity More...
|
|
MoFEMErrorCode | getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const |
| Get the entity data order. More...
|
|
template<EntityType type> |
MoFEMErrorCode | getEntitySense (EntitiesFieldData &data) const |
| Get the entity sense (orientation) More...
|
|
template<EntityType type> |
MoFEMErrorCode | getEntityDataOrder (EntitiesFieldData &data, const FieldSpace space) const |
| Get the entity data order for given space. More...
|
|
MoFEMErrorCode | getFaceNodes (EntitiesFieldData &data) const |
| Get nodes on faces. More...
|
|
MoFEMErrorCode | getSpacesAndBaseOnEntities (EntitiesFieldData &data) const |
| Get field approximation space and base on entities. More...
|
|
virtual int | getRule (int order_row, int order_col, int order_data) |
| another variant of getRule More...
|
|
virtual MoFEMErrorCode | setGaussPts (int order_row, int order_col, int order_data) |
| set user specific integration rule More...
|
|
MoFEMErrorCode | calHierarchicalBaseFunctionsOnElement (const FieldApproximationBase b) |
| Calculate base functions. More...
|
|
MoFEMErrorCode | calHierarchicalBaseFunctionsOnElement () |
| Calculate base functions. More...
|
|
MoFEMErrorCode | calBernsteinBezierBaseFunctionsOnElement () |
| Calculate Bernstein-Bezier base. More...
|
|
MoFEMErrorCode | createDataOnElement (EntityType type) |
| Create a entity data on element object. More...
|
|
MoFEMErrorCode | loopOverOperators () |
| Iterate user data operators. More...
|
|
template<typename EXTRACTOR > |
MoFEMErrorCode | getNodesIndices (const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const |
| get node indices More...
|
|
MoFEMErrorCode | getRowNodesIndices (EntitiesFieldData &data, const int bit_number) const |
| get row node indices from FENumeredDofEntity_multiIndex More...
|
|
MoFEMErrorCode | getColNodesIndices (EntitiesFieldData &data, const int bit_number) const |
| get col node indices from FENumeredDofEntity_multiIndex More...
|
|
template<typename EXTRACTOR > |
MoFEMErrorCode | getEntityIndices (EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const |
|
MoFEMErrorCode | getEntityRowIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const |
|
MoFEMErrorCode | getEntityColIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const |
|
MoFEMErrorCode | getNoFieldIndices (const int bit_number, boost::shared_ptr< FENumeredDofEntity_multiIndex > dofs, VectorInt &nodes_indices) const |
| get NoField indices More...
|
|
MoFEMErrorCode | getNoFieldRowIndices (EntitiesFieldData &data, const int bit_number) const |
| get col NoField indices More...
|
|
MoFEMErrorCode | getNoFieldColIndices (EntitiesFieldData &data, const int bit_number) const |
| get col NoField indices More...
|
|
MoFEMErrorCode | getBitRefLevelOnData () |
|
MoFEMErrorCode | getNoFieldFieldData (const int bit_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const |
| Get field data on nodes. More...
|
|
MoFEMErrorCode | getNoFieldFieldData (EntitiesFieldData &data, const int bit_number) const |
|
MoFEMErrorCode | getNodesFieldData (EntitiesFieldData &data, const int bit_number) const |
| Get data on nodes. More...
|
|
MoFEMErrorCode | getEntityFieldData (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const |
|
MoFEMErrorCode | getProblemNodesIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const |
| get indices of nodal indices which are declared for problem but not this particular element More...
|
|
MoFEMErrorCode | getProblemTypeIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const |
| get indices by type (generic function) which are declared for problem but not this particular element More...
|
|
MoFEMErrorCode | getProblemNodesRowIndices (const std::string &field_name, VectorInt &nodes_indices) const |
|
MoFEMErrorCode | getProblemTypeRowIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const |
|
MoFEMErrorCode | getProblemNodesColIndices (const std::string &field_name, VectorInt &nodes_indices) const |
|
MoFEMErrorCode | getProblemTypeColIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const |
|
virtual int | getRule (int order) |
|
virtual MoFEMErrorCode | setGaussPts (int order) |
|
|
Public Types inherited from MoFEM::ForcesAndSourcesCore |
typedef boost::function< int(int order_row, int order_col, int order_data)> | RuleHookFun |
|
typedef boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> | GaussHookFun |
|
enum | KSPContext { CTX_SETFUNCTION,
CTX_OPERATORS,
CTX_KSPNONE
} |
| pass information about context of KSP/DM for with finite element is computed More...
|
|
enum | DataContext {
CTX_SET_NONE = 0,
CTX_SET_F = 1 << 0,
CTX_SET_A = 1 << 1,
CTX_SET_B = 1 << 2,
CTX_SET_X = 1 << 3,
CTX_SET_X_T = 1 << 4,
CTX_SET_X_TT = 1 << 6,
CTX_SET_TIME = 1 << 7
} |
|
using | Switches = std::bitset< 8 > |
|
enum | SNESContext { CTX_SNESSETFUNCTION,
CTX_SNESSETJACOBIAN,
CTX_SNESNONE
} |
|
enum | TSContext {
CTX_TSSETRHSFUNCTION,
CTX_TSSETRHSJACOBIAN,
CTX_TSSETIFUNCTION,
CTX_TSSETIJACOBIAN,
CTX_TSTSMONITORSET,
CTX_TSNONE
} |
|
static MoFEMErrorCode | getLibVersion (Version &version) |
| Get library version. More...
|
|
static MoFEMErrorCode | getFileVersion (moab::Interface &moab, Version &version) |
| Get database major version. More...
|
|
static MoFEMErrorCode | setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD)) |
| Get database major version. More...
|
|
static MoFEMErrorCode | getInterfaceVersion (Version &version) |
| Get database major version. More...
|
|
static constexpr Switches | CtxSetNone = PetscData::Switches(CTX_SET_NONE) |
|
static constexpr Switches | CtxSetF = PetscData::Switches(CTX_SET_F) |
|
static constexpr Switches | CtxSetA = PetscData::Switches(CTX_SET_A) |
|
static constexpr Switches | CtxSetB = PetscData::Switches(CTX_SET_B) |
|
static constexpr Switches | CtxSetX = PetscData::Switches(CTX_SET_X) |
|
static constexpr Switches | CtxSetX_T = PetscData::Switches(CTX_SET_X_T) |
|
static constexpr Switches | CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT) |
|
static constexpr Switches | CtxSetTime = PetscData::Switches(CTX_SET_TIME) |
|
Face finite element.
User is implementing own operator at Gauss point level, by own object derived from FaceElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to OpPtrVector
- Examples
- approx_sphere.cpp, boundary_marker.cpp, continuity_check_on_skeleton_3d.cpp, eigen_elastic.cpp, elasticity.cpp, EshelbianPlasticity.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hcurl_divergence_operator_2d.cpp, hdiv_divergence_operator.cpp, hello_world.cpp, helmholtz.cpp, MagneticElement.hpp, navier_stokes.cpp, phase.cpp, photon_diffusion.cpp, plot_base.cpp, PoissonOperators.hpp, prism_elements_from_surface.cpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, shallow_wave.cpp, simple_contact.cpp, simple_contact_thermal.cpp, simple_elasticity.cpp, simple_interface.cpp, and test_cache_on_entities.cpp.
Definition at line 23 of file FaceElementForcesAndSourcesCore.hpp.
◆ FaceElementForcesAndSourcesCore()
FaceElementForcesAndSourcesCore::FaceElementForcesAndSourcesCore |
( |
Interface & |
m_field | ) |
|
◆ calculateAreaAndNormal()
MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateAreaAndNormal |
( |
| ) |
|
|
protectedvirtual |
Calculate element area and normal of the face.
Note that at that point is assumed that geometry is exclusively defined by corner nodes.
- Returns
- Error code
Definition at line 113 of file FaceElementForcesAndSourcesCore.cpp.
124 auto calc_normal = [&](
const double *diff_ptr) {
135 diff_ptr, &diff_ptr[1]);
146 for (
int nn = 0; nn !=
num_nodes; ++nn) {
147 t_t1(
i) += t_coords(
i) * t_diff(N0);
148 t_t2(
i) += t_coords(
i) * t_diff(N1);
153 aRea = sqrt(t_normal(
i) * t_normal(
i));
157 const double *diff_ptr;
161 CHKERR calc_normal(diff_ptr);
167 CHKERR calc_normal(diff_ptr);
171 "Element type not implemented");
◆ calculateAreaAndNormalAtIntegrationPts()
MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateAreaAndNormalAtIntegrationPts |
( |
| ) |
|
|
protectedvirtual |
Calculate element area and normal of the face at integration points.
- Returns
- Error code
Definition at line 15 of file FaceElementForcesAndSourcesCore.cpp.
29 auto get_ftensor_n_diff = [&]() {
37 &
m(0, 0), &
m(0, 1), &
m(0, 2));
42 const size_t nb_gauss_pts =
gaussPts.size2();
58 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
67 }
else if (
type == MBQUAD) {
75 const size_t nb_gauss_pts =
gaussPts.size2();
90 auto t_diff = get_ftensor_n_diff();
91 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
92 auto t_coords = get_ftensor_from_vec_3d(
coords);
94 t_t1(
i) += t_coords(
i) * t_diff(N0);
95 t_t2(
i) += t_coords(
i) * t_diff(N1);
107 "Element type not implemented");
◆ calculateCoordinatesAtGaussPts()
MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts |
( |
| ) |
|
|
protectedvirtual |
Calculate coordinate at integration points.
- Returns
- Error code
Definition at line 346 of file FaceElementForcesAndSourcesCore.cpp.
349 const size_t nb_nodes =
351 double *shape_functions =
353 const size_t nb_gauss_pts =
gaussPts.size2();
355 for (
int gg = 0; gg != nb_gauss_pts; ++gg)
356 for (
int dd = 0;
dd != 3; ++
dd)
358 nb_nodes, &shape_functions[nb_nodes * gg], 1, &
coords[
dd], 3);
◆ getSpaceBaseAndOrderOnElement()
MoFEMErrorCode FaceElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement |
( |
| ) |
|
|
protectedvirtual |
Determine approximation space and order of base functions.
- Returns
- Error code
Definition at line 310 of file FaceElementForcesAndSourcesCore.cpp.
317 auto dim_type = CN::Dimension(
type);
319 auto get_data_on_ents = [&](
auto lower_dim,
auto space) {
322 for (
auto dd = dim_type;
dd >= lower_dim; --
dd) {
323 int nb_ents = moab::CN::NumSubEntities(
type,
dd);
324 for (
int ii = 0; ii != nb_ents; ++ii) {
325 auto sub_ent_type = moab::CN::SubEntityType(
type,
dd, ii);
327 auto &data_on_ent = data->dataOnEntities[sub_ent_type];
330 data->spacesOnEntities[sub_ent_type].set(space);
◆ operator()()
◆ setIntegrationPts()
Set integration points.
- Returns
- Error code
Definition at line 177 of file FaceElementForcesAndSourcesCore.cpp.
186 auto get_rule_by_type = [&]() {
189 return getRule(order_row + 1, order_col + 1, order_data + 1);
191 return getRule(order_row, order_col, order_data);
195 const int rule = get_rule_by_type();
197 auto set_integration_pts_for_tri = [&]() {
208 gaussPts.resize(3, nb_gauss_pts,
false);
209 cblas_dcopy(nb_gauss_pts, &
QUAD_2D_TABLE[rule]->points[1], 3,
211 cblas_dcopy(nb_gauss_pts, &
QUAD_2D_TABLE[rule]->points[2], 3,
219 cblas_dcopy(3 * nb_gauss_pts,
QUAD_2D_TABLE[rule]->points, 1, shape_ptr,
232 auto calc_base_for_tri = [&]() {
234 const size_t nb_gauss_pts =
gaussPts.size2();
237 base.resize(nb_gauss_pts, 3,
false);
238 diff_base.resize(3, 2,
false);
247 auto calc_base_for_quad = [&]() {
249 const size_t nb_gauss_pts =
gaussPts.size2();
252 base.resize(nb_gauss_pts, 4,
false);
253 diff_base.resize(nb_gauss_pts, 8,
false);
254 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
276 CHKERR set_integration_pts_for_tri();
281 CHKERR calc_base_for_quad();
285 "Element type not implemented: %d",
type);
291 const size_t nb_gauss_pts =
gaussPts.size2();
295 CHKERR calc_base_for_tri();
298 CHKERR calc_base_for_quad();
302 "Element type not implemented: %d",
type);
◆ OpCopyGeomDataToE
◆ UserDataOperator
- Examples
- boundary_marker.cpp, continuity_check_on_skeleton_3d.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_divergence_operator_2d.cpp, NavierStokesElement.hpp, PoissonOperators.hpp, quad_polynomial_approximation.cpp, reaction_diffusion.cpp, test_cache_on_entities.cpp, and UnsaturatedFlow.hpp.
Definition at line 86 of file FaceElementForcesAndSourcesCore.hpp.
◆ VolumeElementForcesAndSourcesCoreOnSide
◆ aRea
double& MoFEM::FaceElementForcesAndSourcesCore::aRea |
|
protected |
◆ conn
const EntityHandle* MoFEM::FaceElementForcesAndSourcesCore::conn |
|
protected |
◆ coords
◆ meshPositionsFieldName
std::string MoFEM::FaceElementForcesAndSourcesCore::meshPositionsFieldName |
◆ nOrmal
◆ normalsAtGaussPts
MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::normalsAtGaussPts |
|
protected |
◆ num_nodes
int MoFEM::FaceElementForcesAndSourcesCore::num_nodes |
|
protected |
◆ tangentOne
VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOne |
|
protected |
◆ tangentOneAtGaussPts
MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOneAtGaussPts |
|
protected |
◆ tangentTwo
VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwo |
|
protected |
◆ tangentTwoAtGaussPts
MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwoAtGaussPts |
|
protected |
The documentation for this struct was generated from the following files:
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
int getMaxRowOrder() const
Get max order of approximation for field in rows.
std::string meshPositionsFieldName
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
MatrixDouble tangentTwoAtGaussPts
@ L2
field with C-1 continuity
virtual MoFEMErrorCode calculateAreaAndNormal()
Calculate element area and normal of the face.
static QUAD *const QUAD_2D_TABLE[]
UBlasMatrix< double > MatrixDouble
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
int getMaxDataOrder() const
Get max order of approximation for data fields.
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
Get the entity data order.
#define QUAD_2D_TABLE_SIZE
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
int getMaxColOrder() const
Get max order of approximation for field in columns.
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
double zeta
Viscous hardening.
#define CHKERR
Inline error check.
const EntityHandle * conn
virtual moab::Interface & get_moab()=0
MatrixDouble coordsAtGaussPts
coordinated at gauss points
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
ForcesAndSourcesCore(Interface &m_field)
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
FTensor::Index< 'i', SPACE_DIM > i
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
EntitiesFieldData & dataH1
const double v
phase velocity of light in medium (cm/ns)
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
MatrixDouble tangentOneAtGaussPts
FTensor::Index< 'j', 3 > j
@ HCURL
field with continuous tangents
@ MOFEM_DATA_INCONSISTENCY
MatrixDouble gaussPts
Matrix of integration points.
UBlasVector< double > VectorDouble
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
FTensor::Index< 'm', 3 > m
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
MatrixDouble normalsAtGaussPts
FTensor::Index< 'k', 3 > k
virtual MoFEMErrorCode calculateAreaAndNormalAtIntegrationPts()
Calculate element area and normal of the face at integration points.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
@ HDIV
field with continuous normal traction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
get sense (orientation) of entity