|
| v0.14.0
|
Go to the documentation of this file.
13 #ifndef __SIMPLE_HPP__
14 #define __SIMPLE_HPP__
20 template <
int DIM>
struct ParentFiniteElementAdjacencyFunction;
21 template <
int DIM>
struct ParentFiniteElementAdjacencyFunctionSkeleton;
52 const char *file_name,
53 const char *options) {
55 return m_field.
get_moab().load_file(file_name, 0, options);
101 const TagType tag_type = MB_TAG_SPARSE,
122 const TagType tag_type = MB_TAG_SPARSE,
143 const TagType tag_type = MB_TAG_SPARSE,
164 const TagType tag_type = MB_TAG_SPARSE,
184 const TagType tag_type = MB_TAG_SPARSE,
233 const Range *ents = NULL);
515 const std::string col_field)
const;
550 std::list<std::tuple<std::string, int, Range, bool>>
571 template <
int DIM = -1>
576 boost::shared_ptr<ParentFiniteElementAdjacencyFunction<3>>
578 boost::shared_ptr<ParentFiniteElementAdjacencyFunction<2>>
580 boost::shared_ptr<ParentFiniteElementAdjacencyFunction<1>>
583 boost::shared_ptr<ParentFiniteElementAdjacencyFunctionSkeleton<2>>
585 boost::shared_ptr<ParentFiniteElementAdjacencyFunctionSkeleton<1>>
591 #endif // __SIMPLE_HPP__
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.
std::vector< std::string > otherFEs
Other finite elements.
bool addParentAdjacencies
If set to true parent adjacencies are build.
bool addBoundaryFE
Add boundary FE.
MoFEMErrorCode exchangeGhostCells()
std::list< std::tuple< std::string, int, Range, bool > > fieldsOrder
Simple(const MoFEM::Core &core)
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 2 > > parentAdjFunctionDim2
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
std::vector< std::string > noFieldFields
NOFIELD field name.
EntityHandle meshSet
domain meshset
SmartPetscObj< DM > getDM()
Return smart DM object.
auto & getBitAdjParentMask()
bit ref level for parent parent
std::string & getDomainFEName()
Get the Domain FE Name.
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
const std::string getSkeletonFEName() const
Get the Skeleton FE Name.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
int dIm
dimension of problem
EntityHandle & getBoundaryMeshSet()
Get the BoundaryMeshSet object.
EntityHandle & getMeshset()
Get the MeshSet object.
BitRefLevel & getBitRefLevelMask()
Get the BitRefLevel.
bool addSkeletonFE
Add skeleton FE.
std::string nameOfProblem
problem name
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 2 > > parentAdjSkeletonFunctionDim2
Simple interface for fast problem set-up.
auto & getBitAdjEnt()
bit ref level for parent
Deprecated interface functions.
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.
MoFEMErrorCode removeDomainField(const std::string &name)
Remove field form domain.
PetscLogEvent MOFEM_EVENT_SimpleLoadMesh
MoFEMErrorCode deleteDM()
Delete dm.
MoFEMErrorCode getOptions()
get options
BitRefLevel bitLevelMask
BitRefLevel of the problem.
std::vector< std::string > dataFields
Data fields.
FieldSpace
approximation spaces
bool & getAddSkeletonFE()
Get the addSkeletonFE.
BitRefLevel bitAdjEntMask
bit ref level for parent parent
virtual moab::Interface & get_moab()=0
implementation of Data Operators for Forces and Sources
MoFEMErrorCode deleteFiniteElements()
Delete finite elements.
std::string & getProblemName()
Get the Problem Name.
std::string domainFE
domain finite element
DEPRECATED EntityHandle & getMeshSet()
int FieldCoefficientsNumber
Number of field coefficients.
PetscLogEvent MOFEM_EVENT_SimpleBuildFiniteElements
boost::function< MoFEMErrorCode(Interface &, const char *, const char *)> LoadFileFunc
MoFEMErrorCode createSkeletonMeshset()
std::string boundaryFE
boundary finite element
virtual ~Simple()=default
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.
bool & getAddBoundaryFE()
Get the addSkeletonFE.
const std::string getDomainFEName() const
Get the Domain FE Name.
boost::shared_ptr< ParentFiniteElementAdjacencyFunction< 3 > > parentAdjFunctionDim3
BitRefLevel bitAdjParentMask
bit ref level for parent parent
std::string skeletonFE
skeleton finite element
std::string & getBoundaryFEName()
Get the Boundary FE Name.
MoFEMErrorCode removeBoundaryField(const std::string &name)
Remove field form boundary.
EntityHandle & getSkeletonMeshSet()
Get the SkeletonMeshSet object.
BitRefLevel bitAdjEnt
bit ref level for parent
void setDim(int dim)
Set the problem dimension.
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.
static MoFEMErrorCode defaultLoadFileFunc(Interface &m_field, const char *file_name, const char *options)
Set function for loading mesh file.
constexpr auto field_name
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
MoFEMErrorCode buildFiniteElements()
Build finite elements.
auto & getBitAdjParent()
bit ref level for parent
base class for all interface classes
auto & getBitAdjEntMask()
bit ref level for parent parent
std::string & getSkeletonFEName()
Get the Skeleton FE Name.
std::vector< std::string > domainFields
domain fields
bool & getParentAdjacencies()
Get the addParentAdjacencies.
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.
std::vector< std::string > skeletonFields
fields on the skeleton
const std::string getBoundaryFEName() const
Get the Boundary FE Name.
BitRefLevel bitLevel
BitRefLevel of the problem.
MoFEMErrorCode defineFiniteElements()
Define finite elements.
BitRefLevel & getBitRefLevel()
Get the BitRefLevel.
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.
FieldApproximationBase
approximation base
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
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.
MoFEMErrorCode reSetUp(bool only_dm=false)
Rebuild internal MoFEM data structures.
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.
MoFEMErrorCode defineProblem(const PetscBool is_partitioned=PETSC_TRUE)
define problem
SmartPetscObj< DM > dM
Discrete manager (interface to PETSc/MoFEM functions)
MoFEMErrorCode setParentAdjacency()
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
std::vector< std::string > & getOtherFiniteElements()
Get the Other Finite Elements.
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
EntityHandle skeletonMeshset
skeleton meshset with boundary
PetscLogEvent MOFEM_EVENT_SimpleSetUP
PetscLogEvent MOFEM_EVENT_SimpleBuildProblem
MoFEMErrorCode removeSkeletonField(const std::string &name)
Remove field form skeleton.