v0.10.0
quad_polynomial_approximation.cpp

Checking approximation functions for quad

/** \file quad_polynomial_approximation.cpp
\example quad_polynomial_approximation.cpp
\brief Checking approximation functions for quad
*/
/* This file is part of MoFEM.
* MoFEM is free software: you can redistribute it and/or modify it under
* the terms of the GNU Lesser General Public License as published by the
* Free Software Foundation, either version 3 of the License, or (at your
* option) any later version.
*
* MoFEM is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
* License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
static constexpr int approx_order = 6;
static inline double fun(double x, double y) {
double r = 1;
for (int o = 1; o <= approx_order; ++o) {
for (int i = 0; i <= o; ++i) {
int j = o - i;
if (j >= 0)
r += pow(x, i) * pow(y, j);
}
}
return r;
}
static inline VectorDouble3 diff_fun(double x, double y) {
r.clear();
for (int o = 1; o <= approx_order; ++o) {
for (int i = 0; i <= o; ++i) {
int j = o - i;
if (j >= 0) {
r[0] += i > 0 ? i * pow(x, i - 1) * pow(y, j) : 0;
r[1] += j > 0 ? j * pow(x, i) * pow(y, j - 1) : 0;
}
}
}
return r;
}
};
struct QuadOpCheck : public OpEle {
QuadOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
boost::shared_ptr<MatrixDouble> &diff_field_vals);
MoFEMErrorCode doWork(int side, EntityType type,
private:
boost::shared_ptr<VectorDouble> fieldVals;
boost::shared_ptr<MatrixDouble> diffFieldVals;
};
struct QuadOpRhs : public OpEle {
QuadOpRhs(SmartPetscObj<Vec> &f);
MoFEMErrorCode doWork(int side, EntityType type,
private:
SmartPetscObj<Vec> F;
};
struct QuadOpLhs : public OpEle {
QuadOpLhs(SmartPetscObj<Mat> &a);
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
EntityType col_type,
private:
SmartPetscObj<Mat> A;
};
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
// Declare elements
enum bases {
AINSWORTH,
AINSWORTH_LOBATTO,
DEMKOWICZ,
BERNSTEIN,
LASBASETOP
};
const char *list_bases[] = {"ainsworth", "ainsworth_labatto", "demkowicz",
"bernstein"};
PetscBool flg;
PetscInt choice_base_value = AINSWORTH;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
LASBASETOP, &choice_base_value, &flg);
if (flg != PETSC_TRUE)
SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSIBLE_CASE, "base not set");
if (choice_base_value == AINSWORTH)
if (choice_base_value == AINSWORTH_LOBATTO)
else if (choice_base_value == DEMKOWICZ)
else if (choice_base_value == BERNSTEIN)
enum spaces { H1SPACE, L2SPACE, LASBASETSPACE };
const char *list_spaces[] = {"h1", "l2"};
PetscInt choice_space_value = H1SPACE;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
LASBASETSPACE, &choice_space_value, &flg);
if (flg != PETSC_TRUE)
SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSIBLE_CASE, "space not set");
FieldSpace space = H1;
if (choice_space_value == H1SPACE)
space = H1;
else if (choice_space_value == L2SPACE)
space = L2;
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
std::array<double, 12> one_quad_coords = {0, 0, 0,
2, 0, 0,
1, 1, 0,
0, 1, 0};
std::array<EntityHandle, 4> one_quad_nodes;
for (int n = 0; n != 4; ++n)
CHKERR moab.create_vertex(&one_quad_coords[3 * n], one_quad_nodes[n]);
EntityHandle one_quad;
CHKERR moab.create_element(MBQUAD, one_quad_nodes.data(), 4, one_quad);
Range one_quad_range;
one_quad_range.insert(one_quad);
Range one_quad_adj_ents;
CHKERR moab.get_adjacencies(one_quad_range, 1, true, one_quad_adj_ents,
moab::Interface::UNION);
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
BitRefLevel bit_level0 = BitRefLevel().set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 2, bit_level0);
// Fields
CHKERR m_field.add_field("FIELD1", space, base, 1);
CHKERR m_field.add_ents_to_field_by_type(0, MBQUAD, "FIELD1");
CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD1", 1);
CHKERR m_field.set_field_order(0, MBEDGE, "FIELD1", approx_order + 1);
CHKERR m_field.set_field_order(0, MBQUAD, "FIELD1", approx_order + 1);
CHKERR m_field.build_fields();
// FE
CHKERR m_field.add_finite_element("QUAD");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("QUAD", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("QUAD", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("QUAD", "FIELD1");
CHKERR m_field.add_ents_to_finite_element_by_type(0, MBQUAD, "QUAD");
// build finite elemnts
// //build adjacencies
CHKERR m_field.build_adjacencies(bit_level0);
// Problem
CHKERR m_field.add_problem("TEST_PROBLEM");
// set finite elements for problem
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "QUAD");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
// build problem
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
// partition
CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
// what are ghost nodes, see Petsc Manual
CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
// Create matrices
SmartPetscObj<Mat> A;
CHKERR m_field.getInterface<MatrixManager>()
->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("TEST_PROBLEM", A);
SmartPetscObj<Vec> F;
CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
ROW, F);
SmartPetscObj<Vec> D;
CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
COL, D);
auto rule = [&](int, int, int p) { return 2 * (p + 1); };
auto assemble_matrices_and_vectors = [&]() {
Ele fe(m_field);
fe.getRuleHook = rule;
MatrixDouble inv_jac;
fe.getOpPtrVector().push_back(new OpCalculateInvJacForFace(inv_jac));
fe.getOpPtrVector().push_back(new OpSetInvJacH1ForFace(inv_jac));
fe.getOpPtrVector().push_back(new OpSetInvJacL2ForFace(inv_jac));
fe.getOpPtrVector().push_back(new OpMakeHighOrderGeometryWeightsOnFace());
fe.getOpPtrVector().push_back(new QuadOpRhs(F));
fe.getOpPtrVector().push_back(new QuadOpLhs(A));
CHKERR VecZeroEntries(F);
CHKERR MatZeroEntries(A);
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "QUAD", fe);
CHKERR VecAssemblyBegin(F);
CHKERR VecAssemblyEnd(F);
CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
};
auto solve_problem = [&] {
auto solver = createKSP(PETSC_COMM_WORLD);
CHKERR KSPSetOperators(solver, A, A);
CHKERR KSPSetFromOptions(solver);
CHKERR KSPSetUp(solver);
CHKERR KSPSolve(solver, F, D);
CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR m_field.getInterface<VecManager>()->setLocalGhostVector(
"TEST_PROBLEM", COL, D, INSERT_VALUES, SCATTER_REVERSE);
};
auto check_solution = [&] {
Ele fe(m_field);
fe.getRuleHook = rule;
boost::shared_ptr<VectorDouble> field_vals_ptr(new VectorDouble());
boost::shared_ptr<MatrixDouble> diff_field_vals_ptr(new MatrixDouble());
MatrixDouble inv_jac;
fe.getOpPtrVector().push_back(
new OpCalculateScalarFieldValues("FIELD1", field_vals_ptr));
fe.getOpPtrVector().push_back(new OpCalculateInvJacForFace(inv_jac));
fe.getOpPtrVector().push_back(new OpSetInvJacH1ForFace(inv_jac));
fe.getOpPtrVector().push_back(new OpSetInvJacL2ForFace(inv_jac));
fe.getOpPtrVector().push_back(new OpMakeHighOrderGeometryWeightsOnFace());
fe.getOpPtrVector().push_back(new OpCalculateScalarFieldGradient<2>(
"FIELD1", diff_field_vals_ptr, space == L2 ? MBQUAD : MBVERTEX));
fe.getOpPtrVector().push_back(
new QuadOpCheck(field_vals_ptr, diff_field_vals_ptr));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "QUAD", fe);
};
CHKERR assemble_matrices_and_vectors();
CHKERR solve_problem();
CHKERR check_solution();
}
return 0;
}
QuadOpCheck::QuadOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
boost::shared_ptr<MatrixDouble> &diff_field_vals)
: OpEle("FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
fieldVals(field_vals), diffFieldVals(diff_field_vals) {}
if (type == MBQUAD) {
const int nb_gauss_pts = data.getN().size1();
auto t_coords = getFTensor1CoordsAtGaussPts();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
double f = ApproxFunction::fun(t_coords(0), t_coords(1));
constexpr double eps = 1e-6;
std::cout << f - (*fieldVals)[gg] << std::endl;
if (std::abs(f - (*fieldVals)[gg]) > eps)
SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong value %d : %6.4e != %6.4e", gg, f, (*fieldVals)[gg]);
VectorDouble3 diff_f = ApproxFunction::diff_fun(t_coords(0), t_coords(1));
for (auto d : {0, 1})
std::cout << diff_f[d] - (*diffFieldVals)(d, gg) << " ";
std::cout << std::endl;
for (auto d : {0, 1})
if (std::abs(diff_f[d] - (*diffFieldVals)(d, gg)) > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong direvative value (%d) %6.4e != %6.4e", diff_f[d],
(*diffFieldVals)(d, gg));
++t_coords;
}
}
}
QuadOpRhs::QuadOpRhs(SmartPetscObj<Vec> &f)
: OpEle("FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
F(f) {}
MoFEMErrorCode QuadOpRhs::doWork(int side, EntityType type,
const int nb_dofs = data.getIndices().size();
if (nb_dofs) {
const int nb_gauss_pts = data.getN().size1();
VectorDouble nf(nb_dofs);
nf.clear();
auto t_base = data.getFTensor0N();
auto t_coords = getFTensor1CoordsAtGaussPts();
auto a = getMeasure();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
double f = ApproxFunction::fun(t_coords(0), t_coords(1));
double v = a * t_w * f;
double *val = &*nf.begin();
for (int bb = 0; bb != nb_dofs; ++bb) {
*val += v * t_base;
++t_base;
++val;
}
++t_coords;
++t_w;
// ++t_normal;
}
CHKERR VecSetValues(F, data, &*nf.data().begin(), ADD_VALUES);
}
}
QuadOpLhs::QuadOpLhs(SmartPetscObj<Mat> &a)
: OpEle("FIELD1", "FIELD1",
ForcesAndSourcesCore::UserDataOperator::OPROWCOL),
A(a) {
// FIXME: Can be symmetric, is not for simplicity
sYmm = false;
}
MoFEMErrorCode QuadOpLhs::doWork(int row_side, int col_side,
EntityType row_type, EntityType col_type,
const int row_nb_dofs = row_data.getIndices().size();
const int col_nb_dofs = col_data.getIndices().size();
if (row_nb_dofs && col_nb_dofs) {
const int nb_gauss_pts = row_data.getN().size1();
MatrixDouble m(row_nb_dofs, col_nb_dofs);
m.clear();
auto a = getMeasure();
double *row_base_ptr = &*row_data.getN().data().begin();
double *col_base_ptr = &*col_data.getN().data().begin();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
double v = a * t_w;
cblas_dger(CblasRowMajor, row_nb_dofs, col_nb_dofs, v, row_base_ptr, 1,
col_base_ptr, 1, &*m.data().begin(), col_nb_dofs);
row_base_ptr += row_nb_dofs;
col_base_ptr += col_nb_dofs;
++t_w;
}
CHKERR MatSetValues(A, row_data, col_data, &*m.data().begin(), ADD_VALUES);
}
}
static Index< 'm', 3 > m
static Index< 'n', 3 > n
static Index< 'o', 3 > o
static Index< 'p', 3 > p
ForcesAndSourcesCore::UserDataOperator UserDataOperator
@ CblasRowMajor
Definition: cblas.h:10
void cblas_dger(const enum CBLAS_ORDER order, const int M, const int N, const double alpha, const double *X, const int incX, const double *Y, const int incY, double *A, const int lda)
Definition: cblas_dger.c:12
@ COL
Definition: definitions.h:192
@ ROW
Definition: definitions.h:192
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:441
FieldApproximationBase
approximation base
Definition: definitions.h:150
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:154
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:158
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:156
FieldSpace
approximation spaces
Definition: definitions.h:174
@ L2
field with C-1 continuity
Definition: definitions.h:180
@ H1
continuous field
Definition: definitions.h:177
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
@ MOFEM_IMPOSIBLE_CASE
Definition: definitions.h:127
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:132
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define CHKERR
Inline error check.
Definition: definitions.h:604
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string &name_row)=0
set field row which finite element use
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string &name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string &name_row)=0
set field col which finite element use
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
FaceElementForcesAndSourcesCoreSwitch< 0 > FaceElementForcesAndSourcesCore
Face finite element default.
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string &name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
FTensor::Index< 'j', 3 > j
FTensor::Index< 'i', 3 > i
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:74
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:76
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
auto createKSP
OpCalculateInvJacForFaceImpl< 2 > OpCalculateInvJacForFace
CoreTmp< 0 > Core
Definition: Core.hpp:1129
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1943
constexpr auto VecSetValues
constexpr auto MatSetValues
const double D
diffusivity
const double r
rate factor
int main(int argc, char *argv[])
static char help[]
static constexpr int approx_order
DataForcesAndSourcesCore::EntData EntData
static VectorDouble3 diff_fun(double x, double y, double z)
static double fun(double x, double y, double z)
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
Core (interface) class.
Definition: Core.hpp:77
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:60
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:100
Data on single entity (This is passed as argument to DataOperator::doWork)
Deprecated interface functions.
DEPRECATED MoFEMErrorCode loop_finite_elements(const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, MoFEMTypes bh, int verb=DEFAULT_VERBOSITY)
auto getFTensor0IntegrationWeight()
Get integration weights.
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
boost::shared_ptr< VectorDouble > fieldVals
boost::shared_ptr< MatrixDouble > diffFieldVals
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
QuadOpCheck(boost::shared_ptr< VectorDouble > &field_vals, boost::shared_ptr< MatrixDouble > &diff_field_vals)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
QuadOpLhs(SmartPetscObj< Mat > &a)
SmartPetscObj< Mat > A
SmartPetscObj< Vec > F
QuadOpRhs(SmartPetscObj< Vec > &f)
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Operator for linear form, usually to calculate values on right hand side.