v0.14.0
Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
MoFEM::Simple Struct Reference

Simple interface for fast problem set-up. More...

#include <src/interfaces/Simple.hpp>

Inheritance diagram for MoFEM::Simple:
[legend]
Collaboration diagram for MoFEM::Simple:
[legend]

Public Types

typedef boost::function< MoFEMErrorCode(Interface &, const char *, const char *)> LoadFileFunc
 

Public Member Functions

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 EntityHandlegetMeshSet ()
 
EntityHandlegetMeshset ()
 Get the MeshSet object. More...
 
EntityHandlegetBoundaryMeshSet ()
 Get the BoundaryMeshSet object. More...
 
EntityHandlegetSkeletonMeshSet ()
 Get the SkeletonMeshSet object. More...
 
BitRefLevelgetBitRefLevel ()
 Get the BitRefLevel. More...
 
BitRefLevelgetBitRefLevelMask ()
 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...
 
boolgetAddSkeletonFE ()
 Get the addSkeletonFE. More...
 
boolgetAddBoundaryFE ()
 Get the addSkeletonFE. More...
 
boolgetParentAdjacencies ()
 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...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
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
 

Static Public Member Functions

static MoFEMErrorCode defaultLoadFileFunc (Interface &m_field, const char *file_name, const char *options)
 Set function for loading mesh file. More...
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 

Private Member Functions

MoFEMErrorCode createBoundaryMeshset ()
 
MoFEMErrorCode createSkeletonMeshset ()
 
MoFEMErrorCode exchangeGhostCells ()
 
template<int DIM = -1>
MoFEMErrorCode setSkeletonAdjacency (std::string fe_name)
 
template<int DIM = -1>
MoFEMErrorCode setParentAdjacency ()
 
template<>
MoFEMErrorCode setSkeletonAdjacency (std::string fe_name)
 
template<>
MoFEMErrorCode setSkeletonAdjacency (std::string fe_name)
 
template<>
MoFEMErrorCode setSkeletonAdjacency (std::string fe_name)
 
template<>
MoFEMErrorCode setParentAdjacency ()
 
template<>
MoFEMErrorCode setParentAdjacency ()
 
template<>
MoFEMErrorCode setParentAdjacency ()
 

Private Attributes

MoFEM::CorecOre
 
BitRefLevel bitLevel
 BitRefLevel of the problem. More...
 
BitRefLevel bitLevelMask
 BitRefLevel of the problem. More...
 
PetscLogEvent MOFEM_EVENT_SimpleSetUP
 
PetscLogEvent MOFEM_EVENT_SimpleLoadMesh
 
PetscLogEvent MOFEM_EVENT_SimpleBuildFields
 
PetscLogEvent MOFEM_EVENT_SimpleBuildFiniteElements
 
PetscLogEvent MOFEM_EVENT_SimpleBuildProblem
 
PetscLogEvent MOFEM_EVENT_SimpleKSPSolve
 
EntityHandle meshSet
 domain meshset More...
 
EntityHandle boundaryMeshset
 meshset with boundary More...
 
EntityHandle skeletonMeshset
 skeleton meshset with boundary More...
 
bool addSkeletonFE
 Add skeleton FE. More...
 
bool addBoundaryFE
 Add boundary FE. More...
 
bool addParentAdjacencies
 If set to true parent adjacencies are build. More...
 
BitRefLevel bitAdjParent
 bit ref level for parent More...
 
BitRefLevel bitAdjParentMask
 bit ref level for parent parent More...
 
BitRefLevel bitAdjEnt
 bit ref level for parent More...
 
BitRefLevel bitAdjEntMask
 bit ref level for parent parent More...
 
std::vector< std::string > domainFields
 domain fields More...
 
std::vector< std::string > boundaryFields
 boundary fields More...
 
std::vector< std::string > skeletonFields
 fields on the skeleton More...
 
std::vector< std::string > dataFields
 Data fields. More...
 
std::vector< std::string > noFieldFields
 NOFIELD field name. More...
 
std::vector< std::string > noFieldDataFields
 NOFIELD field name. More...
 
std::list< std::tuple< std::string, int, Range, bool > > fieldsOrder
 
std::string nameOfProblem
 problem name More...
 
std::string domainFE
 domain finite element More...
 
std::string boundaryFE
 boundary finite element More...
 
std::string skeletonFE
 skeleton finite element More...
 
std::vector< std::string > otherFEs
 Other finite elements. More...
 
char meshFileName [255]
 mesh file name More...
 
int dIm
 dimension of problem More...
 
SmartPetscObj< DM > dM
 Discrete manager (interface to PETSc/MoFEM functions) More...
 
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 3 > > parentAdjFunctionDim3
 
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 2 > > parentAdjFunctionDim2
 
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 1 > > parentAdjFunctionDim1
 
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 2 > > parentAdjSkeletonFunctionDim2
 
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 1 > > parentAdjSkeletonFunctionDim1
 

Detailed Description

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_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, nonlinear_elastic.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.

Member Typedef Documentation

◆ LoadFileFunc

typedef boost::function<MoFEMErrorCode(Interface &, const char *, const char *)> MoFEM::Simple::LoadFileFunc

Definition at line 43 of file Simple.hpp.

Constructor & Destructor Documentation

◆ Simple()

MoFEM::Simple::Simple ( const MoFEM::Core core)

Definition at line 162 of file Simple.cpp.

163  : cOre(const_cast<Core &>(core)), bitLevel(BitRefLevel().set(0)),
165  skeletonMeshset(0), nameOfProblem("SimpleProblem"), domainFE("dFE"),
166  boundaryFE("bFE"), skeletonFE("sFE"), dIm(-1), addSkeletonFE(false),
167  addBoundaryFE(false), addParentAdjacencies(false),
169  bitAdjEnt(BitRefLevel().set()), bitAdjEntMask(BitRefLevel().set()) {
170  PetscLogEventRegister("SimpSetUp", 0, &MOFEM_EVENT_SimpleSetUP);
171  PetscLogEventRegister("SimpLoadMesh", 0, &MOFEM_EVENT_SimpleLoadMesh);
172  PetscLogEventRegister("SimpBuildFields", 0, &MOFEM_EVENT_SimpleBuildFields);
173  PetscLogEventRegister("SimpBuildFEs", 0,
175  PetscLogEventRegister("SimpSetUp", 0, &MOFEM_EVENT_SimpleBuildProblem);
176  PetscLogEventRegister("SimpKSPSolve", 0, &MOFEM_EVENT_SimpleKSPSolve);
177  strcpy(meshFileName, "mesh.h5m");
178 }

◆ ~Simple()

virtual MoFEM::Simple::~Simple ( )
virtualdefault

Member Function Documentation

◆ addBoundaryField()

MoFEMErrorCode MoFEM::Simple::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.

Parameters
namename of the filed
spacespace (L2,H1,Hdiv,Hcurl)
baseapproximation base, see FieldApproximationBase
nb_of_coefficientsnumber of field coefficients
tag_typetype 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)
bhif MF_EXCL throws error if field exits, MF_ZERO no error if field exist
verbverbosity 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 300 of file Simple.cpp.

304  {
305  Interface &m_field = cOre;
307  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
308  verb);
309  boundaryFields.push_back(name);
310  if (space == NOFIELD)
311  SETERRQ(
312  PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
313  "NOFIELD space for boundary filed not implemented in Simple interface");
315 }

