v0.14.0
SCL-9: Heat method

Table of Contents

Note
Prerequisites of this tutorial include SCL-1: Poisson's equation (homogeneous BC)


Note
Intended learning outcome:
  • general structure of a program developed using MoFEM
  • idea of Simple Interface in MoFEM and how to use it
  • Solving scalar problem embedded in 3d space
  • How to push the developed UDOs to the Pipeline

Source code

/**
* \file heat_method.cpp \example heat_method.cpp
*
* Calculating geodetic distances using heat method. For details see,
*
* K. Crane, C. Weischedel, M. Wardetzky, Geodesics in heat: A new approach to
* computing distance based on heat flow, ACM Transactions on Graphics , vol.
* 32, no. 5, pp. 152:1-152:11, 2013.
*
* Interent resources:
* http://www.cs.cmu.edu/~kmcrane/Projects/HeatMethod/
* http://www.numerical-tours.com/matlab/meshproc_7_geodesic_poisson/
*
*
* Solving two problems in sequence. Show hot to use form integrators and how to
* implement user data operator.
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
double dt = 2; // relative to edge length
#include <MoFEM.hpp>
constexpr int SPACE_DIM = 3;
// Use forms iterators for Grad-Grad term
// Use forms for Mass term
struct Example {
Example(MoFEM::Interface &m_field) : mField(m_field) {}
MoFEMErrorCode runProblem();
private:
Simple *simpleInterface;
Range pinchNodes;
MoFEMErrorCode readMesh();
MoFEMErrorCode setupProblem();
MoFEMErrorCode setIntegrationRules();
MoFEMErrorCode createCommonData();
MoFEMErrorCode boundaryCondition();
MoFEMErrorCode assembleSystem();
MoFEMErrorCode solveSystem();
MoFEMErrorCode outputResults();
MoFEMErrorCode checkResults();
/**
* Use problem specific implementation for second stage of heat methid
*/
struct OpRhs : public DomainEleOp {
OpRhs(boost::shared_ptr<MatrixDouble> u_grad_ptr);
MoFEMErrorCode doWork(int side, EntityType type,
protected:
boost::shared_ptr<MatrixDouble> uGradPtr;
};
};
//! [Run programme]
CHKERR readMesh();
CHKERR setupProblem();
CHKERR createCommonData();
CHKERR boundaryCondition();
CHKERR assembleSystem();
CHKERR setIntegrationRules();
CHKERR solveSystem();
CHKERR outputResults();
CHKERR checkResults();
}
//! [Run programme]
//! [Read mesh]
CHKERR mField.getInterface(simpleInterface);
CHKERR simpleInterface->getOptions();
CHKERR simpleInterface->loadFile();
// FIXME: THis part of algorithm is not working in parallel. If you have bit
// of free time, feel free to generalise code below.
Range nodes;
CHKERR mField.get_moab().get_entities_by_type(0, MBVERTEX, nodes);
// pinchNodes could be get from BLOCKSET
pinchNodes.insert(nodes[0]);
Range edges;
CHKERR mField.get_moab().get_adjacencies(pinchNodes, 1, false, edges,
moab::Interface::UNION);
double l2;
for (auto e : edges) {
Range nodes;
CHKERR mField.get_moab().get_connectivity(Range(e, e), nodes, false);
MatrixDouble coords(nodes.size(), 3);
CHKERR mField.get_moab().get_coords(nodes, &coords(0, 0));
double l2e = 0;
for (int j = 0; j != 3; ++j) {
l2e += pow(coords(0, j) - coords(1, j), 2);
}
l2 = std::max(l2, l2e);
}
dt *= std::sqrt(l2);
}
//! [Read mesh]
//! [Set up problem]
// Only one field
CHKERR simpleInterface->addDomainField("U", H1, AINSWORTH_LEGENDRE_BASE, 1);
constexpr int order = 1;
CHKERR simpleInterface->setFieldOrder("U", order);
CHKERR simpleInterface->setUp();
}
//! [Set up problem]
//! [Set integration rule]
// Set integration order. To 2p is enough to integrate mass matrix exactly.
auto rule = [](int, int, int p) -> int { return 2 * p; };
// Set integration rule for elements assembling matrix and vector. Note that
// rule for vector could be 2p-1, but that not change computation time
// significantly. So shorter code is better.
PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
}
//! [Set integration rule]
//! [Create common data]
}
//! [Create common data]
//! [Boundary condition]
}
//! [Boundary condition]
//! [Push operators to pipeline]
PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
// This will store gradients at integration points on element
auto grad_u_ptr = boost::make_shared<MatrixDouble>();
// Push element from reference configuration to current configuration in 3d
// space
auto set_domain_general = [&](auto &pipeline) {
auto det_ptr = boost::make_shared<VectorDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
pipeline.push_back(new OpCalculateHOJacForFaceEmbeddedIn3DSpace(jac_ptr));
pipeline.push_back(new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
pipeline.push_back(new OpSetInvJacH1ForFaceEmbeddedIn3DSpace(inv_jac_ptr));
};
// Operators for assembling matrix for first stage
auto set_domain_lhs_first = [&](auto &pipeline) {
auto one = [](double, double, double) -> double { return 1.; };
pipeline.push_back(new OpGradGrad("U", "U", one));
auto rho = [](double, double, double) -> double { return 1. / dt; };
pipeline.push_back(new OpMass("U", "U", rho));
};
// Nothing is assembled for vector for first stage of heat method. Later
// simply value of one is set to elements of right hand side vector at pinch
// nodes.
auto set_domain_rhs_first = [&](auto &pipeline) {};
// Operators for assembly of left hand side. This time only Grad-Grad
// operator.
auto set_domain_lhs_second = [&](auto &pipeline) {
auto one = [](double, double, double) { return 1.; };
pipeline.push_back(new OpGradGrad("U", "U", one));
};
// Now apply on the right hand side vector X = Grad/||Grad||
auto set_domain_rhs_second = [&](auto &pipeline) {
pipeline.push_back(new OpCalculateScalarFieldGradient<3>("U", grad_u_ptr));
pipeline.push_back(new OpRhs(grad_u_ptr));
};
// Solver first problem
auto solve_first = [&]() {
auto simple = mField.getInterface<Simple>();
auto pipeline_mng = mField.getInterface<PipelineManager>();
auto solver = pipeline_mng->createKSP();
CHKERR KSPSetFromOptions(solver);
CHKERR KSPSetUp(solver);
auto dm = simpleInterface->getDM();
auto D = createDMVector(dm);
auto F = vectorDuplicate(D);
// Note add one at nodes which are pinch nodes
CHKERR VecZeroEntries(F);
auto problem_ptr = mField.get_problem(simple->getProblemName());
auto bit_number = mField.get_field_bit_number("U");
auto dofs_ptr = problem_ptr->getNumeredRowDofsPtr();
for (auto v : pinchNodes) {
auto dof = dofs_ptr->get<Unique_mi_tag>().find(uid);
if (dof != dofs_ptr->end())
VecSetValue(F, (*dof)->getPetscGlobalDofIdx(), 1, INSERT_VALUES);
}
CHKERR VecAssemblyBegin(F);
CHKERR VecAssemblyEnd(F);
// Solve problem
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);
};
// Solve second problem
auto solve_second = [&]() {
auto simple = mField.getInterface<Simple>();
// Note that DOFs on pinch nodes are removed from the problem
auto prb_mng = mField.getInterface<ProblemsManager>();
CHKERR prb_mng->removeDofsOnEntities(simple->getProblemName(), "U",
pinchNodes, 0, 1);
auto *pipeline_mng = mField.getInterface<PipelineManager>();
auto solver = pipeline_mng->createKSP();
CHKERR KSPSetFromOptions(solver);
CHKERR KSPSetUp(solver);
auto dm = simpleInterface->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);
};
// Post-process results
auto post_proc = [&](const std::string out_name) {
PipelineManager *pipeline_mng = mField.getInterface<PipelineManager>();
auto tmp_lhs_fe = pipeline_mng->getDomainLhsFE();
auto tmp_rhs_fe = pipeline_mng->getDomainRhsFE();
auto det_ptr = boost::make_shared<VectorDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getDomainLhsFE().reset();
pipeline_mng->getDomainRhsFE().reset();
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(mField);
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
new OpInvertMatrix<SPACE_DIM>(jac_ptr, det_ptr, inv_jac_ptr));
post_proc_fe->getOpPtrVector().push_back(
auto u_ptr = boost::make_shared<VectorDouble>();
auto grad_ptr = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues("U", u_ptr));
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(new OpPPMap(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"U", u_ptr}},
{{"GRAD_U", grad_ptr}},
{}, {}));
pipeline_mng->getDomainRhsFE() = post_proc_fe;
CHKERR pipeline_mng->loopFiniteElements();
CHKERR post_proc_fe->writeFile(out_name);
tmp_lhs_fe = pipeline_mng->getDomainLhsFE() = tmp_lhs_fe;
tmp_rhs_fe = pipeline_mng->getDomainRhsFE() = tmp_rhs_fe;
};
set_domain_general(pipeline_mng->getOpDomainLhsPipeline());
set_domain_general(pipeline_mng->getOpDomainRhsPipeline());
set_domain_lhs_first(pipeline_mng->getOpDomainLhsPipeline());
set_domain_rhs_first(pipeline_mng->getOpDomainRhsPipeline());
CHKERR solve_first();
CHKERR post_proc("out_heat_method_first.h5m");
pipeline_mng->getOpDomainLhsPipeline().clear();
pipeline_mng->getOpDomainRhsPipeline().clear();
set_domain_general(pipeline_mng->getOpDomainLhsPipeline());
set_domain_general(pipeline_mng->getOpDomainRhsPipeline());
set_domain_lhs_second(pipeline_mng->getOpDomainLhsPipeline());
set_domain_rhs_second(pipeline_mng->getOpDomainRhsPipeline());
CHKERR solve_second();
CHKERR post_proc("out_heat_method_second.h5m");
};
//! [Push operators to pipeline]
//! [Solve]
}
//! [Solve]
//! [Postprocess results]
}
//! [Postprocess results]
//! [Check results]
}
//! [Check results]
int main(int argc, char *argv[]) {
// Initialisation of MoFEM/PETSc and MOAB data structures
const char param_file[] = "param_file.petsc";
MoFEM::Core::Initialize(&argc, &argv, param_file, help);
// Add logging channel for example
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "EXAMPLE"));
LogManager::setLog("EXAMPLE");
MOFEM_LOG_TAG("EXAMPLE", "example");
try {
//! [Register MoFEM discrete manager in PETSc]
DMType dm_name = "DMMOFEM";
//! [Register MoFEM discrete manager in PETSc
//! [Create MoAB]
moab::Core mb_instance; ///< mesh database
moab::Interface &moab = mb_instance; ///< mesh database interface
//! [Create MoAB]
//! [Create MoFEM]
MoFEM::Core core(moab); ///< finite element database
MoFEM::Interface &m_field = core; ///< finite element database insterface
//! [Create MoFEM]
//! [Example]
Example ex(m_field);
CHKERR ex.runProblem();
//! [Example]
}
}
Example::OpRhs::OpRhs(boost::shared_ptr<MatrixDouble> u_grad_ptr)
: DomainEleOp("U", DomainEleOp::OPROW), uGradPtr(u_grad_ptr) {}
auto nb_dofs = data.getIndices().size();
if (nb_dofs) {
auto t_grad = getFTensor1FromMat<3>(*uGradPtr);
auto nb_base_functions = data.getN().size2();
auto nb_gauss_pts = getGaussPts().size2();
std::array<double, MAX_DOFS_ON_ENTITY> nf;
std::fill(nf.begin(), &nf[nb_dofs], 0);
auto t_diff_base = data.getFTensor1DiffN<3>();
auto t_w = getFTensor0IntegrationWeight();
auto a = getMeasure();
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
double alpha = t_w * a;
const auto l2 = t_grad(i) * t_grad(i);
if (l2 > std::numeric_limits<double>::epsilon())
t_one(i) = t_grad(i) / std::sqrt(l2);
else
t_one(i) = 0;
size_t bb = 0;
for (; bb != nb_dofs; ++bb) {
nf[bb] -= alpha * t_diff_base(i) * t_one(i);
++t_diff_base;
}
for (; bb < nb_base_functions; ++bb) {
++t_diff_base;
}
++t_grad;
}
CHKERR VecSetValues<MoFEM::EssentialBcStorage>(getKSPf(), data, &nf[0],
ADD_VALUES);
}
}
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
Example::checkResults
MoFEMErrorCode checkResults()
[Postprocess results]
Definition: dynamic_first_order_con_law.cpp:1205
MoFEM::PipelineManager::getDomainRhsFE
boost::shared_ptr< FEMethod > & getDomainRhsFE()
Definition: PipelineManager.hpp:405
Example::assembleSystem
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
Definition: dynamic_first_order_con_law.cpp:647
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
FTensor::Tensor1< double, 3 >
MoFEM::OpCalculateHOJacForFaceEmbeddedIn3DSpace
OpCalculateHOJacForFaceImpl< 3 > OpCalculateHOJacForFaceEmbeddedIn3DSpace
Definition: HODataOperators.hpp:265
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
rho
double rho
Definition: plastic.cpp:140
MoFEM::PipelineManager::loopFiniteElements
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
Definition: PipelineManager.cpp:19
OpGradGrad
FormsIntegrators< DomainEleOp >::Assembly< A >::BiLinearForm< I >::OpGradGrad< 1, FIELD_DIM, SPACE_DIM > OpGradGrad
Definition: operators_tests.cpp:44
MoFEM::DofEntity::getUniqueIdCalculate
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
Definition: DofsMultiIndices.hpp:54
help
static char help[]
Definition: activate_deactivate_dofs.cpp:13
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::PipelineManager::getOpDomainRhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
Definition: PipelineManager.hpp:773
MoFEM::CoreInterface::get_field_bit_number
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
MoFEM::CoreInterface::get_problem
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
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::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
order
constexpr int order
Definition: dg_projection.cpp:18
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
Example::readMesh
MoFEMErrorCode readMesh()
[Run problem]
Definition: dynamic_first_order_con_law.cpp:463
MoFEM::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:747
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
Example::OpRhs::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Definition: heat_method.cpp:425
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
a
constexpr double a
Definition: approx_sphere.cpp:30
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
Example
[Example]
Definition: plastic.cpp:177
double
convert.type
type
Definition: convert.py:64
MoFEM::PipelineManager::FaceEle
MoFEM::FaceElementForcesAndSourcesCore FaceEle
Definition: PipelineManager.hpp:35
MoFEM::FieldEntity::getLocalUniqueIdCalculate
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
Definition: FieldEntsMultiIndices.hpp:136
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
Example::solveSystem
MoFEMErrorCode solveSystem()
[Solve]
Definition: dynamic_first_order_con_law.cpp:893
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
Example::OpRhs::uGradPtr
boost::shared_ptr< MatrixDouble > uGradPtr
Definition: heat_method.cpp:77
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::Problem::getNumeredRowDofsPtr
auto & getNumeredRowDofsPtr() const
get access to numeredRowDofsPtr storing DOFs on rows
Definition: ProblemsMultiIndices.hpp:82
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:136
Example::setIntegrationRules
MoFEMErrorCode setIntegrationRules()
[Set up problem]
Definition: plot_base.cpp:229
MoFEM::PipelineManager::setDomainRhsIntegrationRule
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:530
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::PipelineManager::getOpDomainLhsPipeline
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
Definition: PipelineManager.hpp:749
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
BiLinearForm
main
int main(int argc, char *argv[])
Definition: activate_deactivate_dofs.cpp:15
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
FTensor::Index< 'i', 3 >
OpRhs
Definition: approx_sphere.cpp:68
MoFEM::PipelineManager::setDomainLhsIntegrationRule
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Definition: PipelineManager.hpp:503
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
MoFEM::PipelineManager::getDomainLhsFE
boost::shared_ptr< FEMethod > & getDomainLhsFE()
Definition: PipelineManager.hpp:401
DomainEleOp
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::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
Example::OpRhs::OpRhs
OpRhs(boost::shared_ptr< MatrixDouble > u_grad_ptr)
Definition: heat_method.cpp:422
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
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
Example::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: plastic.cpp:233
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
Example::boundaryCondition
MoFEMErrorCode boundaryCondition()
[Set up problem]
Definition: dynamic_first_order_con_law.cpp:536
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::OpInvertMatrix
Definition: UserDataOperators.hpp:3249
Example::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: plastic.cpp:214
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEM::Unique_mi_tag
Definition: TagMultiIndices.hpp:18
Example::createCommonData
MoFEMErrorCode createCommonData()
[Set up problem]
Definition: plastic.cpp:404
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
dt
double dt
Definition: heat_method.cpp:26
Example::outputResults
MoFEMErrorCode outputResults()
[Solve]
Definition: dynamic_first_order_con_law.cpp:1183
MoFEM::OpSetInvJacH1ForFaceEmbeddedIn3DSpace
Definition: UserDataOperators.hpp:2932
convert.int
int
Definition: convert.py:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
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
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698