v0.14.0
EdgePolynomialBase.cpp
Go to the documentation of this file.
1 /** \file EdgePolynomialBase.cpp
2 \brief Implementation of Ainsworth-Cole H1 base on edge
3 */
4 
5 using namespace MoFEM;
6 
8 EdgePolynomialBase::query_interface(boost::typeindex::type_index type_index,
9  UnknownInterface **iface) const {
10  *iface = const_cast<EdgePolynomialBase *>(this);
11  return 0;
12 }
13 
16  boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
18 
20 
21  int nb_gauss_pts = pts.size2();
22  if (!nb_gauss_pts)
24 
25  if (pts.size1() < 1)
26  SETERRQ(
27  PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
28  "Wrong dimension of pts, should be at least 3 rows with coordinates");
29 
30  const FieldApproximationBase base = cTx->bAse;
31  EntitiesFieldData &data = cTx->dAta;
32 
33  if (base != AINSWORTH_BERNSTEIN_BEZIER_BASE) {
34  if (cTx->copyNodeBase == LASTBASE) {
35  data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 2,
36  false);
38  &*data.dataOnEntities[MBVERTEX][0].getN(base).data().begin(),
39  &pts(0, 0), nb_gauss_pts);
40  } else
41  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(base) =
42  data.dataOnEntities[MBVERTEX][0].getNSharedPtr(cTx->copyNodeBase);
43 
44  if (data.dataOnEntities[MBVERTEX][0].getN(base).size1() !=
45  (unsigned int)nb_gauss_pts)
46  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
47  "Base functions or nodes has wrong number of integration points "
48  "for base %s",
50 
51  data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(nb_gauss_pts, 2 * 1,
52  false);
53  for (auto gg = 0; gg != nb_gauss_pts; ++gg)
54  std::copy(Tools::diffShapeFunMBEDGE.begin(),
56  &data.dataOnEntities[MBVERTEX][0].getDiffN(base)(gg, 0));
57  }
58 
59  switch (cTx->sPace) {
60  case H1:
61  CHKERR getValueH1(pts);
62  break;
63  case HDIV:
64  CHKERR getValueHdiv(pts);
65  break;
66  case HCURL:
67  CHKERR getValueHcurl(pts);
68  break;
69  case L2:
70  CHKERR getValueL2(pts);
71  break;
72  default:
73  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
74  }
75 
77 }
78 
81 
82  switch (cTx->bAse) {
86  break;
89  break;
92  break;
93  default:
94  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
95  }
96 
98 }
99 
103 
104  EntitiesFieldData &data = cTx->dAta;
105  const std::string &field_name = cTx->fieldName;
106  const int nb_gauss_pts = pts.size2();
107 
108  auto get_alpha = [field_name](auto &data) -> MatrixInt & {
109  auto &ptr = data.getBBAlphaIndicesSharedPtr(field_name);
110  if (!ptr)
111  ptr.reset(new MatrixInt());
112  return *ptr;
113  };
114 
115  auto get_base = [field_name](auto &data) -> MatrixDouble & {
116  auto &ptr = data.getBBNSharedPtr(field_name);
117  if (!ptr)
118  ptr.reset(new MatrixDouble());
119  return *ptr;
120  };
121 
122  auto get_diff_base = [field_name](auto &data) -> MatrixDouble & {
123  auto &ptr = data.getBBDiffNSharedPtr(field_name);
124  if (!ptr)
125  ptr.reset(new MatrixDouble());
126  return *ptr;
127  };
128 
129  auto &get_n = get_base(data.dataOnEntities[MBVERTEX][0]);
130  auto &get_diff_n = get_diff_base(data.dataOnEntities[MBVERTEX][0]);
131  get_n.resize(nb_gauss_pts, 2, false);
132  get_diff_n.resize(nb_gauss_pts, 2, false);
133  get_n.clear();
134  get_diff_n.clear();
135 
136  if (data.dataOnEntities[MBVERTEX][0].getN(NOBASE).size1() !=
137  (unsigned int)nb_gauss_pts)
138  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
139  "Base functions or nodes has wrong number of integration points "
140  "for base %s",
142  auto &lambda = data.dataOnEntities[MBVERTEX][0].getN(NOBASE);
143 
144  auto &vertex_alpha = get_alpha(data.dataOnEntities[MBVERTEX][0]);
145  vertex_alpha.resize(2, 2, false);
146  vertex_alpha(0, 0) = data.dataOnEntities[MBVERTEX][0].getBBNodeOrder()[0];
147  vertex_alpha(0, 1) = 0;
148  vertex_alpha(1, 0) = 0;
149  vertex_alpha(1, 1) = data.dataOnEntities[MBVERTEX][0].getBBNodeOrder()[1];
150 
152  1, nb_gauss_pts, vertex_alpha.size1(), &vertex_alpha(0, 0), &lambda(0, 0),
153  Tools::diffShapeFunMBEDGE.data(), &get_n(0, 0), &get_diff_n(0, 0));
154  std::array<double, 2> f = {
155  boost::math::factorial<double>(
156  data.dataOnEntities[MBVERTEX][0].getBBNodeOrder()[0]),
157  boost::math::factorial<double>(
158  data.dataOnEntities[MBVERTEX][0].getBBNodeOrder()[1])};
159 
160  for (int g = 0; g != nb_gauss_pts; ++g)
161  for (int n = 0; n != 2; ++n) {
162  get_n(g, n) *= f[n];
163  get_diff_n(g, n) *= f[n];
164  }
165 
166  if (data.spacesOnEntities[MBEDGE].test(H1)) {
167  if (data.dataOnEntities[MBEDGE].size() != 1)
168  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
169  "Wrong size ent of ent data");
170 
171  int order = data.dataOnEntities[MBEDGE][0].getOrder();
172  const int nb_dofs = NBEDGE_H1(order);
173 
174  auto &get_n = get_base(data.dataOnEntities[MBEDGE][0]);
175  auto &get_diff_n = get_diff_base(data.dataOnEntities[MBEDGE][0]);
176  get_n.resize(nb_gauss_pts, nb_dofs, false);
177  get_diff_n.resize(nb_gauss_pts, nb_dofs, false);
178 
179  if (nb_dofs) {
180  auto &edge_alpha = get_alpha(data.dataOnEntities[MBEDGE][0]);
181  edge_alpha.resize(nb_dofs, 2);
184  order, nb_gauss_pts, edge_alpha.size1(), &edge_alpha(0, 0),
185  &lambda(0, 0), Tools::diffShapeFunMBEDGE.data(), &get_n(0, 0),
186  &get_diff_n(0, 0));
187  }
188  } else {
189  get_n.resize(nb_gauss_pts, 0, false);
190  get_diff_n.resize(nb_gauss_pts, 0, false);
191  }
192 
194 }
195 
198 
199  EntitiesFieldData &data = cTx->dAta;
200  const FieldApproximationBase base = cTx->bAse;
201  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
202  double *diffL, const int dim) =
204 
205  int nb_gauss_pts = pts.size2();
206 
207  const int side_number = 0;
208  int sense = data.dataOnEntities[MBEDGE][side_number].getSense();
209  int order = data.dataOnEntities[MBEDGE][side_number].getOrder();
210  data.dataOnEntities[MBEDGE][side_number].getN(base).resize(
211  nb_gauss_pts, NBEDGE_H1(order), false);
212  data.dataOnEntities[MBEDGE][side_number].getDiffN(base).resize(
213  nb_gauss_pts, NBEDGE_H1(order), false);
214 
215  data.dataOnEntities[MBEDGE][side_number].getN(base).clear();
216  data.dataOnEntities[MBEDGE][side_number].getDiffN(base).clear();
217 
218  L.resize(NBEDGE_H1(order), false);
219  diffL.resize(NBEDGE_H1(order), false);
220 
221  if (data.dataOnEntities[MBEDGE][side_number].getOrder() > 1) {
222 
223  double diff_s = 2.; // s = s(xi), ds/dxi = 2., because change of basis
224  for (int gg = 0; gg != nb_gauss_pts; gg++) {
225 
226  double s = 2 * pts(0, gg) - 1; // makes form -1..1
227  if (!sense) {
228  s *= -1;
229  diff_s *= -1;
230  }
231 
232  // calculate Legendre polynomials at integration points
233  CHKERR base_polynomials(NBEDGE_H1(order) - 1, s, &diff_s,
234  &*L.data().begin(), &*diffL.data().begin(), 1);
235 
236  for (unsigned int pp = 0;
237  pp < data.dataOnEntities[MBEDGE][side_number].getN(base).size2();
238  pp++) {
239 
240  // Calculate edge shape functions N0*N1*L(p), where N0 and N1 are nodal
241  // shape functions
242  double v = data.dataOnEntities[MBVERTEX][0].getN(base)(gg, 0) *
243  data.dataOnEntities[MBVERTEX][0].getN(base)(gg, 1);
244  data.dataOnEntities[MBEDGE][side_number].getN(base)(gg, pp) = v * L(pp);
245 
246  // Calculate derivative edge shape functions
247  // dN/dksi = dN0/dxi*N1*L + N0*dN1/ksi*L + N0*N1*dL/dxi
248  data.dataOnEntities[MBEDGE][side_number].getDiffN(base)(gg, pp) =
249  ((+1.) * data.dataOnEntities[MBVERTEX][0].getN(base)(gg, 1) +
250  data.dataOnEntities[MBVERTEX][0].getN(base)(gg, 0) * (-1.)) *
251  L(pp) +
252  v * diffL(pp);
253  }
254  }
255  }
256 
258 }
259 
262 
263  EntitiesFieldData &data = cTx->dAta;
264  const FieldApproximationBase base = cTx->bAse;
265 
266  int nb_gauss_pts = pts.size2();
267 
268  const int side_number = 0;
269  int order = data.dataOnEntities[MBEDGE][side_number].getOrder();
270 
271  data.dataOnEntities[MBEDGE][side_number].getN(base).resize(
272  nb_gauss_pts, NBEDGE_H1(order), false);
273  data.dataOnEntities[MBEDGE][side_number].getDiffN(base).resize(
274  nb_gauss_pts, NBEDGE_H1(order), false);
275 
276  data.dataOnEntities[MBEDGE][side_number].getN(base).clear();
277  data.dataOnEntities[MBEDGE][side_number].getDiffN(base).clear();
278 
279  double diff_shape[] = {-1, 1};
280  MatrixDouble shape(nb_gauss_pts, 2);
281 
282  if (NBEDGE_L2(order)) {
283  for (int gg = 0; gg != nb_gauss_pts; gg++) {
284  const double mu0 = 1.0 - pts(0, gg); // pts ranges over [0, 1]
285  const double mu1 = pts(0, gg);
286  shape(gg, 0) = mu0;
287  shape(gg, 1) = mu1;
288  }
289 
291  order, &*shape.data().begin(), diff_shape,
292  &*data.dataOnEntities[MBEDGE][side_number].getN(base).data().begin(),
293  &*data.dataOnEntities[MBEDGE][side_number]
294  .getDiffN(base)
295  .data()
296  .begin(),
297  nb_gauss_pts);
298  }
299 
301 }
302 
305 
306  switch (cTx->bAse) {
310  break;
313  break;
314  default:
315  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
316  }
317 
319 }
320 
323 
324  EntitiesFieldData &data = cTx->dAta;
326 
327  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
328  double *diffL, const int dim) =
330 
331  if(base == AINSWORTH_LOBATTO_BASE)
332  base_polynomials = Lobatto_polynomials;
333 
334  int nb_gauss_pts = pts.size2();
335 
336  constexpr int side_number = 0;
337  int order = data.dataOnEntities[MBEDGE][side_number].getOrder();
338 
339  data.dataOnEntities[MBEDGE][side_number].getN(base).resize(
340  nb_gauss_pts, NBEDGE_L2(order), false);
341  data.dataOnEntities[MBEDGE][side_number].getDiffN(base).resize(
342  nb_gauss_pts, NBEDGE_L2(order), false);
343 
344  auto *fun_n =
345  &*data.dataOnEntities[MBEDGE][side_number].getN(base).data().begin();
346  auto *diff_fun_n =
347  &*data.dataOnEntities[MBEDGE][side_number].getDiffN(base).data().begin();
348 
349  if (NBEDGE_L2(order)) {
350 
351  double diff_mu = 2;
352  double l[NBEDGE_L2(order)];
353  double diff_l[NBEDGE_L2(order)];
354  for (int gg = 0; gg != nb_gauss_pts; gg++) {
355  double mu = 2 * pts(0, gg) - 1;
356  CHKERR base_polynomials(order, mu, &diff_mu, l, diff_l, 1);
357  int qd_shift = NBEDGE_L2(order) * gg;
358  for (int n = 0; n != NBEDGE_L2(order); n++) {
359  fun_n[qd_shift + n] = l[n];
360  diff_fun_n[qd_shift + n] = diff_l[n];
361  }
362  }
363 
364  }
365 
367 }
368 
370  return getValueL2AinsworthBase(pts);
371 }
372 
374  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
375  "Make no sense, unless problem is 2d (2d not implemented yet)");
376 }
377 
381 
382  EntitiesFieldData &data = cTx->dAta;
383  const FieldApproximationBase base = cTx->bAse;
384 
385  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s, double *L,
386  double *diffL, const int dim) =
388 
389  int nb_gauss_pts = pts.size2();
390  if (data.spacesOnEntities[MBEDGE].test(HCURL)) {
391 
392  if (data.dataOnEntities[MBEDGE].size() != 1)
393  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
394 
395  int sense = data.dataOnEntities[MBEDGE][0].getSense();
396  int order = data.dataOnEntities[MBEDGE][0].getOrder();
397  int nb_dofs =
398  NBEDGE_AINSWORTH_HCURL(data.dataOnEntities[MBEDGE][0].getOrder());
399  data.dataOnEntities[MBEDGE][0].getN(base).resize(nb_gauss_pts, 3 * nb_dofs,
400  false);
401  data.dataOnEntities[MBEDGE][0].getDiffN(base).resize(nb_gauss_pts, 0,
402  false);
403 
405  sense, order, &data.dataOnEntities[MBVERTEX][0].getN(base)(0, 0),
406  &*data.dataOnEntities[MBVERTEX][0].getDiffN(base).data().begin(),
407  &*data.dataOnEntities[MBEDGE][0].getN(base).data().begin(), nullptr,
408  nb_gauss_pts, base_polynomials);
409 
410  } else {
411  data.dataOnEntities[MBEDGE][0].getN(base).resize(nb_gauss_pts, 0, false);
412  data.dataOnEntities[MBEDGE][0].getDiffN(base).resize(nb_gauss_pts, 0,
413  false);
414  }
415 
417 }
418 
422 
423  EntitiesFieldData &data = cTx->dAta;
424  const FieldApproximationBase base = cTx->bAse;
425 
426  int nb_gauss_pts = pts.size2();
427  if (data.spacesOnEntities[MBEDGE].test(HCURL)) {
428 
429  if (data.dataOnEntities[MBEDGE].size() != 1)
430  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
431  "No data structure to store base functions");
432 
433  int sense = data.dataOnEntities[MBEDGE][0].getSense();
434  int order = data.dataOnEntities[MBEDGE][0].getOrder();
435  int nb_dofs =
436  NBEDGE_DEMKOWICZ_HCURL(data.dataOnEntities[MBEDGE][0].getOrder());
437  data.dataOnEntities[MBEDGE][0].getN(base).resize(nb_gauss_pts, 3 * nb_dofs,
438  false);
439  data.dataOnEntities[MBEDGE][0].getDiffN(base).resize(nb_gauss_pts, 0,
440  false);
441  if (nb_dofs) {
443  sense, order, &data.dataOnEntities[MBVERTEX][0].getN(base)(0, 0),
444  &*data.dataOnEntities[MBVERTEX][0].getDiffN(base).data().begin(),
445  &*data.dataOnEntities[MBEDGE][0].getN(base).data().begin(), NULL,
446  nb_gauss_pts);
447  }
448  } else {
449 
450  data.dataOnEntities[MBEDGE][0].getN(base).resize(nb_gauss_pts, 0, false);
451  data.dataOnEntities[MBEDGE][0].getDiffN(base).resize(nb_gauss_pts, 0,
452  false);
453  }
454 
456 }
457 
460 
461  switch (cTx->bAse) {
465  break;
468  break;
469  default:
470  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
471  }
472 
474 }
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
g
constexpr double g
Definition: shallow_wave.cpp:63
MoFEM::Tools::diffShapeFunMBEDGE
static constexpr std::array< double, 2 > diffShapeFunMBEDGE
Definition: Tools.hpp:68
MoFEM::Hcurl_Ainsworth_EdgeBaseFunctions_MBTET_ON_EDGE
MoFEMErrorCode Hcurl_Ainsworth_EdgeBaseFunctions_MBTET_ON_EDGE(int sense, int p, double *N, double *diffN, double *edgeN, double *diff_edgeN, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Edge based H-curl base functions on edge.
Definition: Hcurl.cpp:175
H1
@ H1
continuous field
Definition: definitions.h:85
NBEDGE_H1
#define NBEDGE_H1(P)
Numer of base function on edge for H1 space.
Definition: h1_hdiv_hcurl_l2.h:55
LASTBASE
@ LASTBASE
Definition: definitions.h:69
ShapeMBEDGE
PetscErrorCode ShapeMBEDGE(double *N, const double *G_X, int DIM)
Definition: fem_tools.c:761
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
MoFEM::EdgePolynomialBase::getValueL2
MoFEMErrorCode getValueL2(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:303
MoFEM::EdgePolynomialBase::getValueH1AinsworthBase
MoFEMErrorCode getValueH1AinsworthBase(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:196
NOBASE
@ NOBASE
Definition: definitions.h:59
ApproximationBaseNames
const static char *const ApproximationBaseNames[]
Definition: definitions.h:72
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::Hcurl_Demkowicz_EdgeBaseFunctions_MBEDGE
MoFEMErrorCode Hcurl_Demkowicz_EdgeBaseFunctions_MBEDGE(int sense, int p, double *n, double *diff_n, double *phi, double *diff_phi, int nb_integration_pts)
Edge based H-curl base functions on edge.
Definition: Hcurl.cpp:2144
MoFEM::EntPolynomialBaseCtx::dAta
EntitiesFieldData & dAta
Definition: EntPolynomialBaseCtx.hpp:36
MoFEM::EdgePolynomialBase::getValueH1
MoFEMErrorCode getValueH1(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:79
MoFEM::EntitiesFieldData::spacesOnEntities
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
Definition: EntitiesFieldData.hpp:49
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::EdgePolynomialBase::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: EdgePolynomialBase.cpp:8
Lobatto_polynomials
PetscErrorCode Lobatto_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Lobatto base functions .
Definition: base_functions.c:182
MoFEM::EdgePolynomialBase::L
VectorDouble L
Definition: EdgePolynomialBase.hpp:33
MoFEM::EntPolynomialBaseCtx::copyNodeBase
const FieldApproximationBase copyNodeBase
Definition: EntPolynomialBaseCtx.hpp:41
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::BernsteinBezier::baseFunctionsEdge
static MoFEMErrorCode baseFunctionsEdge(const int N, const int gdim, const int n_alpha, const int *alpha, const double *lambda, const double *grad_lambda, double *base, double *grad_base)
Definition: BernsteinBezier.cpp:419
MoFEM::EdgePolynomialBase::getValueH1DemkowiczBase
MoFEMErrorCode getValueH1DemkowiczBase(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:260
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
NBEDGE_L2
#define NBEDGE_L2(P)
Number of base functions on edge fro L2 space.
Definition: h1_hdiv_hcurl_l2.h:48
MoFEM::Types::MatrixInt
UBlasMatrix< int > MatrixInt
Definition: Types.hpp:76
MoFEM::EdgePolynomialBase::getValueH1BernsteinBezierBase
MoFEMErrorCode getValueH1BernsteinBezierBase(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:101
MoFEM::EntPolynomialBaseCtx::sPace
const FieldSpace sPace
Definition: EntPolynomialBaseCtx.hpp:37
MoFEM::EdgePolynomialBase::getValueL2DemkowiczBase
MoFEMErrorCode getValueL2DemkowiczBase(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:369
NBEDGE_DEMKOWICZ_HCURL
#define NBEDGE_DEMKOWICZ_HCURL(P)
Definition: h1_hdiv_hcurl_l2.h:108
MoFEM::EdgePolynomialBase::diffL
VectorDouble diffL
Definition: EdgePolynomialBase.hpp:33
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
MoFEM::EntPolynomialBaseCtx::bAse
const FieldApproximationBase bAse
Definition: EntPolynomialBaseCtx.hpp:39
MoFEM::EdgePolynomialBase::getValue
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Definition: EdgePolynomialBase.cpp:15
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
convert.n
n
Definition: convert.py:82
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
MoFEM::EdgePolynomialBase::cTx
EntPolynomialBaseCtx * cTx
Definition: EdgePolynomialBase.hpp:31
MoFEM::EdgePolynomialBase::getValueL2AinsworthBase
MoFEMErrorCode getValueL2AinsworthBase(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:321
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
MoFEM::EdgePolynomialBase
Calculate base functions on tetrahedral.
Definition: EdgePolynomialBase.hpp:19
NBEDGE_AINSWORTH_HCURL
#define NBEDGE_AINSWORTH_HCURL(P)
Definition: h1_hdiv_hcurl_l2.h:97
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::EdgePolynomialBase::getValueHcurl
MoFEMErrorCode getValueHcurl(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:458
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::EdgePolynomialBase::getValueHdiv
MoFEMErrorCode getValueHdiv(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:373
mu
double mu
Definition: dynamic_first_order_con_law.cpp:98
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::EntPolynomialBaseCtx::fieldName
const std::string fieldName
Definition: EntPolynomialBaseCtx.hpp:40
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
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
MoFEM::EdgePolynomialBase::getValueHcurlAinsworthBase
MoFEMErrorCode getValueHcurlAinsworthBase(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:379
MoFEM::EdgePolynomialBase::getValueHcurlDemkowiczBase
MoFEMErrorCode getValueHcurlDemkowiczBase(MatrixDouble &pts)
Definition: EdgePolynomialBase.cpp:420
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::BernsteinBezier::generateIndicesEdgeEdge
static MoFEMErrorCode generateIndicesEdgeEdge(const int N, int *alpha)
Definition: BernsteinBezier.cpp:85
MoFEM::DemkowiczHexAndQuad::H1_BubbleShapeFunctions_ONSEGMENT
MoFEMErrorCode H1_BubbleShapeFunctions_ONSEGMENT(int p, double *N, double *diffN, double *bubbleN, double *diff_bubbleN, int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:141
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
MoFEM::EntPolynomialBaseCtx::basePolynomialsType0
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Definition: EntPolynomialBaseCtx.hpp:27