◆ addDataField()

MoFEMErrorCode MoFEM::Simple::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.

Parameters
namename of the filed
spacespace (L2,H1,Hdiv,Hcurl)
baseapproximation base, see FieldApproximationBase
nb_of_coefficientsnumber of field coefficients
tag_typetype 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)
bhif MF_EXCL throws error if field exits, MF_ZERO no error if field exist
verbverbosity level
Returns
error code
Examples
hcurl_curl_operator.cpp, and hdiv_divergence_operator.cpp.

Definition at line 338 of file Simple.cpp.

342  {
343 
344  Interface &m_field = cOre;
346  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
347  verb);
348  if (space == NOFIELD)
349  noFieldDataFields.push_back(name);
350  else
351  dataFields.push_back(name);
353 }

◆ addDomainBrokenField()

MoFEMErrorCode MoFEM::Simple::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.

Parameters
namename of the filed
spacespace (L2,H1,Hdiv,Hcurl)
baseapproximation base, see FieldApproximationBase
nb_of_coefficientsnumber of field coefficients
tag_typetype 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)
bhif MF_EXCL throws error if field exits, MF_ZERO no error if field exist
verbverbosity level
Returns
error code

Definition at line 282 of file Simple.cpp.

286  {
287 
288  Interface &m_field = cOre;
290  CHKERR m_field.add_broken_field(name, space, base, nb_of_cooficients,
291  tag_type, bh, verb);
292  if (space == NOFIELD)
293  noFieldFields.push_back(name);
294  else
295  domainFields.push_back(name);
297 }

◆ addDomainField()

MoFEMErrorCode MoFEM::Simple::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.

Parameters
namename of the filed
spacespace (L2,H1,Hdiv,Hcurl)
baseapproximation base, see FieldApproximationBase
nb_of_coefficientsnumber of field coefficients
tag_typetype 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)
bhif MF_EXCL throws error if field exits, MF_ZERO no error if field exist
verbverbosity 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.

268  {
269 
270  Interface &m_field = cOre;
272  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
273  verb);
274  if (space == NOFIELD)
275  noFieldFields.push_back(name);
276  else
277  domainFields.push_back(name);
279 }

◆ 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_fieldrow filed name
col_fieldcol field name
Returns
MoFEMErrorCode

Definition at line 728 of file Simple.cpp.

729  {
730  Interface &m_field = cOre;
732  CHKERR m_field.getInterface<ProblemsManager>()->addFieldToEmptyFieldBlocks(
733  getProblemName(), row_field, col_field);
735 }

◆ addSkeletonField()

MoFEMErrorCode MoFEM::Simple::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.

Parameters
namename of the filed
spacespace (L2,H1,Hdiv,Hcurl)
baseapproximation base, see FieldApproximationBase
nb_of_coefficientsnumber of field coefficients
tag_typetype 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)
bhif MF_EXCL throws error if field exits, MF_ZERO no error if field exist
verbverbosity 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 318 of file Simple.cpp.

322  {
323 
324  Interface &m_field = cOre;
326  CHKERR m_field.add_field(name, space, base, nb_of_cooficients, tag_type, bh,
327  verb);
328  skeletonFields.push_back(name);
329  if (space == NOFIELD)
330  SETERRQ(
331  PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
332  "NOFIELD space for boundary filed not implemented in Simple interface");
333 
335 }

◆ buildFields()

MoFEMErrorCode MoFEM::Simple::buildFields ( )

Build fields.

Returns
error code
Examples
mesh_smoothing.cpp, and simple_elasticity.cpp.

Definition at line 500 of file Simple.cpp.

500  {
501  Interface &m_field = cOre;
503  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildFields, 0, 0, 0, 0);
504 
505  auto comm_interface_ptr = m_field.getInterface<CommInterface>();
506  auto bit_ref_ptr = m_field.getInterface<BitRefManager>();
507 
508  // Add entities to the fields
509  auto add_ents_to_field = [&](auto meshset, auto dim, auto &fields) {
511  for (auto &field : fields) {
512  CHKERR m_field.add_ents_to_field_by_dim(meshset, dim, field);
513  CHKERR comm_interface_ptr->synchroniseFieldEntities(field, 0);
514  }
516  };
517 
518  auto make_no_field_ents = [&](auto &fields) {
520  for (auto &field : fields) {
521  std::array<double, 6> coords = {0, 0, 0, 0, 0, 0};
522  CHKERR m_field.create_vertices_and_add_to_field(field, coords.data(), 1);
523  CHKERR comm_interface_ptr->makeFieldEntitiesMultishared(field, 0);
524  CHKERR bit_ref_ptr->setFieldEntitiesBitRefLevel(field, bitLevel);
525  }
527  };
528 
529  CHKERR add_ents_to_field(meshSet, dIm, domainFields);
530  CHKERR add_ents_to_field(meshSet, dIm, dataFields);
531  CHKERR add_ents_to_field(boundaryMeshset, dIm - 1, boundaryFields);
532  CHKERR add_ents_to_field(skeletonMeshset, dIm - 1, skeletonFields);
533 
534  std::set<std::string> nofield_fields;
535  for (auto &f : noFieldFields)
536  nofield_fields.insert(f);
537  for (auto &f : noFieldDataFields)
538  nofield_fields.insert(f);
539 
540  CHKERR make_no_field_ents(nofield_fields);
541 
542  // Set order
543 
544  auto get_field_ptr = [&](auto &f) { return m_field.get_field_structure(f); };
545  for (auto &t : fieldsOrder) {
546  const auto f = std::get<0>(t);
547  const auto order = std::get<1>(t);
548 
549  MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "Simple")
550  << "Set order to field " << f << " order " << order;
551  MOFEM_LOG_CHANNEL("SYNC");
552  if (std::get<3>(t)) {
553  MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "Simple")
554  << "To ents: " << std::endl
555  << std::get<2>(t) << std::endl;
556  }
557  MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::verbose);
558 
559  if (std::get<3>(t)) {
560 
561  CHKERR m_field.set_field_order(std::get<2>(t), f, order);
562 
563  } else {
564  auto f_ptr = get_field_ptr(f);
565 
566  if (f_ptr->getSpace() == H1) {
567  if (f_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
568  CHKERR m_field.set_field_order(meshSet, MBVERTEX, f, order);
569  } else {
570  CHKERR m_field.set_field_order(meshSet, MBVERTEX, f, 1);
571  }
572  }
573 
574  for (auto d = 1; d <= dIm; ++d) {
575  for (EntityType t = CN::TypeDimensionMap[d].first;
576  t <= CN::TypeDimensionMap[d].second; ++t) {
577  CHKERR m_field.set_field_order(meshSet, t, f, order);
578  }
579  }
580  }
581  }
582  MOFEM_LOG_CHANNEL("WORLD");
583  // Build fields
584  CHKERR m_field.build_fields();
585  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildFields, 0, 0, 0, 0);
587 }

◆ buildFiniteElements()

MoFEMErrorCode MoFEM::Simple::buildFiniteElements ( )

Build finite elements.

Returns
error code
Examples
mesh_smoothing.cpp, and simple_elasticity.cpp.

Definition at line 589 of file Simple.cpp.

