v0.14.0
FlatPrismPolynomialBase.cpp
Go to the documentation of this file.
1 /** \file FlatPrismPolynomialBase.cpp
2 \brief Implementation of Ainsworth-Cole H1 base on edge
3 */
4 
5 
6 
7 using namespace MoFEM;
8 
10  boost::typeindex::type_index type_index, UnknownInterface **iface) const {
11  *iface = const_cast<FlatPrismPolynomialBaseCtx *>(this);
12  return 0;
13 }
14 
17  const NumeredEntFiniteElement *fe_ptr, const FieldSpace space,
18  const FieldApproximationBase base,
19  const FieldApproximationBase copy_node_base)
20  : EntPolynomialBaseCtx(data, space, CONTINUOUS, base, copy_node_base),
21  mOab(moab), fePtr(fe_ptr) {
22  CHKERR setBase();
23  CHKERRABORT(PETSC_COMM_WORLD, ierr);
24 }
26 
28  boost::typeindex::type_index type_index, UnknownInterface **iface) const {
29  *iface = const_cast<FlatPrismPolynomialBase *>(this);
30  return 0;
31 }
32 
35 
38  boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
40 
42 
43  if (!cTx->fePtr)
44  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
45  "Pointer to element should be given "
46  "when EntPolynomialBaseCtx is constructed "
47  "(use different constructor)");
48 
49  int nb_gauss_pts = pts.size2();
50  if (!nb_gauss_pts)
52 
53  if (pts.size1() < 1)
54  SETERRQ(
55  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
56  "Wrong dimension of pts, should be at least 3 rows with coordinates");
57 
58  const FieldApproximationBase base = cTx->bAse;
59  EntitiesFieldData &data = cTx->dAta;
60 
61  if (cTx->copyNodeBase == LASTBASE)
62  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
63  else
64  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(base) =
65  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(cTx->copyNodeBase);
66 
67  data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 6, false);
68  data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 12,
69  false);
70  if (data.dataOnEntities[MBVERTEX][0].getN(base).size1() !=
71  (unsigned int)nb_gauss_pts) {
72  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
73  "Base functions or nodes has wrong number of integration points "
74  "for base %s",
76  }
77 
78  data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 6, false);
79  data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 12,
80  false);
81  N.resize(nb_gauss_pts, 3, false);
82  diffN.resize(3, 2, false);
83  CHKERR ShapeMBTRI(&*N.data().begin(), &pts(0, 0), &pts(1, 0), nb_gauss_pts);
84  std::copy(Tools::diffShapeFunMBTRI.begin(), Tools::diffShapeFunMBTRI.end(),
85  &*diffN.data().begin());
86 
87  // This is needed to have proper order of nodes on faces
88  CHKERR cTx->mOab.get_connectivity(cTx->fePtr->getEnt(), connPrism, numNodes,
89  true);
90  SideNumber_multiIndex &side_table =
92  SideNumber_multiIndex::nth_index<1>::type::iterator siit3 =
93  side_table.get<1>().find(boost::make_tuple(MBTRI, 3));
94  SideNumber_multiIndex::nth_index<1>::type::iterator siit4 =
95  side_table.get<1>().find(boost::make_tuple(MBTRI, 4));
96  if (siit3 == side_table.get<1>().end())
97  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
98  if (siit4 == side_table.get<1>().end())
99  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
100  CHKERR cTx->mOab.get_connectivity(siit3->get()->ent, connFace3, numNodes,
101  true);
102  CHKERR cTx->mOab.get_connectivity(siit4->get()->ent, connFace4, numNodes,
103  true);
104 
105  for (int nn = 0; nn < 3; nn++) {
106  faceNodes[0][nn] = std::distance(
107  connPrism, std::find(connPrism, connPrism + 3, connFace3[nn]));
108  faceNodes[1][nn] = std::distance(
109  connPrism + 3, std::find(connPrism + 3, connPrism + 6, connFace4[nn]));
110  for (int gg = 0; gg < nb_gauss_pts; gg++) {
111  double val = N(gg, nn);
112  double val_x = diffN(nn, 0);
113  double val_y = diffN(nn, 1);
114  data.dataOnEntities[MBVERTEX][0].getN(base)(gg, nn) = val;
115  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 2 * nn + 0) = val_x;
116  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 2 * nn + 1) = val_y;
117  data.dataOnEntities[MBVERTEX][0].getN(base)(gg, 3 + nn) = val;
118  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 6 + 2 * nn + 0) =
119  val_x;
120  data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 6 + 2 * nn + 1) =
121  val_y;
122  }
123  }
124 
125  switch (cTx->sPace) {
126  case H1:
127  CHKERR getValueH1(pts);
128  break;
129  case HDIV:
130  CHKERR getValueHdiv(pts);
131  break;
132  case HCURL:
133  CHKERR getValueHcurl(pts);
134  break;
135  case L2:
136  CHKERR getValueL2(pts);
137  break;
138  default:
139  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
140  }
141 
143 }
144 
147 
148  EntitiesFieldData &data = cTx->dAta;
149  const FieldApproximationBase base = cTx->bAse;
150  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
151  double *diffL, const int dim) =
153 
154  int nb_gauss_pts = pts.size2();
155 
156  // edges
157  if (data.dataOnEntities[MBEDGE].size() != 9)
158  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
159 
160  int sense[9], order[9];
161  double *H1edgeN[9], *diffH1edgeN[9];
162 
163  auto set_edge_base_data = [&](const int ee) {
165  auto &ent_data = data.dataOnEntities[MBEDGE][ee];
166  if (ent_data.getSense() == 0)
167  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
168  sense[ee] = ent_data.getSense();
169  order[ee] = ent_data.getOrder();
170  int nb_dofs = NBEDGE_H1(ent_data.getOrder());
171  ent_data.getN(base).resize(nb_gauss_pts, nb_dofs, false);
172  ent_data.getDiffN(base).resize(nb_gauss_pts, 2 * nb_dofs, false);
173  H1edgeN[ee] = &*ent_data.getN(base).data().begin();
174  diffH1edgeN[ee] = &*ent_data.getDiffN(base).data().begin();
176  };
177 
178  if ((data.spacesOnEntities[MBEDGE]).test(H1)) {
179 
180  for (int ee = 0; ee != 3; ++ee)
181  CHKERR set_edge_base_data(ee);
182  for (int ee = 6; ee != 9; ++ee)
183  CHKERR set_edge_base_data(ee);
184 
185  // shape functions on face 3
187  &sense[0], &order[0], &*N.data().begin(), &*diffN.data().begin(),
188  &H1edgeN[0], &diffH1edgeN[0], nb_gauss_pts, base_polynomials);
189  // shape functions on face 4
191  &sense[6], &order[6], &*N.data().begin(), &*diffN.data().begin(),
192  &H1edgeN[6], &diffH1edgeN[6], nb_gauss_pts, base_polynomials);
193  }
194 
195  // face
196  if (data.dataOnEntities[MBTRI].size() != 5)
197  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
198 
199  if ((data.spacesOnEntities[MBTRI]).test(H1)) {
200  for (int ff = 3; ff <= 4; ff++) {
201  int nb_dofs = NBFACETRI_H1(data.dataOnEntities[MBTRI][ff].getOrder());
202  data.dataOnEntities[MBTRI][ff].getN(base).resize(nb_gauss_pts, nb_dofs,
203  false);
204  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(nb_gauss_pts,
205  2 * nb_dofs, false);
207  faceNodes[ff - 3], data.dataOnEntities[MBTRI][ff].getOrder(),
208  &*N.data().begin(), &*diffN.data().begin(),
209  &*data.dataOnEntities[MBTRI][ff].getN(base).data().begin(),
210  &*data.dataOnEntities[MBTRI][ff].getDiffN(base).data().begin(),
211  nb_gauss_pts, base_polynomials);
212  }
213  }
214 
216 }
217 
220  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
222 }
223 
225  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
226 }
227 
229  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
230 }
UBlasMatrix< double >
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::FlatPrismPolynomialBase::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: FlatPrismPolynomialBase.cpp:27
MoFEM::FlatPrismPolynomialBase::getValueHdiv
MoFEMErrorCode getValueHdiv(MatrixDouble &pts)
Definition: FlatPrismPolynomialBase.cpp:224
H1
@ H1
continuous field
Definition: definitions.h:85
NBEDGE_H1
#define NBEDGE_H1(P)
Number of base function on edge for H1 space.
Definition: h1_hdiv_hcurl_l2.h:55
LASTBASE
@ LASTBASE
Definition: definitions.h:69
MoFEM::FlatPrismPolynomialBase::connFace3
const EntityHandle * connFace3
Definition: FlatPrismPolynomialBase.hpp:68
MoFEM::FlatPrismPolynomialBase::getValue
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Definition: FlatPrismPolynomialBase.cpp:37
MoFEM::FlatPrismPolynomialBase::~FlatPrismPolynomialBase
~FlatPrismPolynomialBase()
Definition: FlatPrismPolynomialBase.cpp:33
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
MoFEM::FlatPrismPolynomialBase::getValueL2
MoFEMErrorCode getValueL2(MatrixDouble &pts)
Definition: FlatPrismPolynomialBase.cpp:218
H1_EdgeShapeFunctions_MBTRI
PetscErrorCode H1_EdgeShapeFunctions_MBTRI(int *sense, int *p, double *N, double *diffN, double *edgeN[3], double *diff_edgeN[3], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
H1_EdgeShapeFunctions_MBTRI.
Definition: h1.c:17
MoFEM::interface_EntFiniteElement::getSideNumberTable
SideNumber_multiIndex & getSideNumberTable() const
Definition: FEMultiIndices.hpp:705
ApproximationBaseNames
const static char *const ApproximationBaseNames[]
Definition: definitions.h:72
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::FlatPrismPolynomialBase::faceNodes
int faceNodes[2][3]
Definition: FlatPrismPolynomialBase.hpp:70
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::EntPolynomialBaseCtx::dAta
EntitiesFieldData & dAta
Definition: EntPolynomialBaseCtx.hpp:36
MoFEM::interface_RefEntity::getEnt
EntityHandle getEnt() const
Get the entity handle.
Definition: RefEntsMultiIndices.hpp:603
MoFEM::FlatPrismPolynomialBase::cTx
FlatPrismPolynomialBaseCtx * cTx
Definition: FlatPrismPolynomialBase.hpp:56
MoFEM::EntitiesFieldData::spacesOnEntities
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
Definition: EntitiesFieldData.hpp:49
MoFEM::FlatPrismPolynomialBase::FlatPrismPolynomialBase
FlatPrismPolynomialBase()
Definition: FlatPrismPolynomialBase.cpp:34
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
MoFEM::EntPolynomialBaseCtx::copyNodeBase
const FieldApproximationBase copyNodeBase
Definition: EntPolynomialBaseCtx.hpp:41
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::FlatPrismPolynomialBaseCtx::FlatPrismPolynomialBaseCtx
FlatPrismPolynomialBaseCtx(EntitiesFieldData &data, moab::Interface &moab, const NumeredEntFiniteElement *fe_ptr, const FieldSpace space, const FieldApproximationBase base, const FieldApproximationBase copy_node_base=LASTBASE)
Definition: FlatPrismPolynomialBase.cpp:15
MoFEM::FlatPrismPolynomialBaseCtx::fePtr
const NumeredEntFiniteElement * fePtr
Definition: FlatPrismPolynomialBase.hpp:27
MoFEM::FlatPrismPolynomialBase::diffN
MatrixDouble diffN
Definition: FlatPrismPolynomialBase.hpp:72
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
H1_FaceShapeFunctions_MBTRI
PetscErrorCode H1_FaceShapeFunctions_MBTRI(const int *face_nodes, int p, double *N, double *diffN, double *faceN, double *diff_faceN, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition: h1.c:191
SideNumber_multiIndex
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
Definition: RefEntsMultiIndices.hpp:101
MoFEM::FlatPrismPolynomialBase
Calculate base functions on tetrahedral.
Definition: FlatPrismPolynomialBase.hpp:44
MoFEM::EntPolynomialBaseCtx::sPace
const FieldSpace sPace
Definition: EntPolynomialBaseCtx.hpp:37
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
MoFEM::Tools::diffShapeFunMBTRI
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:104
MoFEM::FlatPrismPolynomialBase::connFace4
const EntityHandle * connFace4
Definition: FlatPrismPolynomialBase.hpp:69
MoFEM::EntPolynomialBaseCtx::bAse
const FieldApproximationBase bAse
Definition: EntPolynomialBaseCtx.hpp:39
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
NBFACETRI_H1
#define NBFACETRI_H1(P)
Number of base function on triangle for H1 space.
Definition: h1_hdiv_hcurl_l2.h:60
MoFEM::FlatPrismPolynomialBaseCtx::mOab
moab::Interface & mOab
Definition: FlatPrismPolynomialBase.hpp:26
MoFEM::FlatPrismPolynomialBaseCtx::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: FlatPrismPolynomialBase.cpp:9
MoFEM::FlatPrismPolynomialBase::connPrism
const EntityHandle * connPrism
Definition: FlatPrismPolynomialBase.hpp:67
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::FlatPrismPolynomialBaseCtx
Class used to pass element data to calculate base functions on flat prism.
Definition: FlatPrismPolynomialBase.hpp:21
MoFEM::FlatPrismPolynomialBase::N
MatrixDouble N
Definition: FlatPrismPolynomialBase.hpp:71
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:57
MoFEM::NumeredEntFiniteElement
Partitioned (Indexed) Finite Element in Problem.
Definition: FEMultiIndices.hpp:728
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::FlatPrismPolynomialBaseCtx::~FlatPrismPolynomialBaseCtx
~FlatPrismPolynomialBaseCtx()
Definition: FlatPrismPolynomialBase.cpp:25
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
CONTINUOUS
@ CONTINUOUS
Regular field.
Definition: definitions.h:100
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
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEM::FlatPrismPolynomialBase::getValueHcurl
MoFEMErrorCode getValueHcurl(MatrixDouble &pts)
Definition: FlatPrismPolynomialBase.cpp:228
MoFEM::FlatPrismPolynomialBase::numNodes
int numNodes
Definition: FlatPrismPolynomialBase.hpp:66
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::EntPolynomialBaseCtx::basePolynomialsType0
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Definition: EntPolynomialBaseCtx.hpp:27
ShapeMBTRI
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182
MoFEM::FlatPrismPolynomialBase::getValueH1
MoFEMErrorCode getValueH1(MatrixDouble &pts)
Definition: FlatPrismPolynomialBase.cpp:145
MoFEM::EntPolynomialBaseCtx::setBase
MoFEMErrorCode setBase()
Definition: EntPolynomialBaseCtx.cpp:38