v0.14.0
Public Member Functions | Private Member Functions | Private Attributes | List of all members
SomeUserPolynomialBase Struct Reference

Class used to calculate base functions at integration points. More...

Inheritance diagram for SomeUserPolynomialBase:
[legend]
Collaboration diagram for SomeUserPolynomialBase:
[legend]

Public Member Functions

 SomeUserPolynomialBase ()=default
 
 ~SomeUserPolynomialBase ()=default
 
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 Return interface to this class when one ask for for tetrahedron, otherisw return interface class for generic class. More...
 
MoFEMErrorCode getValue (MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 Calculate base functions at intergeneration points. More...
 
- Public Member Functions inherited from MoFEM::BaseFunction
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, MoFEM::UnknownInterface **iface) const
 
virtual MoFEMErrorCode getValue (MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
virtual MoFEMErrorCode getValue (MatrixDouble &pts_x, MatrixDouble &pts_t, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
- Public Member Functions inherited from MoFEM::BaseFunctionUnknownInterface
virtual ~BaseFunctionUnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 

Private Member Functions

MoFEMErrorCode getValueHdivForCGGBubble (MatrixDouble &pts)
 

Private Attributes

EntPolynomialBaseCtxcTx
 
MatrixDouble shapeFun
 

Additional Inherited Members

- Public Types inherited from MoFEM::BaseFunction
using DofsSideMap = multi_index_container< DofsSideMapData, indexed_by< ordered_non_unique< tag< TypeSide_mi_tag >, composite_key< DofsSideMapData, member< DofsSideMapData, EntityType, &DofsSideMapData::type >, member< DofsSideMapData, int, &DofsSideMapData::side > >>, ordered_unique< tag< EntDofIdx_mi_tag >, member< DofsSideMapData, int, &DofsSideMapData::dof > > > >
 Map entity stype and side to element/entity dof index. More...
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
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...
 

Detailed Description

Class used to calculate base functions at integration points.

Examples
forces_and_sources_testing_users_base.cpp.

Definition at line 37 of file forces_and_sources_testing_users_base.cpp.

Constructor & Destructor Documentation

◆ SomeUserPolynomialBase()

SomeUserPolynomialBase::SomeUserPolynomialBase ( )
default

◆ ~SomeUserPolynomialBase()

SomeUserPolynomialBase::~SomeUserPolynomialBase ( )
default

Member Function Documentation

◆ getValue()

MoFEMErrorCode SomeUserPolynomialBase::getValue ( MatrixDouble &  pts,
boost::shared_ptr< BaseFunctionCtx ctx_ptr 
)
inline

Calculate base functions at intergeneration points.

Parameters
pts
ctx_ptr
Returns
MoFEMErrorCode

Definition at line 63 of file forces_and_sources_testing_users_base.cpp.

64  {
66 
68 
69  int nb_gauss_pts = pts.size2();
70  if (!nb_gauss_pts) {
72  }
73 
74  if (pts.size1() < 3) {
75  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
76  "Wrong dimension of pts, should be at least 3 rows with "
77  "coordinates");
78  }
79 
80  switch (cTx->sPace) {
81  case HDIV:
83  break;
84  default:
85  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
86  }
87 
89  }

◆ getValueHdivForCGGBubble()

MoFEMErrorCode SomeUserPolynomialBase::getValueHdivForCGGBubble ( MatrixDouble &  pts)
inlineprivate

Definition at line 96 of file forces_and_sources_testing_users_base.cpp.

96  {
98 
99  const FieldApproximationBase base = cTx->bAse;
100  // This should be used only in case USER_BASE is selected
101  if (cTx->bAse != USER_BASE) {
102  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
103  "Wrong base, should be USER_BASE");
104  }
105 
106  // This is example, simply use Demkowicz HDiv base to generate base
107  // functions
108 
109  EntitiesFieldData &data = cTx->dAta;
110  int nb_gauss_pts = pts.size2();
111 
112  // calculate shape functions, i.e. barycentric coordinates
113  shapeFun.resize(nb_gauss_pts, 4, false);
114  CHKERR ShapeMBTET(&*shapeFun.data().begin(), &pts(0, 0), &pts(1, 0),
115  &pts(2, 0), nb_gauss_pts);
116  // derivatives of shape functions
117  double diff_shape_fun[12];
118  CHKERR ShapeDiffMBTET(diff_shape_fun);
119 
120  int volume_order = data.dataOnEntities[MBTET][0].getOrder();
121 
122  int p_f[4];
123  double *phi_f[4];
124  double *diff_phi_f[4];
125 
126  // Calculate base function on tet faces
127  for (int ff = 0; ff != 4; ff++) {
128  int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
129  int order = volume_order > face_order ? volume_order : face_order;
130  data.dataOnEntities[MBTRI][ff].getN(base).resize(
131  nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
132  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
133  nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
134  p_f[ff] = order;
135  phi_f[ff] = &*data.dataOnEntities[MBTRI][ff].getN(base).data().begin();
136  diff_phi_f[ff] =
137  &*data.dataOnEntities[MBTRI][ff].getDiffN(base).data().begin();
139  continue;
141  &data.facesNodes(ff, 0), order, &*shapeFun.data().begin(),
142  diff_shape_fun, phi_f[ff], diff_phi_f[ff], nb_gauss_pts, 4);
143  }
144 
145  // Calculate base functions in tet interior
146  if (NBVOLUMETET_DEMKOWICZ_HDIV(volume_order) > 0) {
147  data.dataOnEntities[MBTET][0].getN(base).resize(
148  nb_gauss_pts, 3 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
149  data.dataOnEntities[MBTET][0].getDiffN(base).resize(
150  nb_gauss_pts, 9 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
151  double *phi_v = &*data.dataOnEntities[MBTET][0].getN(base).data().begin();
152  double *diff_phi_v =
153  &*data.dataOnEntities[MBTET][0].getDiffN(base).data().begin();
155  volume_order, &*shapeFun.data().begin(), diff_shape_fun, p_f, phi_f,
156  diff_phi_f, phi_v, diff_phi_v, nb_gauss_pts);
157  }
158 
159  // Set size of face base correctly
160  for (int ff = 0; ff != 4; ff++) {
161  int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
162  data.dataOnEntities[MBTRI][ff].getN(base).resize(
163  nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
164  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
165  nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
166  }
167 
169  }

◆ query_interface()

MoFEMErrorCode SomeUserPolynomialBase::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
inlinevirtual

Return interface to this class when one ask for for tetrahedron, otherisw return interface class for generic class.

Parameters
ifaceinterface class
Returns
MoFEMErrorCode

Implements MoFEM::BaseFunctionUnknownInterface.

Definition at line 49 of file forces_and_sources_testing_users_base.cpp.

50  {
52  *iface = const_cast<SomeUserPolynomialBase *>(this);
54  }

Member Data Documentation

◆ cTx

EntPolynomialBaseCtx* SomeUserPolynomialBase::cTx
private

Definition at line 92 of file forces_and_sources_testing_users_base.cpp.

◆ shapeFun

MatrixDouble SomeUserPolynomialBase::shapeFun
private

Definition at line 94 of file forces_and_sources_testing_users_base.cpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::Hdiv_Demkowicz_Face_MBTET_ON_FACE
MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb)
Definition: Hdiv.cpp:634
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
MoFEM::EntPolynomialBaseCtx::dAta
EntitiesFieldData & dAta
Definition: EntPolynomialBaseCtx.hpp:36
MoFEM::Hdiv_Demkowicz_Interior_MBTET
MoFEMErrorCode Hdiv_Demkowicz_Interior_MBTET(int p, double *N, double *diffN, int p_face[], double *phi_f[4], double *diff_phi_f[4], double *phi_v, double *diff_phi_v, int gdim)
Definition: Hdiv.cpp:780
MoFEM::EntitiesFieldData::facesNodes
MatrixInt facesNodes
nodes on finite element faces
Definition: EntitiesFieldData.hpp:45
USER_BASE
@ USER_BASE
user implemented approximation base
Definition: definitions.h:68
order
constexpr int order
Definition: dg_projection.cpp:18
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
NBFACETRI_DEMKOWICZ_HDIV
#define NBFACETRI_DEMKOWICZ_HDIV(P)
Definition: h1_hdiv_hcurl_l2.h:139
MoFEM::EntPolynomialBaseCtx::sPace
const FieldSpace sPace
Definition: EntPolynomialBaseCtx.hpp:37
SomeUserPolynomialBase::shapeFun
MatrixDouble shapeFun
Definition: forces_and_sources_testing_users_base.cpp:94
NBVOLUMETET_DEMKOWICZ_HDIV
#define NBVOLUMETET_DEMKOWICZ_HDIV(P)
Definition: h1_hdiv_hcurl_l2.h:140
MoFEM::EntPolynomialBaseCtx::bAse
const FieldApproximationBase bAse
Definition: EntPolynomialBaseCtx.hpp:39
ShapeDiffMBTET
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:319
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:57
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
SomeUserPolynomialBase::getValueHdivForCGGBubble
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
Definition: forces_and_sources_testing_users_base.cpp:96
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
SomeUserPolynomialBase
Class used to calculate base functions at integration points.
Definition: forces_and_sources_testing_users_base.cpp:37
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
SomeUserPolynomialBase::cTx
EntPolynomialBaseCtx * cTx
Definition: forces_and_sources_testing_users_base.cpp:92
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
ShapeMBTET
PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition: fem_tools.c:306