|
| v0.14.0
|
Simple interface for fast problem set-up.
More...
#include <src/interfaces/Simple.hpp>
|
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
|
| Simple (const MoFEM::Core &core) |
|
virtual | ~Simple ()=default |
|
MoFEMErrorCode | getOptions () |
| get options More...
|
|
MoFEMErrorCode | loadFile (const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc) |
| Load mesh file. More...
|
|
MoFEMErrorCode | loadFile (const std::string mesh_file_name="") |
| Load mesh file with parallel options if number of cores > 1. More...
|
|
MoFEMErrorCode | addDomainField (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_ZERO, int verb=-1) |
| Add field on domain. More...
|
|
MoFEMErrorCode | addDomainBrokenField (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_ZERO, int verb=-1) |
| Add broken field on domain. More...
|
|
MoFEMErrorCode | addBoundaryField (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_ZERO, int verb=-1) |
| Add field on boundary. More...
|
|
MoFEMErrorCode | addSkeletonField (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_ZERO, int verb=-1) |
| Add field on skeleton. More...
|
|
MoFEMErrorCode | addDataField (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_ZERO, int verb=-1) |
| Add field on domain. More...
|
|
MoFEMErrorCode | removeDomainField (const std::string &name) |
| Remove field form domain. More...
|
|
MoFEMErrorCode | removeBoundaryField (const std::string &name) |
| Remove field form boundary. More...
|
|
MoFEMErrorCode | removeSkeletonField (const std::string &name) |
| Remove field form skeleton. More...
|
|
MoFEMErrorCode | defineFiniteElements () |
| Define finite elements. More...
|
|
MoFEMErrorCode | defineProblem (const PetscBool is_partitioned=PETSC_TRUE) |
| define problem More...
|
|
MoFEMErrorCode | setFieldOrder (const std::string field_name, const int order, const Range *ents=NULL) |
| Set field order. More...
|
|
MoFEMErrorCode | buildFields () |
| Build fields. More...
|
|
MoFEMErrorCode | buildFiniteElements () |
| Build finite elements. More...
|
|
MoFEMErrorCode | setSkeletonAdjacency (int dim=-1, std::string fe_name="") |
| Set the skeleton adjacency object. More...
|
|
MoFEMErrorCode | buildProblem () |
| Build problem. More...
|
|
MoFEMErrorCode | setUp (const PetscBool is_partitioned=PETSC_TRUE) |
| Setup problem. More...
|
|
MoFEMErrorCode | reSetUp (bool only_dm=false) |
| Rebuild internal MoFEM data structures. More...
|
|
MoFEMErrorCode | getDM (DM *dm) |
| Get DM. More...
|
|
SmartPetscObj< DM > | getDM () |
| Return smart DM object. More...
|
|
int | getDim () const |
| Get the problem dimension. More...
|
|
void | setDim (int dim) |
| Set the problem dimension. More...
|
|
DEPRECATED EntityHandle & | getMeshSet () |
|
EntityHandle & | getMeshset () |
| Get the MeshSet object. More...
|
|
EntityHandle & | getBoundaryMeshSet () |
| Get the BoundaryMeshSet object. More...
|
|
EntityHandle & | getSkeletonMeshSet () |
| Get the SkeletonMeshSet object. More...
|
|
BitRefLevel & | getBitRefLevel () |
| Get the BitRefLevel. More...
|
|
BitRefLevel & | getBitRefLevelMask () |
| Get the BitRefLevel. More...
|
|
const std::string | getDomainFEName () const |
| Get the Domain FE Name. More...
|
|
const std::string | getBoundaryFEName () const |
| Get the Boundary FE Name. More...
|
|
const std::string | getSkeletonFEName () const |
| Get the Skeleton FE Name. More...
|
|
const std::string | getProblemName () const |
| Get the Problem Name. More...
|
|
std::string & | getDomainFEName () |
| Get the Domain FE Name. More...
|
|
std::string & | getBoundaryFEName () |
| Get the Boundary FE Name. More...
|
|
std::string & | getSkeletonFEName () |
| Get the Skeleton FE Name. More...
|
|
std::string & | getProblemName () |
| Get the Problem Name. More...
|
|
std::vector< std::string > & | getOtherFiniteElements () |
| Get the Other Finite Elements. More...
|
|
MoFEMErrorCode | deleteDM () |
| Delete dm. More...
|
|
MoFEMErrorCode | deleteFiniteElements () |
| Delete finite elements. More...
|
|
bool & | getAddSkeletonFE () |
| Get the addSkeletonFE. More...
|
|
bool & | getAddBoundaryFE () |
| Get the addSkeletonFE. More...
|
|
bool & | getParentAdjacencies () |
| Get the addParentAdjacencies. More...
|
|
auto & | getBitAdjParent () |
| bit ref level for parent More...
|
|
auto & | getBitAdjParentMask () |
| bit ref level for parent parent More...
|
|
auto & | getBitAdjEnt () |
| bit ref level for parent More...
|
|
auto & | getBitAdjEntMask () |
| bit ref level for parent parent More...
|
|
MoFEMErrorCode | addFieldToEmptyFieldBlocks (const std::string row_field, const std::string col_field) const |
| add empty block to problem More...
|
|
template<class IFACE > |
MoFEMErrorCode | registerInterface (bool error_if_registration_failed=true) |
| Register interface. More...
|
|
template<class IFACE > |
MoFEMErrorCode | getInterface (IFACE *&iface) const |
| Get interface reference to pointer of interface. More...
|
|
template<class IFACE > |
MoFEMErrorCode | getInterface (IFACE **const iface) const |
| Get interface pointer to pointer of interface. More...
|
|
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0> |
IFACE | getInterface () const |
| Get interface pointer to pointer of interface. More...
|
|
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0> |
IFACE | getInterface () const |
| Get reference to interface. More...
|
|
template<class IFACE > |
IFACE * | getInterface () const |
| Function returning pointer to interface. More...
|
|
virtual | ~UnknownInterface ()=default |
|
Simple interface for fast problem set-up.
- Examples
- adolc_plasticity.cpp, approx_sphere.cpp, boundary_marker.cpp, child_and_parent.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, dg_projection.cpp, dynamic_first_order_con_law.cpp, eigen_elastic.cpp, field_evaluator.cpp, free_surface.cpp, hanging_node_approx.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hdiv_check_approx_in_3d.cpp, hdiv_divergence_operator.cpp, heat_equation.cpp, heat_method.cpp, hello_world.cpp, helmholtz.cpp, higher_derivatives.cpp, level_set.cpp, loop_entities.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, operators_tests.cpp, partition_mesh.cpp, phase.cpp, photon_diffusion.cpp, plastic.cpp, plate.cpp, plot_base.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, scalar_check_approximation.cpp, schur_test_diag_mat.cpp, seepage.cpp, shallow_wave.cpp, simple_elasticity.cpp, simple_interface.cpp, simple_l2_only.cpp, tensor_divergence_operator.cpp, test_broken_space.cpp, test_cache_on_entities.cpp, testing_jacobian_of_hook_element.cpp, testing_jacobian_of_hook_scaled_with_density_element.cpp, thermo_elastic.cpp, and wave_equation.cpp.
Definition at line 27 of file Simple.hpp.
◆ LoadFileFunc
◆ Simple()
Definition at line 162 of file Simple.cpp.
173 PetscLogEventRegister(
"SimpBuildFEs", 0,
◆ ~Simple()
virtual MoFEM::Simple::~Simple |
( |
| ) |
|
|
virtualdefault |
◆ addBoundaryField()
Add field on boundary.
- Parameters
-
name | name of the field |
space | space (L2,H1,Hdiv,Hcurl) |
base | approximation base, see FieldApproximationBase |
nb_of_coefficients | number of field coefficients |
tag_type | type of the tag MB_TAG_DENSE or MB_TAG_SPARSE (DENSE is faster and uses less memory, SPARSE is more flexible if you define field on subdomains) |
bh | if MF_EXCL throws error if field exits, MF_ZERO no error if field exist |
verb | verbosity level |
- Returns
- error code
- Examples
- child_and_parent.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, mesh_smoothing.cpp, and simple_interface.cpp.
Definition at line 354 of file Simple.cpp.
361 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
367 "NOFIELD space for boundary filed not implemented in Simple interface");
◆ addDataField()
Add field on domain.
- Parameters
-
name | name of the field |
space | space (L2,H1,Hdiv,Hcurl) |
base | approximation base, see FieldApproximationBase |
nb_of_coefficients | number of field coefficients |
tag_type | type of the tag MB_TAG_DENSE or MB_TAG_SPARSE (DENSE is faster and uses less memory, SPARSE is more flexible if you define field on subdomains) |
bh | if MF_EXCL throws error if field exits, MF_ZERO no error if field exist |
verb | verbosity level |
- Returns
- error code
- Examples
- hcurl_curl_operator.cpp, and hdiv_divergence_operator.cpp.
Definition at line 392 of file Simple.cpp.
400 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
◆ addDomainBrokenField()
Add broken field on domain.
- Parameters
-
name | name of the field |
space | space (L2,H1,Hdiv,Hcurl) |
base | approximation base, see FieldApproximationBase |
nb_of_coefficients | number of field coefficients |
tag_type | type of the tag MB_TAG_DENSE or MB_TAG_SPARSE (DENSE is faster and uses less memory, SPARSE is more flexible if you define field on subdomains) |
bh | if MF_EXCL throws error if field exits, MF_ZERO no error if field exist |
verb | verbosity level |
- Returns
- error code
Definition at line 282 of file Simple.cpp.
291 auto get_side_map_hcurl = [&]() {
294 std::pair<EntityType,
308 auto get_side_map_hdiv = [&]() {
311 std::pair<EntityType,
325 auto get_side_map = [&](
auto space) {
328 return get_side_map_hcurl();
330 return get_side_map_hdiv();
334 <<
"Side dof map for broken space <" << space <<
"> not implemented";
338 std::pair<EntityType,
344 CHKERR m_field.add_broken_field(name, space, base, nb_of_coefficients,
345 get_side_map(space), tag_type, bh, verb);
◆ addDomainField()
Add field on domain.
- Parameters
-
name | name of the field |
space | space (L2,H1,Hdiv,Hcurl) |
base | approximation base, see FieldApproximationBase |
nb_of_coefficients | number of field coefficients |
tag_type | type of the tag MB_TAG_DENSE or MB_TAG_SPARSE (DENSE is faster and uses less memory, SPARSE is more flexible if you define field on subdomains) |
bh | if MF_EXCL throws error if field exits, MF_ZERO no error if field exist |
verb | verbosity level |
- Returns
- error code
- Examples
- child_and_parent.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, loop_entities.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, simple_elasticity.cpp, simple_interface.cpp, simple_l2_only.cpp, test_cache_on_entities.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.
Definition at line 264 of file Simple.cpp.
272 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
◆ addFieldToEmptyFieldBlocks()
MoFEMErrorCode MoFEM::Simple::addFieldToEmptyFieldBlocks |
( |
const std::string |
row_field, |
|
|
const std::string |
col_field |
|
) |
| const |
add empty block to problem
MatrixManager assumes that all blocks, i.e. all fields combinations are non zero. This is not always the case, to optimise code and reduce memory usage user can specifi which blocks are empty.
- Parameters
-
row_field | row field name |
col_field | col field name |
- Returns
- MoFEMErrorCode
Definition at line 782 of file Simple.cpp.
◆ addSkeletonField()
Add field on skeleton.
- Parameters
-
name | name of the field |
space | space (L2,H1,Hdiv,Hcurl) |
base | approximation base, see FieldApproximationBase |
nb_of_coefficients | number of field coefficients |
tag_type | type of the tag MB_TAG_DENSE or MB_TAG_SPARSE (DENSE is faster and uses less memory, SPARSE is more flexible if you define field on subdomains) |
bh | if MF_EXCL throws error if field exits, MF_ZERO no error if field exist |
verb | verbosity level |
- Returns
- error code
- Examples
- continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, and simple_interface.cpp.
Definition at line 372 of file Simple.cpp.
380 CHKERR m_field.add_field(name, space, base, nb_of_coefficients, tag_type, bh,
386 "NOFIELD space for boundary filed not implemented in Simple interface");
◆ buildFields()
Build fields.
- Returns
- error code
- Examples
- mesh_smoothing.cpp, and simple_elasticity.cpp.
Definition at line 554 of file Simple.cpp.
559 auto comm_interface_ptr = m_field.getInterface<CommInterface>();
560 auto bit_ref_ptr = m_field.getInterface<BitRefManager>();
563 auto add_ents_to_field = [&](
auto meshset,
auto dim,
auto &fields) {
565 for (
auto &field : fields) {
566 CHKERR m_field.add_ents_to_field_by_dim(meshset, dim, field);
567 CHKERR comm_interface_ptr->synchroniseFieldEntities(field, 0);
572 auto make_no_field_ents = [&](
auto &fields) {
574 for (
auto &field : fields) {
575 std::array<double, 6> coords = {0, 0, 0, 0, 0, 0};
576 CHKERR m_field.create_vertices_and_add_to_field(field, coords.data(), 1);
577 CHKERR comm_interface_ptr->makeFieldEntitiesMultishared(field, 0);
588 std::set<std::string> nofield_fields;
590 nofield_fields.insert(
f);
592 nofield_fields.insert(
f);
594 CHKERR make_no_field_ents(nofield_fields);
598 auto get_field_ptr = [&](
auto &
f) {
return m_field.get_field_structure(
f); };
600 const auto f = std::get<0>(
t);
601 const auto order = std::get<1>(
t);
604 <<
"Set order to field " <<
f <<
" order " <<
order;
606 if (std::get<3>(
t)) {
608 <<
"To ents: " << std::endl
609 << std::get<2>(
t) << std::endl;
613 if (std::get<3>(
t)) {
618 auto f_ptr = get_field_ptr(
f);
620 if (f_ptr->getSpace() ==
H1) {
628 for (
auto d = 1;
d <=
dIm; ++
d) {
629 for (EntityType
t = CN::TypeDimensionMap[
d].first;
630 t <= CN::TypeDimensionMap[
d].second; ++
t) {
638 CHKERR m_field.build_fields();
◆ buildFiniteElements()
◆ buildProblem()
Build problem.
- Returns
- error code
- Examples
- mesh_smoothing.cpp, and simple_elasticity.cpp.
Definition at line 669 of file Simple.cpp.
676 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
678 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
◆ createBoundaryMeshset()
Definition at line 791 of file Simple.cpp.
794 ParallelComm *pcomm =
797 auto get_skin = [&](
auto meshset) {
801 CHKERR m_field.get_moab().get_entities_by_dimension(meshset,
dIm,
803 CHKERR pcomm->filter_pstatus(domain_ents,
804 PSTATUS_SHARED | PSTATUS_MULTISHARED,
805 PSTATUS_NOT, -1,
nullptr);
807 Skinner skin(&m_field.get_moab());
809 CHKERR skin.find_skin(0, domain_ents,
false, domain_skin);
810 CHKERR pcomm->filter_pstatus(domain_skin,
811 PSTATUS_SHARED | PSTATUS_MULTISHARED,
812 PSTATUS_NOT, -1,
nullptr);
816 auto create_boundary_meshset = [&](
auto &&domain_skin) {
826 CHKERR m_field.get_moab().get_adjacencies(domain_skin,
dd,
false, adj,
827 moab::Interface::UNION);
◆ createSkeletonMeshset()
Definition at line 838 of file Simple.cpp.
842 ParallelComm *pcomm =
845 auto create_skeleton_meshset = [&](
auto meshset) {
851 Range boundary_ents, skeleton_ents;
853 dIm - 1, boundary_ents);
855 CHKERR m_field.get_moab().get_entities_by_dimension(meshset,
dIm,
857 CHKERR m_field.get_moab().get_adjacencies(
858 domain_ents,
dIm - 1,
false, skeleton_ents, moab::Interface::UNION);
859 skeleton_ents = subtract(skeleton_ents, boundary_ents);
860 CHKERR pcomm->filter_pstatus(skeleton_ents, PSTATUS_NOT_OWNED, PSTATUS_NOT,
◆ defaultLoadFileFunc()
static MoFEMErrorCode MoFEM::Simple::defaultLoadFileFunc |
( |
Interface & |
m_field, |
|
|
const char * |
file_name, |
|
|
const char * |
options |
|
) |
| |
|
inlinestatic |
Set function for loading mesh file.
- Parameters
-
m_field | interface |
file_name | file name |
options | loader options |
- Returns
- error code
Definition at line 51 of file Simple.hpp.
55 return m_field.get_moab().load_file(file_name, 0, options);
◆ defineFiniteElements()
Define finite elements.
- Returns
- Error code
- Examples
- mesh_smoothing.cpp, and simple_elasticity.cpp.
Definition at line 452 of file Simple.cpp.
456 auto clear_rows_and_cols = [&](
auto &fe_name) {
458 auto fe_ptr = m_field.get_finite_elements();
460 ->get<FiniteElement_name_mi_tag>();
461 auto it_fe = fe_by_name.find(fe_name);
462 if (it_fe != fe_by_name.end()) {
464 if (!fe_by_name.modify(it_fe, FiniteElement_row_change_bit_reset()))
466 "modification unsuccessful");
468 if (!fe_by_name.modify(it_fe, FiniteElement_col_change_bit_reset()))
470 "modification unsuccessful");
481 auto add_fields = [&](
auto &fe_name,
auto &fields) {
483 for (
auto &field : fields) {
484 CHKERR m_field.modify_finite_element_add_field_row(fe_name, field);
485 CHKERR m_field.modify_finite_element_add_field_col(fe_name, field);
486 CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
491 auto add_data_fields = [&](
auto &fe_name,
auto &fields) {
493 for (
auto &field : fields)
494 CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
◆ defineProblem()
MoFEMErrorCode MoFEM::Simple::defineProblem |
( |
const PetscBool |
is_partitioned = PETSC_TRUE | ) |
|
◆ deleteDM()
Delete dm.
- Returns
- MoFEMErrorCode
Definition at line 759 of file Simple.cpp.
◆ deleteFiniteElements()
Delete finite elements.
- Returns
- MoFEMErrorCode
Definition at line 770 of file Simple.cpp.
774 if (m_field.check_finite_element(fe)) {
775 CHKERR m_field.delete_finite_element(fe);
◆ exchangeGhostCells()
get all adjacent ghosted entities
Definition at line 872 of file Simple.cpp.
875 MOFEM_LOG(
"WORLD", Sev::verbose) <<
"Exchange ghost cells";
877 ParallelComm *pcomm =
880 pcomm =
new ParallelComm(&m_field.get_moab(), m_field.get_comm());
887 CHKERR m_field.get_moab().get_entities_by_dimension(0,
dIm, shared);
888 for (
auto d =
dIm - 1;
d >= 0; --
d) {
889 CHKERR m_field.get_moab().get_adjacencies(shared,
d,
false, shared,
890 moab::Interface::UNION);
892 CHKERR pcomm->filter_pstatus(shared, PSTATUS_SHARED | PSTATUS_MULTISHARED,
893 PSTATUS_OR, -1, &shared);
894 Tag part_tag = pcomm->part_tag();
895 CHKERR pcomm->exchange_tags(part_tag, shared);
896 CHKERR m_field.getInterface<CommInterface>()->resolveParentEntities(shared,
◆ getAddBoundaryFE()
bool& MoFEM::Simple::getAddBoundaryFE |
( |
| ) |
|
|
inline |
◆ getAddSkeletonFE()
bool& MoFEM::Simple::getAddSkeletonFE |
( |
| ) |
|
|
inline |
◆ getBitAdjEnt()
auto& MoFEM::Simple::getBitAdjEnt |
( |
| ) |
|
|
inline |
bit ref level for parent
- Returns
- auto&
Definition at line 494 of file Simple.hpp.
◆ getBitAdjEntMask()
auto& MoFEM::Simple::getBitAdjEntMask |
( |
| ) |
|
|
inline |
bit ref level for parent parent
- Returns
- auto&
Definition at line 501 of file Simple.hpp.
◆ getBitAdjParent()
auto& MoFEM::Simple::getBitAdjParent |
( |
| ) |
|
|
inline |
bit ref level for parent
- Returns
- auto&
Definition at line 480 of file Simple.hpp.
◆ getBitAdjParentMask()
auto& MoFEM::Simple::getBitAdjParentMask |
( |
| ) |
|
|
inline |
bit ref level for parent parent
- Returns
- auto&
Definition at line 487 of file Simple.hpp.
◆ getBitRefLevel()
◆ getBitRefLevelMask()
◆ getBoundaryFEName() [1/2]
std::string& MoFEM::Simple::getBoundaryFEName |
( |
| ) |
|
|
inline |
Get the Boundary FE Name.
- Returns
- std::string&
Definition at line 403 of file Simple.hpp.
◆ getBoundaryFEName() [2/2]
const std::string MoFEM::Simple::getBoundaryFEName |
( |
| ) |
const |
|
inline |
◆ getBoundaryMeshSet()
Get the BoundaryMeshSet object.
- Returns
- EntityHandle&
Definition at line 340 of file Simple.hpp.
◆ getDim()
int MoFEM::Simple::getDim |
( |
| ) |
const |
|
inline |
◆ getDM() [1/2]
◆ getDM() [2/2]
Get DM.
- Parameters
-
- Returns
- error code
- Examples
- child_and_parent.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, helmholtz.cpp, loop_entities.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, simple_elasticity.cpp, simple_interface.cpp, test_cache_on_entities.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.
Definition at line 747 of file Simple.cpp.
◆ getDomainFEName() [1/2]
std::string& MoFEM::Simple::getDomainFEName |
( |
| ) |
|
|
inline |
Get the Domain FE Name.
- Returns
- std::string&
Definition at line 396 of file Simple.hpp.
◆ getDomainFEName() [2/2]
const std::string MoFEM::Simple::getDomainFEName |
( |
| ) |
const |
|
inline |
◆ getMeshSet()
◆ getMeshset()
Get the MeshSet object.
- Returns
- EntityHandle&
Definition at line 333 of file Simple.hpp.
◆ getOptions()
get options
- Returns
- error code
- Examples
- child_and_parent.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, loop_entities.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, simple_elasticity.cpp, simple_interface.cpp, simple_l2_only.cpp, test_cache_on_entities.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.
Definition at line 180 of file Simple.cpp.
181 PetscBool flg = PETSC_TRUE;
183 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,
"",
"Simple interface options",
186 ierr = PetscOptionsString(
"-file_name",
"file name",
"",
"mesh.h5m",
189 ierr = PetscOptionsEnd();
◆ getOtherFiniteElements()
std::vector<std::string>& MoFEM::Simple::getOtherFiniteElements |
( |
| ) |
|
|
inline |
Get the Other Finite Elements.
User can create finite elements using directly core interface and and them to simple problem by this function
- Returns
- std::vector<std::string>&
- Examples
- mesh_smoothing.cpp.
Definition at line 427 of file Simple.hpp.
◆ getParentAdjacencies()
bool& MoFEM::Simple::getParentAdjacencies |
( |
| ) |
|
|
inline |
Get the addParentAdjacencies.
If set true add parent adjacencies
- Returns
- true
-
false
Definition at line 473 of file Simple.hpp.
◆ getProblemName() [1/2]
std::string& MoFEM::Simple::getProblemName |
( |
| ) |
|
|
inline |
◆ getProblemName() [2/2]
const std::string MoFEM::Simple::getProblemName |
( |
| ) |
const |
|
inline |
◆ getSkeletonFEName() [1/2]
std::string& MoFEM::Simple::getSkeletonFEName |
( |
| ) |
|
|
inline |
Get the Skeleton FE Name.
- Returns
- std::string&
Definition at line 410 of file Simple.hpp.
◆ getSkeletonFEName() [2/2]
const std::string MoFEM::Simple::getSkeletonFEName |
( |
| ) |
const |
|
inline |
◆ getSkeletonMeshSet()
Get the SkeletonMeshSet object.
- Returns
- EntityHandle&
Definition at line 347 of file Simple.hpp.
◆ loadFile() [1/2]
MoFEMErrorCode MoFEM::Simple::loadFile |
( |
const std::string |
mesh_file_name = "" | ) |
|
Load mesh file with parallel options if number of cores > 1.
- Parameters
-
mesh_file_name | file name if not set default or set by command line is used. |
- Returns
- error code
Definition at line 250 of file Simple.cpp.
253 if (m_field.get_comm_size() == 1)
257 "PARALLEL_RESOLVE_SHARED_ENTS;"
258 "PARTITION=PARALLEL_PARTITION;",
◆ loadFile() [2/2]
Load mesh file.
- Parameters
-
options | file load options |
mesh_file_name | file name if not set default or set by command line |
loadFunc | function for loading mesh file is used. |
- Note
- If bitRefLevel is set to any, bit ref level of loaded entities is not changed. After mesh is load, bit ref level should be set to create problem. Default setting of bit ref level is on first bit, and if is set all mesh entities on load are set to set level.
- Returns
- error code
- Examples
- child_and_parent.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, loop_entities.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, simple_elasticity.cpp, simple_interface.cpp, simple_l2_only.cpp, test_cache_on_entities.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.
Definition at line 194 of file Simple.cpp.
206 CHKERR m_field.rebuild_database();
211 CHKERR m_field.get_moab().get_number_entities_by_dimension(
213 if (nb_ents_3d > 0) {
217 CHKERR m_field.get_moab().get_number_entities_by_dimension(
219 if (nb_ents_2d > 0) {
238 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(ents,
bitLevel,
241 MOFEM_LOG(
"WORLD", Sev::warning) <<
"BitRefLevel is none and not set";
242 CHKERR m_field.getInterface<BitRefManager>()->addToDatabaseBitRefLevelByDim(
◆ query_interface()
◆ removeBoundaryField()
MoFEMErrorCode MoFEM::Simple::removeBoundaryField |
( |
const std::string & |
name | ) |
|
Remove field form boundary.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 424 of file Simple.cpp.
427 auto remove_field_from_list = [&](
auto &vec) {
428 auto it = std::find(vec.begin(), vec.end(), name);
◆ removeDomainField()
MoFEMErrorCode MoFEM::Simple::removeDomainField |
( |
const std::string & |
name | ) |
|
Remove field form domain.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 409 of file Simple.cpp.
412 auto remove_field_from_list = [&](
auto &vec) {
413 auto it = std::find(vec.begin(), vec.end(), name);
◆ removeSkeletonField()
MoFEMErrorCode MoFEM::Simple::removeSkeletonField |
( |
const std::string & |
name | ) |
|
Remove field form skeleton.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 438 of file Simple.cpp.
441 auto remove_field_from_list = [&](
auto &vec) {
442 auto it = std::find(vec.begin(), vec.end(), name);
◆ reSetUp()
Rebuild internal MoFEM data structures.
Call this function after you add field or remove it.
- Note
- If you add field, or remove it, finite element and problem needs to be rebuild. However DM can remain the same.
- Returns
- MoFEMErrorCode
- Examples
- child_and_parent.cpp.
Definition at line 709 of file Simple.cpp.
730 const Problem *problem_ptr;
732 const auto problem_name = problem_ptr->getName();
733 CHKERR m_field.modify_problem_ref_level_set_bit(problem_name,
bitLevel);
734 CHKERR m_field.modify_problem_mask_ref_level_set_bit(problem_name,
739 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
741 m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
◆ setDim()
void MoFEM::Simple::setDim |
( |
int |
dim | ) |
|
|
inline |
Set the problem dimension.
- Returns
- int
Definition at line 319 of file Simple.hpp.
◆ setFieldOrder()
MoFEMErrorCode MoFEM::Simple::setFieldOrder |
( |
const std::string |
field_name, |
|
|
const int |
order, |
|
|
const Range * |
ents = NULL |
|
) |
| |
Set field order.
- Parameters
-
std::field_name | field name |
order | order |
range | of entities to which order is set (If null it sat to all entities) |
- Returns
- error code
- Examples
- child_and_parent.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, loop_entities.cpp, mesh_smoothing.cpp, mixed_poisson.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, simple_elasticity.cpp, simple_interface.cpp, simple_l2_only.cpp, test_cache_on_entities.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.
Definition at line 545 of file Simple.cpp.
550 ents == NULL ?
false :
true);
◆ setParentAdjacency() [1/4]
Definition at line 72 of file Simple.cpp.
77 boost::make_shared<ParentFiniteElementAdjacencyFunction<3>>(
80 boost::make_shared<ParentFiniteElementAdjacencyFunction<2>>(
83 CHKERR m_field.modify_finite_element_adjacency_table(
domainFE, MBTET,
85 CHKERR m_field.modify_finite_element_adjacency_table(
domainFE, MBHEX,
88 CHKERR m_field.modify_finite_element_adjacency_table(
90 CHKERR m_field.modify_finite_element_adjacency_table(
94 CHKERR m_field.modify_finite_element_adjacency_table(
96 CHKERR m_field.modify_finite_element_adjacency_table(
◆ setParentAdjacency() [2/4]
Definition at line 103 of file Simple.cpp.
108 boost::make_shared<ParentFiniteElementAdjacencyFunction<2>>(
111 boost::make_shared<ParentFiniteElementAdjacencyFunction<1>>(
114 CHKERR m_field.modify_finite_element_adjacency_table(
domainFE, MBTRI,
116 CHKERR m_field.modify_finite_element_adjacency_table(
domainFE, MBQUAD,
119 CHKERR m_field.modify_finite_element_adjacency_table(
122 CHKERR m_field.modify_finite_element_adjacency_table(
◆ setParentAdjacency() [3/4]
Definition at line 128 of file Simple.cpp.
134 return setParentAdjacency<2>();
136 return setParentAdjacency<3>();
◆ setParentAdjacency() [4/4]
Definition at line 67 of file Simple.cpp.
68 static_assert(DIM == 2 || DIM == 3,
"not implemented");
◆ setSkeletonAdjacency() [1/5]
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency |
( |
int |
dim = -1 , |
|
|
std::string |
fe_name = "" |
|
) |
| |
Set the skeleton adjacency object.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 143 of file Simple.cpp.
153 return setSkeletonAdjacency<2>(fe_name);
155 return setSkeletonAdjacency<3>(fe_name);
◆ setSkeletonAdjacency() [2/5]
template<>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency |
( |
std::string |
fe_name | ) |
|
|
private |
Definition at line 21 of file Simple.cpp.
26 boost::make_shared<ParentFiniteElementAdjacencyFunctionSkeleton<1>>(
28 CHKERR m_field.modify_finite_element_adjacency_table(
◆ setSkeletonAdjacency() [3/5]
template<>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency |
( |
std::string |
fe_name | ) |
|
|
private |
Definition at line 35 of file Simple.cpp.
40 boost::make_shared<ParentFiniteElementAdjacencyFunctionSkeleton<2>>(
43 CHKERR m_field.modify_finite_element_adjacency_table(
45 CHKERR m_field.modify_finite_element_adjacency_table(
◆ setSkeletonAdjacency() [4/5]
template<>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency |
( |
std::string |
fe_name | ) |
|
|
private |
Definition at line 52 of file Simple.cpp.
59 return setSkeletonAdjacency<2>(fe_name);
61 return setSkeletonAdjacency<3>(fe_name);
◆ setSkeletonAdjacency() [5/5]
template<int DIM>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency |
( |
std::string |
fe_name | ) |
|
|
private |
Definition at line 15 of file Simple.cpp.
16 static_assert(DIM == 2 || DIM == 3,
"not implemented");
◆ setUp()
MoFEMErrorCode MoFEM::Simple::setUp |
( |
const PetscBool |
is_partitioned = PETSC_TRUE | ) |
|
Setup problem.
- Returns
- error code
- Examples
- child_and_parent.cpp, continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, field_evaluator.cpp, hcurl_check_approx_in_2d.cpp, hcurl_curl_operator.cpp, hdiv_divergence_operator.cpp, loop_entities.cpp, mixed_poisson.cpp, poisson_2d_dis_galerkin.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, simple_interface.cpp, simple_l2_only.cpp, test_cache_on_entities.cpp, testing_jacobian_of_hook_element.cpp, and testing_jacobian_of_hook_scaled_with_density_element.cpp.
Definition at line 683 of file Simple.cpp.
◆ addBoundaryFE
bool MoFEM::Simple::addBoundaryFE |
|
private |
◆ addParentAdjacencies
bool MoFEM::Simple::addParentAdjacencies |
|
private |
If set to true parent adjacencies are build.
Definition at line 536 of file Simple.hpp.
◆ addSkeletonFE
bool MoFEM::Simple::addSkeletonFE |
|
private |
◆ bitAdjEnt
bit ref level for parent
Definition at line 540 of file Simple.hpp.
◆ bitAdjEntMask
bit ref level for parent parent
Definition at line 541 of file Simple.hpp.
◆ bitAdjParent
bit ref level for parent
Definition at line 538 of file Simple.hpp.
◆ bitAdjParentMask
bit ref level for parent parent
Definition at line 539 of file Simple.hpp.
◆ bitLevel
BitRefLevel of the problem.
Definition at line 520 of file Simple.hpp.
◆ bitLevelMask
BitRefLevel of the problem.
Definition at line 521 of file Simple.hpp.
◆ boundaryFE
std::string MoFEM::Simple::boundaryFE |
|
private |
boundary finite element
Definition at line 556 of file Simple.hpp.
◆ boundaryFields
std::vector<std::string> MoFEM::Simple::boundaryFields |
|
private |
◆ boundaryMeshset
◆ cOre
◆ dataFields
std::vector<std::string> MoFEM::Simple::dataFields |
|
private |
◆ dIm
◆ dM
Discrete manager (interface to PETSc/MoFEM functions)
Definition at line 565 of file Simple.hpp.
◆ domainFE
std::string MoFEM::Simple::domainFE |
|
private |
◆ domainFields
std::vector<std::string> MoFEM::Simple::domainFields |
|
private |
◆ fieldsOrder
std::list<std::tuple<std::string, int, Range, bool> > MoFEM::Simple::fieldsOrder |
|
private |
fields order. 1: field name, order, range, set by range if true
Definition at line 551 of file Simple.hpp.
◆ meshFileName
char MoFEM::Simple::meshFileName[255] |
|
private |
◆ meshSet
◆ MOFEM_EVENT_SimpleBuildFields
PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleBuildFields |
|
private |
◆ MOFEM_EVENT_SimpleBuildFiniteElements
PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleBuildFiniteElements |
|
private |
◆ MOFEM_EVENT_SimpleBuildProblem
PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleBuildProblem |
|
private |
◆ MOFEM_EVENT_SimpleKSPSolve
PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleKSPSolve |
|
private |
◆ MOFEM_EVENT_SimpleLoadMesh
PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleLoadMesh |
|
private |
◆ MOFEM_EVENT_SimpleSetUP
PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleSetUP |
|
private |
◆ nameOfProblem
std::string MoFEM::Simple::nameOfProblem |
|
private |
◆ noFieldDataFields
std::vector<std::string> MoFEM::Simple::noFieldDataFields |
|
private |
◆ noFieldFields
std::vector<std::string> MoFEM::Simple::noFieldFields |
|
private |
◆ otherFEs
std::vector<std::string> MoFEM::Simple::otherFEs |
|
private |
Other finite elements.
Definition at line 559 of file Simple.hpp.
◆ parentAdjFunctionDim1
◆ parentAdjFunctionDim2
◆ parentAdjFunctionDim3
◆ parentAdjSkeletonFunctionDim1
◆ parentAdjSkeletonFunctionDim2
◆ skeletonFE
std::string MoFEM::Simple::skeletonFE |
|
private |
skeleton finite element
Definition at line 557 of file Simple.hpp.
◆ skeletonFields
std::vector<std::string> MoFEM::Simple::skeletonFields |
|
private |
fields on the skeleton
Definition at line 545 of file Simple.hpp.
◆ skeletonMeshset
skeleton meshset with boundary
Definition at line 532 of file Simple.hpp.
The documentation for this struct was generated from the following files:
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
std::vector< std::string > boundaryFields
boundary fields
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 1 > > parentAdjFunctionDim1
MoFEMErrorCode setSkeletonAdjacency(int dim=-1, std::string fe_name="")
Set the skeleton adjacency object.
#define MYPCOMM_INDEX
default communicator number PCOMM
std::vector< std::string > otherFEs
Other finite elements.
bool addParentAdjacencies
If set to true parent adjacencies are build.
bool addBoundaryFE
Add boundary FE.
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
MoFEMErrorCode exchangeGhostCells()
std::list< std::tuple< std::string, int, Range, bool > > fieldsOrder
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Simple(const MoFEM::Core &core)
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 2 > > parentAdjFunctionDim2
std::vector< std::string > noFieldFields
NOFIELD field name.
multi_index_container< DofsSideMapData, indexed_by< ordered_non_unique< tag< TypeSide_mi_tag >, composite_key< DofsSideMapData, member< DofsSideMapData, EntityType, &DofsSideMapData::type >, member< DofsSideMapData, int, &DofsSideMapData::side > >>, ordered_unique< tag< EntDofIdx_mi_tag >, member< DofsSideMapData, int, &DofsSideMapData::dof > > > > DofsSideMap
Map entity stype and side to element/entity dof index.
EntityHandle meshSet
domain meshset
MoFEMErrorCode createBoundaryMeshset()
char meshFileName[255]
mesh file name
PetscLogEvent MOFEM_EVENT_SimpleKSPSolve
MoFEMErrorCode buildProblem()
Build problem.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
PetscLogEvent MOFEM_EVENT_SimpleBuildFields
PetscObject getPetscObject(T obj)
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
PetscErrorCode DMSetUp_MoFEM(DM dm)
int dIm
dimension of problem
bool addSkeletonFE
Add skeleton FE.
std::string nameOfProblem
problem name
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 2 > > parentAdjSkeletonFunctionDim2
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
DeprecatedCoreInterface Interface
PetscLogEvent MOFEM_EVENT_SimpleLoadMesh
BitRefLevel bitLevelMask
BitRefLevel of the problem.
std::vector< std::string > dataFields
Data fields.
#define CHKERR
Inline error check.
BitRefLevel bitAdjEntMask
bit ref level for parent parent
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
std::string domainFE
domain finite element
PetscLogEvent MOFEM_EVENT_SimpleBuildFiniteElements
MoFEMErrorCode createSkeletonMeshset()
static MoFEMErrorCode setDofsSideMap(const FieldSpace space, const FieldContinuity continuity, const FieldApproximationBase base, DofsSideMap &)
Set map of dof to side number.
std::string boundaryFE
boundary finite element
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 3 > > parentAdjFunctionDim3
BitRefLevel bitAdjParentMask
bit ref level for parent parent
std::string skeletonFE
skeleton finite element
@ MOFEM_OPERATION_UNSUCCESSFUL
BitRefLevel bitAdjEnt
bit ref level for parent
constexpr double t
plate stiffness
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
constexpr auto field_name
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
MoFEMErrorCode buildFiniteElements()
Build finite elements.
static MoFEMErrorCode setDofsSideMap(const FieldSpace space, const FieldContinuity continuity, const FieldApproximationBase base, DofsSideMap &)
Set map of dof to side number.
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define MOFEM_LOG(channel, severity)
Log.
std::vector< std::string > domainFields
domain fields
@ DISCONTINUOUS
Broken continuity (No effect on L2 space)
std::vector< std::string > skeletonFields
fields on the skeleton
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
@ HCURL
field with continuous tangents
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
BitRefLevel bitLevel
BitRefLevel of the problem.
MoFEMErrorCode defineFiniteElements()
Define finite elements.
int getDim() const
Get the problem dimension.
MoFEMErrorCode buildFields()
Build fields.
BitRefLevel bitAdjParent
bit ref level for parent
std::vector< std::string > noFieldDataFields
NOFIELD field name.
multi_index_container< boost::shared_ptr< FiniteElement >, indexed_by< hashed_unique< tag< FiniteElement_Meshset_mi_tag >, member< FiniteElement, EntityHandle, &FiniteElement::meshset > >, hashed_unique< tag< BitFEId_mi_tag >, const_mem_fun< FiniteElement, BitFEId, &FiniteElement::getId >, HashBit< BitFEId >, EqBit< BitFEId > >, ordered_unique< tag< FiniteElement_name_mi_tag >, const_mem_fun< FiniteElement, boost::string_ref, &FiniteElement::getNameRef > > > > FiniteElement_multiIndex
MultiIndex for entities for FiniteElement.
EntityHandle boundaryMeshset
meshset with boundary
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 1 > > parentAdjSkeletonFunctionDim1
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
SmartPetscObj< DM > dM
Discrete manager (interface to PETSc/MoFEM functions)
MoFEMErrorCode setParentAdjacency()
const std::string getProblemName() const
Get the Problem Name.
MoFEMErrorCode addFieldToEmptyFieldBlocks(const std::string row_field, const std::string col_field) const
add empty block to problem
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
EntityHandle skeletonMeshset
skeleton meshset with boundary
PetscLogEvent MOFEM_EVENT_SimpleSetUP
@ HDIV
field with continuous normal traction
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ NOFIELD
scalar or vector of scalars describe (no true field)
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
PetscLogEvent MOFEM_EVENT_SimpleBuildProblem