589  {
590  Interface &m_field = cOre;
592  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildFiniteElements, 0, 0, 0, 0);
593  // Add finite elements
594  CHKERR m_field.add_ents_to_finite_element_by_dim(meshSet, dIm, domainFE,
595  true);
596  CHKERR m_field.build_finite_elements(domainFE);
597  if (addBoundaryFE || boundaryFields.size()) {
598  CHKERR m_field.add_ents_to_finite_element_by_dim(boundaryMeshset, dIm - 1,
599  boundaryFE, true);
600  CHKERR m_field.build_finite_elements(boundaryFE);
601  }
602  if (addSkeletonFE || skeletonFields.size()) {
603  CHKERR m_field.add_ents_to_finite_element_by_dim(skeletonMeshset, dIm - 1,
604  skeletonFE, true);
605  CHKERR m_field.build_finite_elements(skeletonFE);
606  }
607  for (std::vector<std::string>::iterator fit = otherFEs.begin();
608  fit != otherFEs.end(); ++fit) {
609  CHKERR m_field.build_finite_elements(*fit);
610  }
611  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildFiniteElements, 0, 0, 0, 0);
613 }

◆ buildProblem()

MoFEMErrorCode MoFEM::Simple::buildProblem ( )

Build problem.

Returns
error code
Examples
mesh_smoothing.cpp, and simple_elasticity.cpp.

Definition at line 615 of file Simple.cpp.

615  {
616  Interface &m_field = cOre;
618  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
619  CHKERR m_field.build_adjacencies(bitLevel, bitLevelMask);
620  // Set problem by the DOFs on the fields rather that by adding DOFs on the
621  // elements
622  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
624  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
625  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
627 }

◆ createBoundaryMeshset()

MoFEMErrorCode MoFEM::Simple::createBoundaryMeshset ( )
private

Definition at line 737 of file Simple.cpp.

737  {
738  Interface &m_field = cOre;
740  ParallelComm *pcomm =
741  ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
742 
743  auto get_skin = [&](auto meshset) {
744  // filter not owned entities, those are not on boundary
745 
746  Range domain_ents;
747  CHKERR m_field.get_moab().get_entities_by_dimension(meshset, dIm,
748  domain_ents, true);
749  CHKERR pcomm->filter_pstatus(domain_ents,
750  PSTATUS_SHARED | PSTATUS_MULTISHARED,
751  PSTATUS_NOT, -1, nullptr);
752 
753  Skinner skin(&m_field.get_moab());
754  Range domain_skin;
755  CHKERR skin.find_skin(0, domain_ents, false, domain_skin);
756  CHKERR pcomm->filter_pstatus(domain_skin,
757  PSTATUS_SHARED | PSTATUS_MULTISHARED,
758  PSTATUS_NOT, -1, nullptr);
759  return domain_skin;
760  };
761 
762  auto create_boundary_meshset = [&](auto &&domain_skin) {
764  // create boundary meshset
765  if (boundaryMeshset != 0) {
767  }
768  CHKERR m_field.get_moab().create_meshset(MESHSET_SET, boundaryMeshset);
769  CHKERR m_field.get_moab().add_entities(boundaryMeshset, domain_skin);
770  for (int dd = 0; dd != dIm - 1; dd++) {
771  Range adj;
772  CHKERR m_field.get_moab().get_adjacencies(domain_skin, dd, false, adj,
773  moab::Interface::UNION);
774  CHKERR m_field.get_moab().add_entities(boundaryMeshset, adj);
775  }
777  };
778 
779  CHKERR create_boundary_meshset(get_skin(meshSet));
780 
782 }

◆ createSkeletonMeshset()

MoFEMErrorCode MoFEM::Simple::createSkeletonMeshset ( )
private

Definition at line 784 of file Simple.cpp.

784  {
785  Interface &m_field = cOre;
787 
788  ParallelComm *pcomm =
789  ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
790 
791  auto create_skeleton_meshset = [&](auto meshset) {
793  // create boundary meshset
794  if (skeletonMeshset != 0) {
796  }
797  Range boundary_ents, skeleton_ents;
798  CHKERR m_field.get_moab().get_entities_by_dimension(boundaryMeshset,
799  dIm - 1, boundary_ents);
800  Range domain_ents;
801  CHKERR m_field.get_moab().get_entities_by_dimension(meshset, dIm,
802  domain_ents, true);
803  CHKERR m_field.get_moab().get_adjacencies(
804  domain_ents, dIm - 1, false, skeleton_ents, moab::Interface::UNION);
805  skeleton_ents = subtract(skeleton_ents, boundary_ents);
806  CHKERR pcomm->filter_pstatus(skeleton_ents, PSTATUS_NOT_OWNED, PSTATUS_NOT,
807  -1, nullptr);
808  CHKERR m_field.get_moab().create_meshset(MESHSET_SET, skeletonMeshset);
809  CHKERR m_field.get_moab().add_entities(skeletonMeshset, skeleton_ents);
811  };
812 
813  CHKERR create_skeleton_meshset(meshSet);
814 
816 }

◆ 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_fieldinterface
file_namefile name
optionsloader options
Returns
error code

Definition at line 51 of file Simple.hpp.

53  {
54  // Default behavior using the member moab.load_file
55  return m_field.get_moab().load_file(file_name, 0, options);
56  }

◆ defineFiniteElements()

MoFEMErrorCode MoFEM::Simple::defineFiniteElements ( )

Define finite elements.

Returns
Error code
Examples
mesh_smoothing.cpp, and simple_elasticity.cpp.

Definition at line 398 of file Simple.cpp.

398  {
399  Interface &m_field = cOre;
401 
402  auto clear_rows_and_cols = [&](auto &fe_name) {
404  auto fe_ptr = m_field.get_finite_elements();
405  auto &fe_by_name = const_cast<FiniteElement_multiIndex *>(fe_ptr)
406  ->get<FiniteElement_name_mi_tag>();
407  auto it_fe = fe_by_name.find(fe_name);
408  if (it_fe != fe_by_name.end()) {
409 
410  if (!fe_by_name.modify(it_fe, FiniteElement_row_change_bit_reset()))
411  SETERRQ(m_field.get_comm(), MOFEM_OPERATION_UNSUCCESSFUL,
412  "modification unsuccessful");
413 
414  if (!fe_by_name.modify(it_fe, FiniteElement_col_change_bit_reset()))
415  SETERRQ(m_field.get_comm(), MOFEM_OPERATION_UNSUCCESSFUL,
416  "modification unsuccessful");
417  }
419  };
420  CHKERR clear_rows_and_cols(domainFE);
421  CHKERR clear_rows_and_cols(boundaryFE);
422  CHKERR clear_rows_and_cols(skeletonFE);
423 
424  // Define finite elements
425  CHKERR m_field.add_finite_element(domainFE, MF_ZERO);
426 
427  auto add_fields = [&](auto &fe_name, auto &fields) {
429  for (auto &field : fields) {
430  CHKERR m_field.modify_finite_element_add_field_row(fe_name, field);
431  CHKERR m_field.modify_finite_element_add_field_col(fe_name, field);
432  CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
433  }
435  };
436 
437  auto add_data_fields = [&](auto &fe_name, auto &fields) {
439  for (auto &field : fields)
440  CHKERR m_field.modify_finite_element_add_field_data(fe_name, field);
442  };
443 
444  CHKERR add_fields(domainFE, domainFields);
445  CHKERR add_data_fields(domainFE, dataFields);
446  CHKERR add_fields(domainFE, noFieldFields);
447  CHKERR add_data_fields(domainFE, noFieldDataFields);
448 
449  if (addBoundaryFE || !boundaryFields.empty()) {
450  CHKERR m_field.add_finite_element(boundaryFE, MF_ZERO);
451  CHKERR add_fields(boundaryFE, domainFields);
452  if (!boundaryFields.empty())
453  CHKERR add_fields(boundaryFE, boundaryFields);
454  CHKERR add_data_fields(boundaryFE, dataFields);
455  CHKERR add_data_fields(boundaryFE, noFieldDataFields);
456  CHKERR add_fields(boundaryFE, noFieldFields);
457  }
458  if (addSkeletonFE || !skeletonFields.empty()) {
459  CHKERR m_field.add_finite_element(skeletonFE, MF_ZERO);
460  CHKERR add_fields(skeletonFE, domainFields);
461  if (!skeletonFields.empty())
462  CHKERR add_fields(skeletonFE, skeletonFields);
463  CHKERR add_data_fields(skeletonFE, dataFields);
464  CHKERR add_data_fields(skeletonFE, noFieldDataFields);
465  CHKERR add_fields(skeletonFE, noFieldFields);
466  }
468 }

