|
| v0.14.0
|
Go to the source code of this file.
|
int | main (int argc, char *argv[]) |
|
◆ main()
int main |
( |
int |
argc, |
|
|
char * |
argv[] |
|
) |
| |
- Examples
- remove_entities_from_problem_not_partitioned.cpp.
Definition at line 14 of file remove_entities_from_problem_not_partitioned.cpp.
23 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
25 PetscBool flg = PETSC_TRUE;
27 #if PETSC_VERSION_GE(3, 6, 4)
35 SETERRQ(PETSC_COMM_WORLD, 1,
"*** ERROR -my_file (MESH FILE NEEDED)");
47 BOOST_LOG_SCOPED_THREAD_ATTR(
"Timeline", attrs::timer());
69 for (EntityType
t = CN::TypeDimensionMap[1].first;
73 for (EntityType
t = CN::TypeDimensionMap[2].first;
108 auto get_triangles_on_skin = [&](
Range &tets_skin) {
114 CHKERR skin.find_skin(0, tets,
false, tets_skin);
117 adj, moab::Interface::UNION);
118 tets_skin.merge(adj);
127 "P1",
ROW,
"F1", 0, 1, is_before_remove, &tets);
129 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Remove dofs";
132 CHKERR get_triangles_on_skin(tets_skin);
134 "P1",
"F1", tets_skin, 0, 1, 0, 100,
VERBOSE,
true);
136 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Check consistency";
140 "P1",
ROW,
"F1", 0, 1, is_after_remove, &tets);
142 auto test_is = [&]() {
146 auto sub_ao = sub_data->getSmartRowMap();
148 CHKERR AOApplicationToPetscIS(sub_ao, is_before_remove);
149 PetscBool is_the_same;
150 CHKERR ISEqual(is_before_remove, is_after_remove, &is_the_same);
151 if (is_the_same == PETSC_FALSE) {
153 "IS should be the same if map is correctly implemented");
155 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Sub data map is correct";
159 "AO map should exist");
163 "Sub DM should exist");
174 ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>(
"P1", -2,
178 MOFEM_LOG(
"WORLD", Sev::inform) <<
"done";
◆ help
◆ SPACE_DIM
constexpr int SPACE_DIM = 3 |
|
constexpr |
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Problem manager is used to build and partition problems.
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
virtual MoFEMErrorCode add_ents_to_field_by_dim(const Range &ents, const int dim, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
Deprecated interface functions.
DeprecatedCoreInterface Interface
#define CHKERR
Inline error check.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual moab::Interface & get_moab()=0
MoFEMErrorCode removeDofsOnEntitiesNotDistributed(const std::string problem_name, const std::string field_name, const Range ents, const int lo_coeff=0, const int hi_coeff=MAX_DOFS_ON_ENTITY, const int lo_order=0, const int hi_order=100, int verb=VERBOSE, const bool debug=false)
Remove DOFs from problem.
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
Section manager is used to create indexes and sections.
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
boost::shared_ptr< SubProblemData > & getSubData() const
Get main problem of sub-problem is.
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to finite element
constexpr double t
plate stiffness
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
Vector manager is used to create vectors \mofem_vectors.
Matrix manager is used to build and partition problems.
const double v
phase velocity of light in medium (cm/ns)
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define MOFEM_LOG(channel, severity)
Log.
#define CATCH_ERRORS
Catch errors.
@ HCURL
field with continuous tangents
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
@ MOFEM_ATOM_TEST_INVALID
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
keeps basic data about problem
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
virtual MoFEMErrorCode add_field(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_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)