v0.14.0
scalar_check_approximation.cpp

Approximate polynomial in 2D and check derivatives

/**
* \file scalar_check_approximation.cpp
* \example scalar_check_approximation.cpp
*
* Approximate polynomial in 2D and check derivatives
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
static int approx_order = 4;
template <int DIM> struct ApproxFunctionsImpl {};
template <int DIM> struct ElementsAndOps {};
constexpr int SPACE_DIM =
EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
using BoundaryEle =
template <> struct ApproxFunctionsImpl<2> {
static double fUn(const double x, const double y, double z) {
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 FTensor::Tensor1<double, 2> diffFun(const double x, const double y,
double z) {
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;
}
};
template <> struct ApproxFunctionsImpl<3> {
static double fUn(const double x, const 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 FTensor::Tensor1<double, 3> diffFun(const double x, const double y,
double z) {
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;
}
};
struct OpValsDiffVals : public DomainEleOp {
VectorDouble &vAls;
MatrixDouble &diffVals;
const bool checkGradients;
OpValsDiffVals(VectorDouble &vals, MatrixDouble &diff_vals, bool check_grads)
: DomainEleOp("FIELD1", OPROW), vAls(vals), diffVals(diff_vals),
checkGradients(check_grads) {}
MoFEMErrorCode doWork(int side, EntityType type,
const int nb_gauss_pts = getGaussPts().size2();
if (type == MBVERTEX) {
vAls.resize(nb_gauss_pts, false);
diffVals.resize(SPACE_DIM, nb_gauss_pts, false);
vAls.clear();
diffVals.clear();
}
const int nb_dofs = data.getIndices().size();
if (nb_dofs) {
MOFEM_LOG("AT", Sev::noisy)
<< "Type " << moab::CN::EntityTypeName(type) << " side " << side;
MOFEM_LOG("AT", Sev::noisy) << data.getN();
MOFEM_LOG("AT", Sev::noisy) << data.getDiffN();
auto t_vals = getFTensor0FromVec(vAls);
auto t_base_fun = data.getFTensor0N();
for (int gg = 0; gg != nb_gauss_pts; gg++) {
auto t_data = data.getFTensor0FieldData();
for (int bb = 0; bb != nb_dofs; bb++) {
t_vals += t_base_fun * t_data;
++t_base_fun;
++t_data;
}
const double v = t_vals;
if (!std::isnormal(v))
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Not a number");
++t_vals;
}
if (checkGradients) {
auto t_diff_vals = getFTensor1FromMat<SPACE_DIM>(diffVals);
auto t_diff_base_fun = data.getFTensor1DiffN<SPACE_DIM>();
for (int gg = 0; gg != nb_gauss_pts; gg++) {
auto t_data = data.getFTensor0FieldData();
for (int bb = 0; bb != nb_dofs; bb++) {
t_diff_vals(i) += t_diff_base_fun(i) * t_data;
++t_diff_base_fun;
++t_data;
}
for (int d = 0; d != SPACE_DIM; ++d)
if (!std::isnormal(t_diff_vals(d)))
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Not a number");
++t_diff_vals;
}
}
}
}
};
VectorDouble &vAls;
MatrixDouble &diffVals;
boost::shared_ptr<VectorDouble> ptrVals;
boost::shared_ptr<MatrixDouble> ptrDiffVals;
const bool checkGradients;
boost::shared_ptr<VectorDouble> &ptr_vals,
boost::shared_ptr<MatrixDouble> &ptr_diff_vals,
bool check_grads)
: DomainEleOp("FIELD1", OPROW), vAls(vals), diffVals(diff_vals),
ptrVals(ptr_vals), ptrDiffVals(ptr_diff_vals),
checkGradients(check_grads) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type,
const double eps = 1e-6;
const int nb_gauss_pts = data.getN().size1();
auto t_vals = getFTensor0FromVec(vAls);
auto t_ptr_vals = getFTensor0FromVec(*ptrVals);
for (int gg = 0; gg != nb_gauss_pts; gg++) {
double err_val;
// Check user data operators
err_val = std::abs(t_vals - t_ptr_vals);
MOFEM_LOG("AT", Sev::noisy) << "Val op error " << err_val;
if (err_val > eps)
SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong value from operator %4.3e", err_val);
const double x = getCoordsAtGaussPts()(gg, 0);
const double y = getCoordsAtGaussPts()(gg, 1);
const double z = getCoordsAtGaussPts()(gg, 2);
// Check approximation
const double delta_val = t_vals - ApproxFunctions::fUn(x, y, z);
err_val = std::fabs(delta_val * delta_val);
MOFEM_LOG("AT", Sev::verbose) << err_val << " : " << t_vals;
if (err_val > eps)
SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Wrong value %4.3e",
err_val);
++t_vals;
++t_ptr_vals;
}
if (checkGradients) {
auto t_diff_vals = getFTensor1FromMat<SPACE_DIM>(diffVals);
auto t_ptr_diff_vals = getFTensor1FromMat<SPACE_DIM>(*ptrDiffVals);
for (int gg = 0; gg != nb_gauss_pts; gg++) {
double err_diff_val;
t_delta_diff_val(i) = t_diff_vals(i) - t_ptr_diff_vals(i);
err_diff_val = sqrt(t_delta_diff_val(i) * t_delta_diff_val(i));
MOFEM_LOG("AT", Sev::noisy) << "Diff val op error " << err_diff_val;
if (err_diff_val > eps)
SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong derivatives from operator %4.3e", err_diff_val);
const double x = getCoordsAtGaussPts()(gg, 0);
const double y = getCoordsAtGaussPts()(gg, 1);
const double z = getCoordsAtGaussPts()(gg, 2);
// Check approximation
auto t_diff_anal = ApproxFunctions::diffFun(x, y, z);
t_delta_diff_val(i) = t_diff_vals(i) - t_diff_anal(i);
err_diff_val = sqrt(t_delta_diff_val(i) * t_delta_diff_val(i));
if (SPACE_DIM == 3)
MOFEM_LOG("AT", Sev::noisy)
<< "Diff val " << err_diff_val << " : "
<< sqrt(t_diff_vals(i) * t_diff_vals(i)) << " : "
<< t_diff_vals(0) << " (" << t_diff_anal(0) << ") "
<< t_diff_vals(1) << " (" << t_diff_anal(1) << ") "
<< t_diff_vals(2) << " (" << t_diff_anal(2) << ")";
else
MOFEM_LOG("AT", Sev::noisy)
<< "Diff val " << err_diff_val << " : "
<< sqrt(t_diff_vals(i) * t_diff_vals(i)) << " : "
<< t_diff_vals(0) << " (" << t_diff_anal(0) << ") "
<< t_diff_vals(1) << " (" << t_diff_anal(1) << ")";
MOFEM_LOG("AT", Sev::verbose)
<< getCoords()(3 * 1 + 0) - getCoords()(3 * 0 + 0);
MOFEM_LOG("AT", Sev::verbose)
<< getCoords()(3 * 1 + 1) - getCoords()(3 * 0 + 1);
MOFEM_LOG("AT", Sev::verbose)
<< getCoords()(3 * 1 + 2) - getCoords()(3 * 0 + 2);
MOFEM_LOG("AT", Sev::verbose) << "Diff val error " << err_diff_val;
if (err_diff_val > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong derivative of value %4.3e %4.3e", err_diff_val,
t_diff_anal.l2());
++t_diff_vals;
++t_ptr_diff_vals;
}
}
}
};
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
DMType dm_name = "DMMOFEM";
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
// Add logging channel for example
auto core_log = logging::core::get();
core_log->add_sink(
MOFEM_LOG_TAG("AT", "atom_test");
// Create MoFEM instance
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
PipelineManager *pipeline_mng = m_field.getInterface<PipelineManager>();
CHKERR simple->getOptions();
simple->getAddBoundaryFE() = true;
CHKERR simple->loadFile("", "");
// Declare elements
enum bases {
AINSWORTH,
AINSWORTH_LOBATTO,
DEMKOWICZ,
BERNSTEIN,
LASBASETOP
};
const char *list_bases[] = {"ainsworth", "ainsworth_lobatto", "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_IMPOSSIBLE_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_IMPOSSIBLE_CASE, "space not set");
FieldSpace space = H1;
if (choice_space_value == H1SPACE)
space = H1;
else if (choice_space_value == L2SPACE)
space = L2;
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &approx_order,
PETSC_NULL);
CHKERR simple->addDomainField("FIELD1", space, base, 1);
CHKERR simple->setFieldOrder("FIELD1", approx_order);
Range edges, faces;
CHKERR moab.get_entities_by_dimension(0, 1, edges);
CHKERR moab.get_entities_by_dimension(0, 2, faces);
if (choice_base_value != BERNSTEIN) {
Range rise_order;
int i = 0;
for (auto e : edges) {
if (!(i % 2)) {
rise_order.insert(e);
}
++i;
}
for (auto f : faces) {
if (!(i % 3)) {
rise_order.insert(f);
}
++i;
}
Range rise_twice;
for (auto e : rise_order) {
if (!(i % 2)) {
rise_twice.insert(e);
}
++i;
}
CHKERR simple->setFieldOrder("FIELD1", approx_order + 1, &rise_order);
CHKERR simple->setFieldOrder("FIELD1", approx_order + 2, &rise_twice);
}
CHKERR simple->defineFiniteElements();
auto volume_adj = [](moab::Interface &moab, const Field &field,
const EntFiniteElement &fe,
std::vector<EntityHandle> &adjacency) {
EntityHandle fe_ent = fe.getEnt();
switch (field.getSpace()) {
case H1:
CHKERR moab.get_connectivity(&fe_ent, 1, adjacency, true);
case HCURL:
CHKERR moab.get_adjacencies(&fe_ent, 1, 1, false, adjacency,
moab::Interface::UNION);
case HDIV:
case L2:
CHKERR moab.get_adjacencies(&fe_ent, 1, 2, false, adjacency,
moab::Interface::UNION);
adjacency.push_back(fe_ent);
// build side table
for (auto ent : adjacency)
fe.getSideNumberPtr(ent);
break;
case NOFIELD: {
CHKERR moab.get_entities_by_handle(field.getMeshset(), adjacency,
false);
for (auto ent : adjacency) {
const_cast<SideNumber_multiIndex &>(fe.getSideNumberTable())
.insert(
boost::shared_ptr<SideNumber>(new SideNumber(ent, -1, 0, 0)));
}
} break;
default:
SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
"this field is not implemented for TRI finite element");
}
};
simple->getDomainFEName(), MBTET, volume_adj);
simple->getDomainFEName(), MBHEX, volume_adj);
CHKERR simple->defineProblem(PETSC_TRUE);
CHKERR simple->buildFields();
CHKERR simple->buildFiniteElements();
CHKERR simple->buildProblem();
auto dm = simple->getDM();
MatrixDouble diff_vals;
auto assemble_matrices_and_vectors = [&]() {
PETSC>::LinearForm<GAUSS>::OpSource<1, 1>;
pipeline_mng->getOpDomainLhsPipeline(), {space});
pipeline_mng->getOpDomainRhsPipeline(), {space});
pipeline_mng->getOpDomainLhsPipeline().push_back(
pipeline_mng->getOpDomainLhsPipeline().push_back(new OpMass(
"FIELD1", "FIELD1", [](double, double, double) { return 1.; }));
pipeline_mng->getOpDomainLhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpSource("FIELD1", ApproxFunctions::fUn));
auto integration_rule = [](int, int, int p_data) {
return 2 * p_data + 1;
};
CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
};
auto solve_problem = [&] {
auto solver = pipeline_mng->createKSP();
CHKERR KSPSetFromOptions(solver);
CHKERR KSPSetUp(solver);
auto dm = simple->getDM();
auto D = createDMVector(dm);
auto F = vectorDuplicate(D);
CHKERR KSPSolve(solver, F, D);
CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
};
auto check_solution = [&] {
auto ptr_values = boost::make_shared<VectorDouble>();
auto ptr_diff_vals = boost::make_shared<MatrixDouble>();
pipeline_mng->getDomainLhsFE().reset();
pipeline_mng->getOpDomainRhsPipeline().clear();
pipeline_mng->getOpDomainRhsPipeline(), {space});
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpValsDiffVals(vals, diff_vals, true));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateScalarFieldValues("FIELD1", ptr_values));
pipeline_mng->getOpDomainRhsPipeline().push_back(
ptr_diff_vals));
pipeline_mng->getOpDomainRhsPipeline().push_back(new OpCheckValsDiffVals(
vals, diff_vals, ptr_values, ptr_diff_vals, true));
CHKERR pipeline_mng->loopFiniteElements();
};
auto post_proc = [&] {
auto get_domain_post_proc_fe = [&](auto post_proc_mesh_ptr) {
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoabBaseCont<DomainEle>>(
m_field, post_proc_mesh_ptr);
post_proc_fe->getOpPtrVector(), {space});
auto ptr_values = boost::make_shared<VectorDouble>();
auto ptr_diff_vals = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues("FIELD1", ptr_values));
post_proc_fe->getOpPtrVector().push_back(
ptr_diff_vals));
post_proc_fe->getOpPtrVector().push_back(
new OpPPMap(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"FIELD1", ptr_values}},
{{"FIELD1_GRAD", ptr_diff_vals}},
{}, {})
);
return post_proc_fe;
};
auto get_bdy_post_proc_fe = [&](auto post_proc_mesh_ptr) {
auto bdy_post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoabBaseCont<BoundaryEle>>(
m_field, post_proc_mesh_ptr);
auto op_loop_side = new OpLoopSide<EleOnSide>(
m_field, simple->getDomainFEName(), SPACE_DIM, Sev::noisy,
boost::make_shared<
auto ptr_values = boost::make_shared<VectorDouble>();
auto ptr_diff_vals = boost::make_shared<MatrixDouble>();
// push operators to side element
op_loop_side->getOpPtrVector(), {space});
op_loop_side->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues("FIELD1", ptr_values));
op_loop_side->getOpPtrVector().push_back(
ptr_diff_vals));
// push op to boundary element
bdy_post_proc_fe->getOpPtrVector().push_back(op_loop_side);
bdy_post_proc_fe->getOpPtrVector().push_back(
new OpPPMap(
bdy_post_proc_fe->getPostProcMesh(),
bdy_post_proc_fe->getMapGaussPts(),
{{"FIELD1", ptr_values}},
{{"FIELD1_GRAD", ptr_diff_vals}},
{}, {})
);
return bdy_post_proc_fe;
};
auto post_proc_mesh_ptr = boost::make_shared<moab::Core>();
auto post_proc_begin =
boost::make_shared<PostProcBrokenMeshInMoabBaseBegin>(
m_field, post_proc_mesh_ptr);
auto post_proc_end =
boost::make_shared<PostProcBrokenMeshInMoabBaseEnd>(
m_field, post_proc_mesh_ptr);
auto domain_post_proc_fe = get_domain_post_proc_fe(post_proc_mesh_ptr);
auto bdy_post_proc_fe = get_bdy_post_proc_fe(post_proc_mesh_ptr);
CHKERR DMoFEMPreProcessFiniteElements(dm, post_proc_begin->getFEMethod());
CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
domain_post_proc_fe);
CHKERR DMoFEMLoopFiniteElements(dm, simple->getBoundaryFEName(),
bdy_post_proc_fe);
CHKERR DMoFEMPostProcessFiniteElements(dm, post_proc_end->getFEMethod());
CHKERR post_proc_end->writeFile("out_post_proc.h5m");
};
CHKERR assemble_matrices_and_vectors();
CHKERR solve_problem();
CHKERR check_solution();
CHKERR post_proc();
}
}
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreInterface::modify_finite_element_adjacency_table
virtual MoFEMErrorCode modify_finite_element_adjacency_table(const std::string &fe_name, const EntityType type, ElementAdjacencyFunct function)=0
modify finite element table, only for advanced user
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
help
static char help[]
Definition: scalar_check_approximation.cpp:13
MoFEM::ForcesAndSourcesCore::UserDataOperator::AdjCache
std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >> AdjCache
Definition: ForcesAndSourcesCore.hpp:904
EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION
Definition: plastic.cpp:13
BoundaryEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Definition: child_and_parent.cpp:39
ApproxFunctions
Definition: hcurl_check_approx_in_2d.cpp:44
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
SPACE_DIM
constexpr int SPACE_DIM
Definition: scalar_check_approximation.cpp:21
OpValsDiffVals
Definition: scalar_check_approximation.cpp:98
EleOnSide
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::FaceSideEle EleOnSide
Definition: scalar_check_approximation.cpp:27
MoFEM::PipelineManager::ElementsAndOpsByDim
Definition: PipelineManager.hpp:38
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::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:523
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::EntFiniteElement
Finite element data for entity.
Definition: FEMultiIndices.hpp:501
MoFEM::EntitiesFieldData::EntData::getDiffN
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
Definition: EntitiesFieldData.hpp:1329
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
ApproxFunctions::fUn
static FTensor::Tensor1< double, BASE_DIM > fUn(const double x, const double y, double z)
Definition: hcurl_check_approx_in_2d.cpp:46
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
sdf.r
int r
Definition: sdf.py:8
MoFEM::Field
Provide data structure for (tensor) field approximation.
Definition: FieldMultiIndices.hpp:51
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1293
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1502
main
int main(int argc, char *argv[])
Definition: scalar_check_approximation.cpp:284
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
MoFEM::createOpSchurAssembleEnd
OpSchurAssembleBase * createOpSchurAssembleEnd(std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range >> field_ents, SmartPetscObj< AO > ao, SmartPetscObj< Mat > schur, bool sym_schur, bool symm_op)
Construct a new Op Schur Assemble End object.
Definition: Schur.cpp:2186
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
OpMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition: seepage.cpp:57
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
convert.type
type
Definition: convert.py:64
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
SideNumber_multiIndex
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
Definition: RefEntsMultiIndices.hpp:101
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::LogManager::getStrmWorld
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::createOpSchurAssembleBegin
OpSchurAssembleBase * createOpSchurAssembleBegin()
Definition: Schur.cpp:2181
AINSWORTH_LOBATTO_BASE
@ AINSWORTH_LOBATTO_BASE
Definition: definitions.h:62
MoFEM::DMoFEMPreProcessFiniteElements
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:546
MoFEM::DMoFEMPostProcessFiniteElements
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:556
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
BiLinearForm
FTensor::Index< 'i', SPACE_DIM >
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
ElementsAndOps
Definition: child_and_parent.cpp:18
Range
MoFEM::EntitiesFieldData::EntData::getFTensor0FieldData
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Resturn scalar files as a FTensor of rank 0.
Definition: EntitiesFieldData.cpp:506
DomainEleOp
MoFEM::SideNumber
keeps information about side number for the finite element
Definition: RefEntsMultiIndices.hpp:57
MoFEM::CoreTmp< 0 >::Initialize
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
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
Definition: EntitiesFieldData.cpp:526
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
ApproxFunctionsImpl
Definition: scalar_check_approximation.cpp:17
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1318
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
ApproxFunctions::diffFun
static FTensor::Tensor2< double, BASE_DIM, SPACE_DIM > diffFun(const double x, const double y)
Definition: hcurl_check_approx_in_2d.cpp:68
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
EntData
EntitiesFieldData::EntData EntData
Definition: scalar_check_approximation.cpp:29
approx_order
int approx_order
Definition: test_broken_space.cpp:50
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:203
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
OpCheckValsDiffVals
Definition: hcurl_check_approx_in_2d.cpp:161
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:586
convert.int
int
Definition: convert.py:64
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
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::SCHUR
@ SCHUR
Definition: FormsIntegrators.hpp:106
MoFEM::OpLoopSide
Element used to execute operators on side of the element.
Definition: ForcesAndSourcesCore.hpp:1290
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
F
@ F
Definition: free_surface.cpp:394
NOFIELD
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition: definitions.h:84
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698