Testing edge and bubble shape functions on a quad face.
static char help[] =
"...\n\n";
double s = 0;
for (
unsigned int ii = 0; ii <
m.size1(); ii++) {
for (
unsigned int jj = 0; jj <
m.size2(); jj++) {
s += std::abs(
m(ii, jj));
}
}
return s;
}
int main(
int argc,
char *argv[]) {
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
PetscBool flg = PETSC_TRUE;
#if PETSC_VERSION_GE(3, 6, 4)
255, &flg);
#else
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL,
"-my_file",
#endif
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, 1, "Error -my_file (mesh file needed)");
}
const char *option;
option = "";
CHKERR moab.get_entities_by_type(0, MBVERTEX, verts,
true);
CHKERR moab.get_coords(verts, &coords(0, 0));
for (
int i = 0;
i < verts.size(); ++
i) {
}
{
std::array<int, 4> faces_nodes = {0, 1, 2, 3};
constexpr double eps = 1e-4;
double quad_bubbles_sum = 3.275491e+01;
double quad_diff_bubbles_sum = 5.838131e+02;
faces_nodes.data(),
p, &*
N.data().begin(), &*diffN.data().begin(),
&*quad_bubbles.data().begin(), &*quad_diff_bubbles.data().begin(),
if (std::abs(quad_bubbles_sum -
sum_matrix(quad_bubbles)) >
eps) {
"wrong result = %8.6e",
sum_matrix(quad_bubbles));
}
if (std::abs(quad_diff_bubbles_sum -
sum_matrix(quad_diff_bubbles)) >
"wrong result = %8.6e",
sum_matrix(quad_diff_bubbles));
}
}
{
int sense[] = {1, -1, 1, -1};
int P[4];
int ee = 0;
for (; ee < 4; ee++) {
}
double *quad_edges_ptr[4];
std::array<MatrixDouble, 4> quad_edges;
for (auto ee : {0, 1, 2, 3}) {
quad_edges_ptr[ee] = &quad_edges[ee](0, 0);
}
double *quad_diff_edges_ptr[4];
std::array<MatrixDouble, 4> quad_diff_edges;
for (auto ee : {0, 1, 2, 3}) {
quad_diff_edges[ee] =
MatrixDouble(verts.size(), P[ee] * 2);
quad_diff_edges_ptr[ee] = &quad_diff_edges[ee](0, 0);
}
double quad_edges_sum[] = {1.237500e+01, 1.535050e+01, 1.781890e+01,
2.015678e+01};
double quad_diff_edges_sum[] = {8.470000e+01, 1.192510e+02, 1.580128e+02,
1.976378e+02};
sense,
p, &*
N.data().begin(), &*diffN.data().begin(), quad_edges_ptr,
for (int ee = 0; ee < 4; ++ee) {
if (std::abs(quad_edges_sum[ee] -
sum_matrix(quad_edges[ee])) >
eps)
"edge %d wrong result = %8.6e", ee,
if (std::abs(quad_diff_edges_sum[ee] -
"edge %d wrong result = %8.6e", ee,
}
}
}
}
#define CATCH_ERRORS
Catch errors.
@ MOFEM_ATOM_TEST_INVALID
#define CHKERR
Inline error check.
static double sum_matrix(MatrixDouble &m)
FTensor::Index< 'm', SPACE_DIM > m
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
#define NBFACEQUAD_H1(P)
Number of base functions on quad for H1 space.
PetscErrorCode H1_QuadShapeFunctions_MBQUAD(int *faces_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))
#define NBEDGE_H1(P)
Numer of base function on edge for H1 space.
PetscErrorCode H1_EdgeShapeFunctions_MBQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *diff_edgeN[4], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
FTensor::Index< 'i', SPACE_DIM > i
UBlasMatrix< double > MatrixDouble
implementation of Data Operators for Forces and Sources
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.