◆ defineProblem()

MoFEMErrorCode MoFEM::Simple::defineProblem ( const PetscBool  is_partitioned = PETSC_TRUE)

define problem

Returns
error code
Examples
mesh_smoothing.cpp, and simple_elasticity.cpp.

Definition at line 470 of file Simple.cpp.

470  {
471  Interface &m_field = cOre;
473  // Create dm instance
474  dM = createDM(m_field.get_comm(), "DMMOFEM");
475  // set dm data structure which created mofem data structures
476  CHKERR DMMoFEMCreateMoFEM(dM, &m_field, nameOfProblem.c_str(), bitLevel,
477  bitLevelMask);
478  CHKERR DMSetFromOptions(dM);
480  if (addBoundaryFE || !boundaryFields.empty()) {
482  }
483  if (addSkeletonFE || !skeletonFields.empty()) {
485  }
487  CHKERR DMMoFEMSetIsPartitioned(dM, is_partitioned);
489 }

◆ deleteDM()

MoFEMErrorCode MoFEM::Simple::deleteDM ( )

Delete dm.

Returns
MoFEMErrorCode

Definition at line 705 of file Simple.cpp.

705  {
707  dM.reset();
709 }

◆ deleteFiniteElements()

MoFEMErrorCode MoFEM::Simple::deleteFiniteElements ( )

Delete finite elements.

Returns
MoFEMErrorCode

Definition at line 716 of file Simple.cpp.

716  {
717  Interface &m_field = cOre;
719  for (auto fe : {domainFE, boundaryFE, skeletonFE}) {
720  if (m_field.check_finite_element(fe)) {
721  CHKERR m_field.delete_finite_element(fe);
722  }
723  }
725 }

◆ exchangeGhostCells()

MoFEMErrorCode MoFEM::Simple::exchangeGhostCells ( )
private

get all adjacent ghosted entities

Definition at line 818 of file Simple.cpp.

818  {
819  Interface &m_field = cOre;
821  MOFEM_LOG("WORLD", Sev::verbose) << "Exchange ghost cells";
822 
823  ParallelComm *pcomm =
824  ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
825  if (pcomm == NULL)
826  pcomm = new ParallelComm(&m_field.get_moab(), m_field.get_comm());
827 
828  CHKERR pcomm->exchange_ghost_cells(getDim(), getDim() - 1, 1,
829  3 /**get all adjacent ghosted entities */,
830  true, true, meshSet ? &meshSet : nullptr);
831 
832  Range shared;
833  CHKERR m_field.get_moab().get_entities_by_dimension(0, dIm, shared);
834  for (auto d = dIm - 1; d >= 0; --d) {
835  CHKERR m_field.get_moab().get_adjacencies(shared, d, false, shared,
836  moab::Interface::UNION);
837  }
838  CHKERR pcomm->filter_pstatus(shared, PSTATUS_SHARED | PSTATUS_MULTISHARED,
839  PSTATUS_OR, -1, &shared);
840  Tag part_tag = pcomm->part_tag();
841  CHKERR pcomm->exchange_tags(part_tag, shared);
842  CHKERR m_field.getInterface<CommInterface>()->resolveParentEntities(shared,
843  VERBOSE);
844 
846 }

◆ getAddBoundaryFE()

bool& MoFEM::Simple::getAddBoundaryFE ( )
inline

Get the addSkeletonFE.

If variable set to true, boundary element is created regardless field on skelton is added or not.

Returns
true
false
Examples
poisson_2d_dis_galerkin.cpp, and simple_l2_only.cpp.

Definition at line 463 of file Simple.hpp.

463 { return addBoundaryFE; }

◆ getAddSkeletonFE()

bool& MoFEM::Simple::getAddSkeletonFE ( )
inline

Get the addSkeletonFE.

If variable set to true, skeleton element is created regardless field on skelton is added or not.

Returns
true
false
Examples
continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, poisson_2d_dis_galerkin.cpp, and simple_l2_only.cpp.

Definition at line 452 of file Simple.hpp.

452 { return addSkeletonFE; }

◆ getBitAdjEnt()

auto& MoFEM::Simple::getBitAdjEnt ( )
inline

bit ref level for parent

Returns
auto&

Definition at line 494 of file Simple.hpp.

494 { return bitAdjEnt; }

◆ getBitAdjEntMask()

auto& MoFEM::Simple::getBitAdjEntMask ( )
inline

bit ref level for parent parent

Returns
auto&

Definition at line 501 of file Simple.hpp.

501 { return bitAdjEntMask; }

◆ getBitAdjParent()

auto& MoFEM::Simple::getBitAdjParent ( )
inline

bit ref level for parent

Returns
auto&

Definition at line 480 of file Simple.hpp.

480 { return bitAdjParent; }

◆ getBitAdjParentMask()

auto& MoFEM::Simple::getBitAdjParentMask ( )
inline

bit ref level for parent parent

Returns
auto&

Definition at line 487 of file Simple.hpp.

487 { return bitAdjParentMask; }

◆ getBitRefLevel()

BitRefLevel& MoFEM::Simple::getBitRefLevel ( )
inline

Get the BitRefLevel.

Returns
BitRefLevel
Examples
child_and_parent.cpp, level_set.cpp, and mixed_poisson.cpp.

Definition at line 354 of file Simple.hpp.

354 { return bitLevel; }

◆ getBitRefLevelMask()

BitRefLevel& MoFEM::Simple::getBitRefLevelMask ( )
inline

Get the BitRefLevel.

Returns
BitRefLevel
Examples
child_and_parent.cpp.

Definition at line 361 of file Simple.hpp.

361 { return bitLevelMask; }

◆ 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.

403 { return boundaryFE; }

◆ getBoundaryFEName() [2/2]

const std::string MoFEM::Simple::getBoundaryFEName ( ) const
inline

Get the Boundary FE Name.

