|
| v0.14.0
|
Go to the source code of this file.
|
int | main (int argc, char *argv[]) |
|
|
static char | help [] = "...\n\n" |
|
◆ main()
int main |
( |
int |
argc, |
|
|
char * |
argv[] |
|
) |
| |
- Examples
- loop_entities.cpp.
Definition at line 47 of file loop_entities.cpp.
61 DMType dm_name =
"DMMOFEM";
83 CHKERR m_field.
get_moab().get_entities_by_handle(root_set, all_ents);
84 all_ents = subtract(all_ents, all_ents.subset_by_type(MBENTITYSET));
86 auto testingEntitiesInDatabase = [&]() {
93 CHKERR testingEntitiesInDatabase();
95 auto testingEntitiesInDatabaseSubRange = [&]() {
97 Range edges = all_ents.subset_by_type(MBEDGE);
103 CHKERR testingEntitiesInDatabaseSubRange();
105 auto dm = simple_interface->
getDM();
107 auto testingEntitiesInDM = [&]() {
117 CHKERR testingEntitiesInDM();
◆ help
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Simple interface for fast problem set-up.
Deprecated interface functions.
DeprecatedCoreInterface Interface
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
#define CHKERR
Inline error check.
virtual moab::Interface & get_moab()=0
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.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
virtual int get_comm_size() const =0
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
keeps basic data about problem
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
virtual MoFEMErrorCode loop_entities(const std::string field_name, EntityMethod &method, Range const *const ents=nullptr, int verb=DEFAULT_VERBOSITY)=0
Loop over field entities.