v0.10.0
prism_polynomial_approximation.cpp

Checking approximation functions on prism

/** \file prism_polynomial_approximation.cpp
\example prism_polynomial_approximation.cpp
\brief Checking approximation functions on prism
*/
/* 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 z) {
double r = 1;
for (int o = 1; o <= approx_order; ++o) {
for (int i = 0; i <= o; ++i) {
for (int j = 0; j <= (o - i); ++j) {
int k = o - i - j;
if (k >= 0) {
r += pow(x, i) * pow(y, j) * pow(z, k);
}
}
}
}
return r;
}
static inline VectorDouble3 diff_fun(double x, double y, double z) {
r.clear();
for (int o = 1; o <= approx_order; ++o) {
for (int i = 0; i <= o; ++i) {
for (int j = 0; j <= (o - i); ++j) {
int k = o - i - j;
if (k >= 0) {
r[0] += i > 0 ? i * pow(x, i - 1) * pow(y, j) * pow(z, k) : 0;
r[1] += j > 0 ? j * pow(x, i) * pow(y, j - 1) * pow(z, k) : 0;
r[2] += k > 0 ? k * pow(x, i) * pow(y, j) * pow(z, k - 1) : 0;
}
}
}
}
return r;
}
};
PrismOpCheck(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 PrismOpRhs
PrismOpRhs(SmartPetscObj<Vec> &f);
MoFEMErrorCode doWork(int side, EntityType type,
private:
SmartPetscObj<Vec> F;
};
struct PrismOpLhs
PrismOpLhs(SmartPetscObj<Mat> &a);
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
EntityType col_type,
private:
SmartPetscObj<Mat> A;
};
int getRuleTrianglesOnly(int order);
int getRuleThroughThickness(int order);
};
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
if (pcomm == NULL)
pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
0, 0, 1, 1, 0, 1, 0, 1, 1};
std::array<EntityHandle, 6> one_prism_nodes;
for (int n = 0; n != 6; ++n)
CHKERR moab.create_vertex(&one_prism_coords[3 * n], one_prism_nodes[n]);
EntityHandle one_prism;
CHKERR moab.create_element(MBPRISM, one_prism_nodes.data(), 6, one_prism);
Range one_prism_range;
one_prism_range.insert(one_prism);
Range one_prism_adj_ents;
for (int d = 1; d != 3; ++d)
CHKERR moab.get_adjacencies(one_prism_range, d, true, one_prism_adj_ents,
moab::Interface::UNION);
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
PrismsFromSurfaceInterface *prisms_from_surface_interface;
CHKERR m_field.getInterface(prisms_from_surface_interface);
BitRefLevel bit_level0;
bit_level0.set(0);
CHKERR m_field.getInterface<BitRefManager>()->setEntitiesBitRefLevel(
one_prism_range, bit_level0);
CHKERR prisms_from_surface_interface->seedPrismsEntities(one_prism_range,
bit_level0);
// Fields
CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
CHKERR m_field.add_ents_to_field_by_type(0, MBPRISM, "FIELD1");
CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD1", 1);
CHKERR m_field.set_field_order(0, MBEDGE, "FIELD1", approx_order);
CHKERR m_field.set_field_order(0, MBTRI, "FIELD1", approx_order);
CHKERR m_field.set_field_order(0, MBQUAD, "FIELD1", approx_order);
CHKERR m_field.set_field_order(0, MBPRISM, "FIELD1", approx_order);
CHKERR m_field.build_fields();
// FE
CHKERR m_field.add_finite_element("PRISM");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("PRISM", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("PRISM", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("PRISM", "FIELD1");
CHKERR m_field.add_ents_to_finite_element_by_type(0, MBPRISM, "PRISM");
// 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", "PRISM");
// 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 assemble_matrices_and_vectors = [&]() {
PrismFE fe(m_field);
MatrixDouble inv_jac;
fe.getOpPtrVector().push_back(
new OpMultiplyDeterminantOfJacobianAndWeightsForFatPrisms());
fe.getOpPtrVector().push_back(
fe.getOpPtrVector().push_back(
fe.getOpPtrVector().push_back(new PrismOpRhs(F));
fe.getOpPtrVector().push_back(new PrismOpLhs(A));
CHKERR VecZeroEntries(F);
CHKERR MatZeroEntries(A);
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", 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 = [&] {
PrismFE fe(m_field);
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 OpCalculateInvJacForFatPrism(inv_jac));
fe.getOpPtrVector().push_back(
new OpSetInvJacH1ForFatPrism(inv_jac));
fe.getOpPtrVector().push_back(
new OpCalculateScalarFieldValues("FIELD1", field_vals_ptr));
fe.getOpPtrVector().push_back(
new OpCalculateScalarFieldGradient<3>("FIELD1", diff_field_vals_ptr));
fe.getOpPtrVector().push_back(
new PrismOpCheck(field_vals_ptr, diff_field_vals_ptr));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", fe);
};
CHKERR assemble_matrices_and_vectors();
CHKERR solve_problem();
CHKERR check_solution();
}
return 0;
}
PrismOpCheck::PrismOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
boost::shared_ptr<MatrixDouble> &diff_field_vals)
: FatPrismElementForcesAndSourcesCore::UserDataOperator(
"FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
fieldVals(field_vals), diffFieldVals(diff_field_vals) {}
if (type == MBVERTEX) {
const int nb_gauss_pts = data.getN().size2();
auto t_coords = getFTensor1CoordsAtGaussPts();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
double f = ApproxFunction::fun(t_coords(0), t_coords(1), t_coords(2));
VectorDouble3 diff_f =
ApproxFunction::diff_fun(t_coords(0), t_coords(1), t_coords(2));
std::cout << f - (*fieldVals)[gg] << " : ";
for (auto d : {0, 1, 2})
std::cout << diff_f[d] - (*diffFieldVals)(d, gg) << " ";
std::cout << std::endl;
constexpr double eps = 1e-6;
if (std::abs(f - (*fieldVals)[gg]) > eps ||
!std::isnormal((*fieldVals)[gg]))
SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong value %6.4e != %6.4e (%6.4e)", f, (*fieldVals)[gg],
f - (*fieldVals)[gg]);
for (auto d : {0, 1, 2})
if (std::abs(diff_f[d] - (*diffFieldVals)(d, gg)) > eps ||
!std::isnormal((*diffFieldVals)(d, gg)))
SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong diff value %6.4e != %6.4e (%6.4e)", diff_f[d],
(*diffFieldVals)(d, gg),
diff_f[d] - (*diffFieldVals)(d, gg));
++t_coords;
}
}
}
PrismOpRhs::PrismOpRhs(SmartPetscObj<Vec> &f)
: FatPrismElementForcesAndSourcesCore::UserDataOperator(
"FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
F(f) {}
MoFEMErrorCode PrismOpRhs::doWork(int side, EntityType type,
const int nb_dofs = data.getN().size2();
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();
double vol = getVolume();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
double f = ApproxFunction::fun(t_coords(0), t_coords(1), t_coords(2));
double v = t_w * vol * f;
double *val = &*nf.begin();
for (int bb = 0; bb != nb_dofs; ++bb) {
*val += v * t_base;
++t_base;
++val;
}
++t_coords;
++t_w;
}
CHKERR VecSetValues(F, data, &*nf.data().begin(), ADD_VALUES);
}
}
PrismOpLhs::PrismOpLhs(SmartPetscObj<Mat> &a)
: FatPrismElementForcesAndSourcesCore::UserDataOperator(
"FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROWCOL),
A(a) {
// FIXME: Can be symmetric, is not for simplicity
sYmm = false;
}
MoFEMErrorCode PrismOpLhs::doWork(int row_side, int col_side,
EntityType row_type, EntityType col_type,
const int row_nb_dofs = row_data.getN().size2();
const int col_nb_dofs = col_data.getN().size2();
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();
double vol = getVolume();
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 = t_w * vol;
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);
}
}
int PrismFE::getRuleTrianglesOnly(int order) { return 2 * (order + 1); };
int PrismFE::getRuleThroughThickness(int order) { return 2 * (order + 1); };
static Index< 'm', 3 > m
static Index< 'n', 3 > n
static Index< 'o', 3 > o
@ 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
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
@ H1
continuous field
Definition: definitions.h:177
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:292
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
@ 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.
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
FTensor::Index< 'k', 3 > k
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
CoreTmp< 0 > Core
Definition: Core.hpp:1129
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
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.
Calculate inverse of jacobian for face element.
Transform local reference derivatives of shape functions to global derivatives.
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
int getRuleThroughThickness(int order)
int getRuleTrianglesOnly(int order)
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > diffFieldVals
boost::shared_ptr< VectorDouble > fieldVals
PrismOpCheck(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.
SmartPetscObj< Mat > A
PrismOpLhs(SmartPetscObj< Mat > &a)
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
SmartPetscObj< Vec > F
PrismOpRhs(SmartPetscObj< Vec > &f)