Returns
const std::string
Examples
mesh_smoothing.cpp, and simple_interface.cpp.

Definition at line 375 of file Simple.hpp.

375 { return boundaryFE; }

◆ getBoundaryMeshSet()

EntityHandle& MoFEM::Simple::getBoundaryMeshSet ( )
inline

Get the BoundaryMeshSet object.

Returns
EntityHandle&

Definition at line 340 of file Simple.hpp.

340 { return boundaryMeshset; }

◆ getDim()

int MoFEM::Simple::getDim ( ) const
inline

Get the problem dimension.

Problem dimension is determined by highest dimension of entities on the mesh.

Returns
int
Examples
child_and_parent.cpp, and field_evaluator.cpp.

Definition at line 312 of file Simple.hpp.

312 { return dIm; }

◆ getDM() [1/2]

SmartPetscObj<DM> MoFEM::Simple::getDM ( )
inline

Return smart DM object.

{
auto dm = simple_interface->getDM();
// ...
// dm is automatically destroyed when out of the scope
}
Returns
SmartPetscObj<DM>

Definition at line 302 of file Simple.hpp.

302 { return dM; }

◆ getDM() [2/2]

MoFEMErrorCode MoFEM::Simple::getDM ( DM *  dm)

◆ 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.

396 { return domainFE; }

◆ getDomainFEName() [2/2]

const std::string MoFEM::Simple::getDomainFEName ( ) const
inline

◆ getMeshSet()

DEPRECATED EntityHandle& MoFEM::Simple::getMeshSet ( )
inline
Deprecated:
Use getMeshset
Returns
EntityHandle&

Definition at line 326 of file Simple.hpp.

326 { return meshSet; }

◆ getMeshset()

EntityHandle& MoFEM::Simple::getMeshset ( )
inline

Get the MeshSet object.

Returns
EntityHandle&

Definition at line 333 of file Simple.hpp.

333 { return meshSet; }

◆ getOptions()

MoFEMErrorCode MoFEM::Simple::getOptions ( )

◆ 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.

427 { return otherFEs; }

◆ 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.

473 { return addParentAdjacencies; }

◆ getProblemName() [1/2]

std::string& MoFEM::Simple::getProblemName ( )
inline

Get the Problem Name.

Returns
std::string&

Definition at line 417 of file Simple.hpp.

417 { return nameOfProblem; }

◆ getProblemName() [2/2]

const std::string MoFEM::Simple::getProblemName ( ) const
inline

Get the Problem Name.

Returns
const std::string
Examples
child_and_parent.cpp, poisson_2d_homogeneous.cpp, reaction_diffusion.cpp, and test_cache_on_entities.cpp.

Definition at line 389 of file Simple.hpp.

389 { return nameOfProblem; }

◆ 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.

410 { return skeletonFE; }

◆ getSkeletonFEName() [2/2]

const std::string MoFEM::Simple::getSkeletonFEName ( ) const
inline

◆ getSkeletonMeshSet()

EntityHandle& MoFEM::Simple::getSkeletonMeshSet ( )
inline

Get the SkeletonMeshSet object.

Returns
EntityHandle&

Definition at line 347 of file Simple.hpp.

347 { return skeletonMeshset; }

◆ 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_namefile name if not set default or set by command line is used.
Returns
error code

Definition at line 250 of file Simple.cpp.

250  {
252  Interface &m_field = cOre;
253  if (m_field.get_comm_size() == 1)
255  else
256  CHKERR loadFile("PARALLEL=READ_PART;"
257  "PARALLEL_RESOLVE_SHARED_ENTS;"
258  "PARTITION=PARALLEL_PARTITION;",
259  meshFileName);
261 }

◆ loadFile() [2/2]

MoFEMErrorCode MoFEM::Simple::loadFile ( const std::string  options,
const std::string  mesh_file_name,
LoadFileFunc  loadFunc = defaultLoadFileFunc 
)

Load mesh file.

Parameters
optionsfile load options
mesh_file_namefile name if not set default or set by command line
loadFuncfunction 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.

196  {
197  Interface &m_field = cOre;
199  PetscLogEventBegin(MOFEM_EVENT_SimpleLoadMesh, 0, 0, 0, 0);
200 
201  if (!mesh_file_name.empty())
202  strcpy(meshFileName, mesh_file_name.c_str());
203 
204  // Invoke the boost function, passing in required arguments
205  CHKERR loadFunc(m_field, meshFileName, options.c_str());
206  CHKERR m_field.rebuild_database();
207 
208  // determine problem dimension
209  if (dIm == -1) {
210  int nb_ents_3d;
211  CHKERR m_field.get_moab().get_number_entities_by_dimension(
212  meshSet, 3, nb_ents_3d, true);
213  if (nb_ents_3d > 0) {
214  dIm = 3;
215  } else {
216  int nb_ents_2d;
217  CHKERR m_field.get_moab().get_number_entities_by_dimension(
218  meshSet, 2, nb_ents_2d, true);
219  if (nb_ents_2d > 0) {
220  dIm = 2;
221  } else {
222  dIm = 1;
223  }
224  }
225  }
226 
227  if (!boundaryMeshset)
229  if (!skeletonMeshset)
231  if (addSkeletonFE)
233 
234  if (bitLevel.any()) {
235  Range ents;
236  CHKERR m_field.get_moab().get_entities_by_dimension(meshSet, dIm, ents,
237  true);
238  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(ents, bitLevel,
239  false);
240  } else {
241  MOFEM_LOG("WORLD", Sev::warning) << "BitRefLevel is none and not set";
242  CHKERR m_field.getInterface<BitRefManager>()->addToDatabaseBitRefLevelByDim(
243  dIm, BitRefLevel().set(), BitRefLevel().set());
244  }
245 
246  PetscLogEventEnd(MOFEM_EVENT_SimpleLoadMesh, 0, 0, 0, 0);
248 }

◆ query_interface()

MoFEMErrorCode MoFEM::Simple::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Implements MoFEM::UnknownInterface.

Definition at line 8 of file Simple.cpp.

9  {
10  *iface = const_cast<Simple *>(this);
11  return 0;
12 }

◆ removeBoundaryField()

MoFEMErrorCode MoFEM::Simple::removeBoundaryField ( const std::string &  name)

Remove field form boundary.

Parameters
name
Returns
MoFEMErrorCode

Definition at line 370 of file Simple.cpp.

370  {
372 
373  auto remove_field_from_list = [&](auto &vec) {
374  auto it = std::find(vec.begin(), vec.end(), name);
375  if (it != vec.end())
376  vec.erase(it);
377  };
378 
379  remove_field_from_list(boundaryFields);
380 
382 }

◆ removeDomainField()

MoFEMErrorCode MoFEM::Simple::removeDomainField ( const std::string &  name)

Remove field form domain.

Parameters
name
Returns
MoFEMErrorCode

Definition at line 355 of file Simple.cpp.

355  {
357 
358  auto remove_field_from_list = [&](auto &vec) {
359  auto it = std::find(vec.begin(), vec.end(), name);
360  if (it != vec.end())
361  vec.erase(it);
362  };
363 
364  remove_field_from_list(noFieldFields);
365  remove_field_from_list(domainFields);
366 
368 }

◆ removeSkeletonField()

MoFEMErrorCode MoFEM::Simple::removeSkeletonField ( const std::string &  name)

