v0.14.0
Loading...
Searching...
No Matches
prism_polynomial_approximation.cpp

Checking approximation functions on prism.

Checking approximation functions on prism

/** \file prism_polynomial_approximation.cpp
\example prism_polynomial_approximation.cpp
\brief Checking approximation functions on prism
*/
#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);
private:
boost::shared_ptr<VectorDouble> fieldVals;
boost::shared_ptr<MatrixDouble> diffFieldVals;
};
struct PrismOpRhs
private:
};
struct PrismOpLhs
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
EntityType col_type,
private:
};
};
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);
auto moab_comm_wrap =
boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
if (pcomm == NULL)
pcomm = new ParallelComm(&moab, moab_comm_wrap->get_comm());
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
->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("TEST_PROBLEM", A);
CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
ROW, F);
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(
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(
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)
"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;
}
}
}
"FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
F(f) {}
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);
}
}
"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< 'o', 3 > o
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static char help[]
int main()
Definition: adol-c_atom.cpp:46
constexpr double a
static const double eps
@ COL
Definition: definitions.h:123
@ ROW
Definition: definitions.h:123
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
auto fun
Function to approximate.
constexpr int order
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'm', SPACE_DIM > m
@ F
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 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 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_row(const std::string &fe_name, const std::string name_row)=0
set field row 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 loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
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
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
FTensor::Index< 'i', SPACE_DIM > i
auto diff_fun
Function derivative.
double D
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'j', 3 > j
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:56
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
auto createKSP(MPI_Comm comm)
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
constexpr AssemblyType A
static constexpr int approx_order
static constexpr int approx_order
static VectorDouble3 diff_fun(double x, double y, double z)
static double fun(double x, double y, double z)
Managing BitRefLevels.
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:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
structure to get information form mofem into EntitiesFieldData
Matrix manager is used to build and partition problems.
Calculate inverse of jacobian for face element.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Get value at integration points for scalar field.
Operator for fat prism element updating integration weights in the volume.
Transform local reference derivatives of shape functions to global derivatives.
MoFEMErrorCode seedPrismsEntities(Range &prisms, const BitRefLevel &bit, int verb=-1)
Seed prism entities by bit level.
Problem manager is used to build and partition problems.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
int getRuleThroughThickness(int order)
int getRuleTrianglesOnly(int order)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::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)
SmartPetscObj< Mat > A
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
PrismOpLhs(SmartPetscObj< Mat > &a)
SmartPetscObj< Vec > F
PrismOpRhs(SmartPetscObj< Vec > &f)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.