Remove field form skeleton.

Parameters
name
Returns
MoFEMErrorCode

Definition at line 384 of file Simple.cpp.

384  {
386 
387  auto remove_field_from_list = [&](auto &vec) {
388  auto it = std::find(vec.begin(), vec.end(), name);
389  if (it != vec.end())
390  vec.erase(it);
391  };
392 
393  remove_field_from_list(skeletonFields);
394 
396 }

◆ reSetUp()

MoFEMErrorCode MoFEM::Simple::reSetUp ( bool  only_dm = false)

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 655 of file Simple.cpp.

655  {
656  Interface &m_field = cOre;
658  PetscLogEventBegin(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
659 
660  if (!only_dm) {
663  if (addSkeletonFE || !skeletonFields.empty()) {
665  if (addBoundaryFE || !boundaryFields.empty()) {
667  }
668  }
672  }
673 
674  CHKERR m_field.build_adjacencies(bitLevel, bitLevelMask);
675 
676  const Problem *problem_ptr;
677  CHKERR DMMoFEMGetProblemPtr(dM, &problem_ptr);
678  const auto problem_name = problem_ptr->getName();
679  CHKERR m_field.modify_problem_ref_level_set_bit(problem_name, bitLevel);
680  CHKERR m_field.modify_problem_mask_ref_level_set_bit(problem_name,
681  bitLevelMask);
682 
683  // Set problem by the DOFs on the fields rather that by adding DOFs on the
684  // elements
685  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_TRUE;
687  m_field.getInterface<ProblemsManager>()->buildProblemFromFields = PETSC_FALSE;
688 
689  PetscLogEventEnd(MOFEM_EVENT_SimpleBuildProblem, 0, 0, 0, 0);
691 }

◆ setDim()

void MoFEM::Simple::setDim ( int  dim)
inline

Set the problem dimension.

Returns
int

Definition at line 319 of file Simple.hpp.

319 { dIm = dim; };

◆ setFieldOrder()

MoFEMErrorCode MoFEM::Simple::setFieldOrder ( const std::string  field_name,
const int  order,
const Range ents = NULL 
)

◆ setParentAdjacency() [1/4]

template<>
MoFEMErrorCode MoFEM::Simple::setParentAdjacency ( )
private

Definition at line 72 of file Simple.cpp.

72  {
73  Interface &m_field = cOre;
75 
77  boost::make_shared<ParentFiniteElementAdjacencyFunction<3>>(
80  boost::make_shared<ParentFiniteElementAdjacencyFunction<2>>(
82 
83  CHKERR m_field.modify_finite_element_adjacency_table(domainFE, MBTET,
85  CHKERR m_field.modify_finite_element_adjacency_table(domainFE, MBHEX,
87  if (addBoundaryFE || !boundaryFields.empty()) {
88  CHKERR m_field.modify_finite_element_adjacency_table(
90  CHKERR m_field.modify_finite_element_adjacency_table(
92  }
93  if (addSkeletonFE || !skeletonFields.empty()) {
94  CHKERR m_field.modify_finite_element_adjacency_table(
96  CHKERR m_field.modify_finite_element_adjacency_table(
98  }
99 
101 }

◆ setParentAdjacency() [2/4]

template<>
MoFEMErrorCode MoFEM::Simple::setParentAdjacency ( )
private

Definition at line 103 of file Simple.cpp.

103  {
104  Interface &m_field = cOre;
106 
108  boost::make_shared<ParentFiniteElementAdjacencyFunction<2>>(
111  boost::make_shared<ParentFiniteElementAdjacencyFunction<1>>(
113 
114  CHKERR m_field.modify_finite_element_adjacency_table(domainFE, MBTRI,
116  CHKERR m_field.modify_finite_element_adjacency_table(domainFE, MBQUAD,
118  if (addBoundaryFE || !boundaryFields.empty())
119  CHKERR m_field.modify_finite_element_adjacency_table(
121  if (addSkeletonFE || !skeletonFields.empty())
122  CHKERR m_field.modify_finite_element_adjacency_table(
124 
126 }

◆ setParentAdjacency() [3/4]

template<>
MoFEMErrorCode MoFEM::Simple::setParentAdjacency ( )
private

Definition at line 128 of file Simple.cpp.

128  {
130  switch (getDim()) {
131  case 1:
132  THROW_MESSAGE("Not implemented");
133  case 2:
134  return setParentAdjacency<2>();
135  case 3:
136  return setParentAdjacency<3>();
137  default:
138  THROW_MESSAGE("Not implemented");
139  }
141 }

◆ setParentAdjacency() [4/4]

template<int DIM>
MoFEMErrorCode MoFEM::Simple::setParentAdjacency
private

Definition at line 67 of file Simple.cpp.

67  {
68  static_assert(DIM == 2 || DIM == 3, "not implemented");
69  return MOFEM_NOT_IMPLEMENTED;
70 }

◆ setSkeletonAdjacency() [1/5]

MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency ( int  dim = -1,
std::string  fe_name = "" 
)

Set the skeleton adjacency object.

Parameters
dim
fe_name
Returns
MoFEMErrorCode

Definition at line 143 of file Simple.cpp.

143  {
145  if (dim == -1)
146  dim = getDim();
147 
148  if (fe_name.empty())
149  fe_name = skeletonFE;
150 
151  switch (dim) {
152  case 2:
153  return setSkeletonAdjacency<2>(fe_name);
154  case 3:
155  return setSkeletonAdjacency<3>(fe_name);
156  default:
157  SETERRQ(PETSC_COMM_WORLD, MOFEM_NOT_IMPLEMENTED, "Not implemented");
158  }
160 }

◆ setSkeletonAdjacency() [2/5]

template<>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency ( std::string  fe_name)
private

Definition at line 21 of file Simple.cpp.

21  {
22  Interface &m_field = cOre;
24 
26  boost::make_shared<ParentFiniteElementAdjacencyFunctionSkeleton<1>>(
28  CHKERR m_field.modify_finite_element_adjacency_table(
29  fe_name, MBEDGE, *parentAdjSkeletonFunctionDim1);
30 
32 }

◆ setSkeletonAdjacency() [3/5]

template<>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency ( std::string  fe_name)
private

Definition at line 35 of file Simple.cpp.

35  {
36  Interface &m_field = cOre;
38 
40  boost::make_shared<ParentFiniteElementAdjacencyFunctionSkeleton<2>>(
42 
43  CHKERR m_field.modify_finite_element_adjacency_table(
44  fe_name, MBTRI, *parentAdjSkeletonFunctionDim2);
45  CHKERR m_field.modify_finite_element_adjacency_table(
46  fe_name, MBQUAD, *parentAdjSkeletonFunctionDim2);
47 
49 }

◆ setSkeletonAdjacency() [4/5]

template<>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency ( std::string  fe_name)
private

Definition at line 52 of file Simple.cpp.

52  {
54 
55  switch (getDim()) {
56  case 1:
57  THROW_MESSAGE("Not implemented");
58  case 2:
59  return setSkeletonAdjacency<2>(fe_name);
60  case 3:
61  return setSkeletonAdjacency<3>(fe_name);
62  default:
63  THROW_MESSAGE("Not implemented");
64  }
66 }

◆ setSkeletonAdjacency() [5/5]

template<int DIM>
MoFEMErrorCode MoFEM::Simple::setSkeletonAdjacency ( std::string  fe_name)
private

Definition at line 15 of file Simple.cpp.

15  {
16  static_assert(DIM == 2 || DIM == 3, "not implemented");
17  return MOFEM_NOT_IMPLEMENTED;
18 }

◆ setUp()

MoFEMErrorCode MoFEM::Simple::setUp ( const PetscBool  is_partitioned = PETSC_TRUE)

Member Data Documentation

◆ addBoundaryFE

bool MoFEM::Simple::addBoundaryFE
private

Add boundary FE.

Definition at line 535 of file Simple.hpp.

◆ 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

Add skeleton FE.

Definition at line 534 of file Simple.hpp.

◆ bitAdjEnt

BitRefLevel MoFEM::Simple::bitAdjEnt
private

bit ref level for parent

Definition at line 540 of file Simple.hpp.

◆ bitAdjEntMask

BitRefLevel MoFEM::Simple::bitAdjEntMask
private

bit ref level for parent parent

Definition at line 541 of file Simple.hpp.

◆ bitAdjParent

BitRefLevel MoFEM::Simple::bitAdjParent
private

bit ref level for parent

Definition at line 538 of file Simple.hpp.

◆ bitAdjParentMask

BitRefLevel MoFEM::Simple::bitAdjParentMask
private

bit ref level for parent parent

Definition at line 539 of file Simple.hpp.

◆ bitLevel

BitRefLevel MoFEM::Simple::bitLevel
private

BitRefLevel of the problem.

Definition at line 520 of file Simple.hpp.

◆ bitLevelMask

BitRefLevel MoFEM::Simple::bitLevelMask
private

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

boundary fields

Definition at line 544 of file Simple.hpp.

◆ boundaryMeshset

EntityHandle MoFEM::Simple::boundaryMeshset
private

meshset with boundary

Definition at line 531 of file Simple.hpp.

◆ cOre

MoFEM::Core& MoFEM::Simple::cOre
private

Definition at line 518 of file Simple.hpp.

◆ dataFields

std::vector<std::string> MoFEM::Simple::dataFields
private

Data fields.

Definition at line 546 of file Simple.hpp.

◆ dIm

int MoFEM::Simple::dIm
private

dimension of problem

Definition at line 562 of file Simple.hpp.

◆ dM

SmartPetscObj<DM> MoFEM::Simple::dM
private

Discrete manager (interface to PETSc/MoFEM functions)

Definition at line 565 of file Simple.hpp.

◆ domainFE

std::string MoFEM::Simple::domainFE
private

domain finite element

Definition at line 555 of file Simple.hpp.

◆ domainFields

std::vector<std::string> MoFEM::Simple::domainFields
private

domain fields

Definition at line 543 of file Simple.hpp.

◆ 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

mesh file name

Definition at line 561 of file Simple.hpp.

◆ meshSet

EntityHandle MoFEM::Simple::meshSet
private

domain meshset

Definition at line 530 of file Simple.hpp.

◆ MOFEM_EVENT_SimpleBuildFields

PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleBuildFields
private

Definition at line 525 of file Simple.hpp.

◆ MOFEM_EVENT_SimpleBuildFiniteElements

PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleBuildFiniteElements
private

Definition at line 526 of file Simple.hpp.

◆ MOFEM_EVENT_SimpleBuildProblem

PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleBuildProblem
private

Definition at line 527 of file Simple.hpp.

◆ MOFEM_EVENT_SimpleKSPSolve

PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleKSPSolve
private

Definition at line 528 of file Simple.hpp.

◆ MOFEM_EVENT_SimpleLoadMesh

PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleLoadMesh
private

Definition at line 524 of file Simple.hpp.

◆ MOFEM_EVENT_SimpleSetUP

PetscLogEvent MoFEM::Simple::MOFEM_EVENT_SimpleSetUP
private

Definition at line 523 of file Simple.hpp.

◆ nameOfProblem

std::string MoFEM::Simple::nameOfProblem
private

problem name

Definition at line 554 of file Simple.hpp.

◆ noFieldDataFields

std::vector<std::string> MoFEM::Simple::noFieldDataFields
private

NOFIELD field name.

Definition at line 548 of file Simple.hpp.

◆ noFieldFields

std::vector<std::string> MoFEM::Simple::noFieldFields
private

NOFIELD field name.

Definition at line 547 of file Simple.hpp.

◆ otherFEs

std::vector<std::string> MoFEM::Simple::otherFEs
private

Other finite elements.

Definition at line 559 of file Simple.hpp.

◆ parentAdjFunctionDim1

boost::shared_ptr<ParentFiniteElementAdjacencyFunction<1> > MoFEM::Simple::parentAdjFunctionDim1
private

Definition at line 581 of file Simple.hpp.

◆ parentAdjFunctionDim2

boost::shared_ptr<ParentFiniteElementAdjacencyFunction<2> > MoFEM::Simple::parentAdjFunctionDim2
private

Definition at line 579 of file Simple.hpp.

◆ parentAdjFunctionDim3

boost::shared_ptr<ParentFiniteElementAdjacencyFunction<3> > MoFEM::Simple::parentAdjFunctionDim3
private

Definition at line 577 of file Simple.hpp.

◆ parentAdjSkeletonFunctionDim1

boost::shared_ptr<ParentFiniteElementAdjacencyFunctionSkeleton<1> > MoFEM::Simple::parentAdjSkeletonFunctionDim1
private

Definition at line 586 of file Simple.hpp.

◆ parentAdjSkeletonFunctionDim2

boost::shared_ptr<ParentFiniteElementAdjacencyFunctionSkeleton<2> > MoFEM::Simple::parentAdjSkeletonFunctionDim2
private

Definition at line 584 of file Simple.hpp.

◆ 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

EntityHandle MoFEM::Simple::skeletonMeshset
private

skeleton meshset with boundary

Definition at line 532 of file Simple.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::Simple::boundaryFields
std::vector< std::string > boundaryFields
boundary fields
Definition: Simple.hpp:544
MoFEM::Simple::parentAdjFunctionDim1
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 1 > > parentAdjFunctionDim1
Definition: Simple.hpp:581
MoFEM::Simple::setSkeletonAdjacency
MoFEMErrorCode setSkeletonAdjacency(int dim=-1, std::string fe_name="")
Set the skeleton adjacency object.
Definition: Simple.cpp:143
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::Simple::otherFEs
std::vector< std::string > otherFEs
Other finite elements.
Definition: Simple.hpp:559
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::Simple::addParentAdjacencies
bool addParentAdjacencies
If set to true parent adjacencies are build.
Definition: Simple.hpp:536
MoFEM::Simple::addBoundaryFE
bool addBoundaryFE
Add boundary FE.
Definition: Simple.hpp:535
MOFEM_LOG_SEVERITY_SYNC
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
Definition: LogManager.hpp:352
MoFEM::Simple::exchangeGhostCells
MoFEMErrorCode exchangeGhostCells()
Definition: Simple.cpp:818
MoFEM::Simple::fieldsOrder
std::list< std::tuple< std::string, int, Range, bool > > fieldsOrder
Definition: Simple.hpp:551
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
MoFEM::Simple::Simple
Simple(const MoFEM::Core &core)
Definition: Simple.cpp:162
MoFEM::Simple::parentAdjFunctionDim2
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 2 > > parentAdjFunctionDim2
Definition: Simple.hpp:579
MoFEM::Simple::noFieldFields
std::vector< std::string > noFieldFields
NOFIELD field name.
Definition: Simple.hpp:547
MoFEM::Simple::meshSet
EntityHandle meshSet
domain meshset
Definition: Simple.hpp:530
MoFEM::Simple::createBoundaryMeshset
MoFEMErrorCode createBoundaryMeshset()
Definition: Simple.cpp:737
MoFEM::Simple::meshFileName
char meshFileName[255]
mesh file name
Definition: Simple.hpp:561
MoFEM::Simple::MOFEM_EVENT_SimpleKSPSolve
PetscLogEvent MOFEM_EVENT_SimpleKSPSolve
Definition: Simple.hpp:528
MoFEM::Simple::buildProblem
MoFEMErrorCode buildProblem()
Build problem.
Definition: Simple.cpp:615
MoFEM::Simple::MOFEM_EVENT_SimpleBuildFields
PetscLogEvent MOFEM_EVENT_SimpleBuildFields
Definition: Simple.hpp:525
MoFEM::getPetscObject
PetscObject getPetscObject(T obj)
Definition: PetscSmartObj.hpp:9
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
MoFEM::DMSetUp_MoFEM
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition: DMMoFEM.cpp:1303
MoFEM::Simple::dIm
int dIm
dimension of problem
Definition: Simple.hpp:562
MoFEM::Simple::addSkeletonFE
bool addSkeletonFE
Add skeleton FE.
Definition: Simple.hpp:534
MoFEM::Simple::nameOfProblem
std::string nameOfProblem
problem name
Definition: Simple.hpp:554
MoFEM::Simple::parentAdjSkeletonFunctionDim2
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 2 > > parentAdjSkeletonFunctionDim2
Definition: Simple.hpp:584
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:574
MoFEM::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:497
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2002
MoFEM::Simple::MOFEM_EVENT_SimpleLoadMesh
PetscLogEvent MOFEM_EVENT_SimpleLoadMesh
Definition: Simple.hpp:524
MoFEM::Simple::bitLevelMask
BitRefLevel bitLevelMask
BitRefLevel of the problem.
Definition: Simple.hpp:521
MoFEM::Simple::dataFields
std::vector< std::string > dataFields
Data fields.
Definition: Simple.hpp:546
VERBOSE
@ VERBOSE
Definition: definitions.h:222
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::Simple::bitAdjEntMask
BitRefLevel bitAdjEntMask
bit ref level for parent parent
Definition: Simple.hpp:541
MoFEM::createDM
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
Definition: PetscSmartObj.hpp:141
MoFEM::Simple::domainFE
std::string domainFE
domain finite element
Definition: Simple.hpp:555
MoFEM::Simple::MOFEM_EVENT_SimpleBuildFiniteElements
PetscLogEvent MOFEM_EVENT_SimpleBuildFiniteElements
Definition: Simple.hpp:526
MoFEM::Simple::createSkeletonMeshset
MoFEMErrorCode createSkeletonMeshset()
Definition: Simple.cpp:784
MoFEM::Simple::boundaryFE
std::string boundaryFE
boundary finite element
Definition: Simple.hpp:556
MoFEM::Simple::parentAdjFunctionDim3
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 3 > > parentAdjFunctionDim3
Definition: Simple.hpp:577
MoFEM::Simple::bitAdjParentMask
BitRefLevel bitAdjParentMask
bit ref level for parent parent
Definition: Simple.hpp:539
MoFEM::Simple::skeletonFE
std::string skeletonFE
skeleton finite element
Definition: Simple.hpp:557
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
MoFEM::Simple::bitAdjEnt
BitRefLevel bitAdjEnt
bit ref level for parent
Definition: Simple.hpp:540
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
MoFEM::DMMoFEMCreateMoFEM
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.
Definition: DMMoFEM.cpp:114
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::Simple::buildFiniteElements
MoFEMErrorCode buildFiniteElements()
Build finite elements.
Definition: Simple.cpp:589
Range
FTensor::dd
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)
Definition: ddTensor0.hpp:33
MF_ZERO
@ MF_ZERO
Definition: definitions.h:111
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::Simple::domainFields
std::vector< std::string > domainFields
domain fields
Definition: Simple.hpp:543
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1140
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::Simple::skeletonFields
std::vector< std::string > skeletonFields
fields on the skeleton
Definition: Simple.hpp:545
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::DMMoFEMGetProblemPtr
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:426
MoFEM::Simple::bitLevel
BitRefLevel bitLevel
BitRefLevel of the problem.
Definition: Simple.hpp:520
MoFEM::Simple::defineFiniteElements
MoFEMErrorCode defineFiniteElements()
Define finite elements.
Definition: Simple.cpp:398
MoFEM::Simple::getDim
int getDim() const
Get the problem dimension.
Definition: Simple.hpp:312
MoFEM::Simple::buildFields
MoFEMErrorCode buildFields()
Build fields.
Definition: Simple.cpp:500
MoFEM::Simple::bitAdjParent
BitRefLevel bitAdjParent
bit ref level for parent
Definition: Simple.hpp:538
MoFEM::Simple::noFieldDataFields
std::vector< std::string > noFieldDataFields
NOFIELD field name.
Definition: Simple.hpp:548
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
FiniteElement_multiIndex
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.
Definition: FEMultiIndices.hpp:849
MoFEM::Simple::boundaryMeshset
EntityHandle boundaryMeshset
meshset with boundary
Definition: Simple.hpp:531
MoFEM::Simple::parentAdjSkeletonFunctionDim1
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 1 > > parentAdjSkeletonFunctionDim1
Definition: Simple.hpp:586
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::Simple::defineProblem
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
Definition: Simple.cpp:470
MoFEM::Simple::dM
SmartPetscObj< DM > dM
Discrete manager (interface to PETSc/MoFEM functions)
Definition: Simple.hpp:565
MoFEM::Simple::setParentAdjacency
MoFEMErrorCode setParentAdjacency()
Definition: Simple.cpp:67
MoFEM::Simple::getProblemName
const std::string getProblemName() const
Get the Problem Name.
Definition: Simple.hpp:389
MoFEM::Simple::addFieldToEmptyFieldBlocks
MoFEMErrorCode addFieldToEmptyFieldBlocks(const std::string row_field, const std::string col_field) const
add empty block to problem
Definition: Simple.cpp:728
MoFEM::Simple::cOre
MoFEM::Core & cOre
Definition: Simple.hpp:518
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::Simple::skeletonMeshset
EntityHandle skeletonMeshset
skeleton meshset with boundary
Definition: Simple.hpp:532
MoFEM::Simple::MOFEM_EVENT_SimpleSetUP
PetscLogEvent MOFEM_EVENT_SimpleSetUP
Definition: Simple.hpp:523
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
NOFIELD
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition: definitions.h:84
MoFEM::DMMoFEMSetIsPartitioned
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMoFEM.cpp:1123
MoFEM::Simple::MOFEM_EVENT_SimpleBuildProblem
PetscLogEvent MOFEM_EVENT_SimpleBuildProblem
Definition: Simple.hpp:527