v0.14.0
Loading...
Searching...
No Matches
remove_entities_from_problem.cpp

Remove field entities from the problem.

Remove field entities from the problem

/** \file remove_entities_from_problem.cpp
\example remove_entities_from_problem.cpp
\brief Remove field entities from the problem
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
constexpr int SPACE_DIM = 3;
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
#if PETSC_VERSION_GE(3, 6, 4)
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
#else
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
mesh_file_name, 255, &flg);
#endif
if (!flg)
SETERRQ(PETSC_COMM_WORLD, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
if (pcomm == NULL)
pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
const char *option;
option = "PARALLEL=BCAST_DELETE;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION="
"PARALLEL_PARTITION;";
CHKERR moab.load_file(mesh_file_name, 0, option);
// Create MoFEM database
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
// set entities bit level
const auto bit_level = BitRefLevel().set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_level);
// Fields
// meshset consisting all entities in mesh
EntityHandle root_set = moab.get_root_set();
// add entities to field
CHKERR m_field.add_ents_to_field_by_dim(root_set, SPACE_DIM, "F1");
CHKERR m_field.add_ents_to_field_by_dim(root_set, SPACE_DIM, "F2");
// set app. order
// see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
// (Mark Ainsworth & Joe Coyle)
int order = 2;
for (EntityType t = CN::TypeDimensionMap[1].first;
t <= CN::TypeDimensionMap[SPACE_DIM].second; ++t) {
CHKERR m_field.set_field_order(root_set, t, "F1", order);
}
for (EntityType t = CN::TypeDimensionMap[2].first;
t <= CN::TypeDimensionMap[SPACE_DIM].second; ++t) {
CHKERR m_field.set_field_order(root_set, t, "F2", order);
}
CHKERR m_field.build_fields();
// add elements
CHKERR m_field.add_finite_element("E1");
CHKERR m_field.build_adjacencies(bit_level);
// Problems
CHKERR m_field.add_problem("P1");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_level);
// build problems
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblemOnDistributedMesh("P1", true);
CHKERR prb_mng_ptr->partitionFiniteElements("P1");
CHKERR prb_mng_ptr->partitionGhostDofs("P1");
auto get_triangles_on_skin = [&](Range &tets_skin) {
Range tets;
CHKERR m_field.get_moab().get_entities_by_dimension(root_set, SPACE_DIM,
tets);
Skinner skin(&m_field.get_moab());
CHKERR skin.find_skin(0, tets, false, tets_skin);
Range adj;
CHKERR m_field.get_moab().get_adjacencies(tets_skin, SPACE_DIM - 2, false,
adj, moab::Interface::UNION);
tets_skin.merge(adj);
};
Range tets;
CHKERR m_field.get_moab().get_entities_by_dimension(root_set, SPACE_DIM,
tets);
SmartPetscObj<IS> is_before_remove;
CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
"P1", ROW, "F1", 0, 1, is_before_remove, &tets);
Range tets_skin;
CHKERR get_triangles_on_skin(tets_skin);
CHKERR prb_mng_ptr->removeDofsOnEntities("P1", "F1", tets_skin, 0, 1, 0,
100, NOISY, true);
SmartPetscObj<IS> is_after_remove;
CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
"P1", ROW, "F1", 0, 1, is_after_remove, &tets);
auto test_is = [&](auto prb_name, auto is_before, auto is_after) {
const Problem *prb_ptr = m_field.get_problem(prb_name);
if (auto sub_data = prb_ptr->getSubData()) {
auto sub_ao = sub_data->getSmartRowMap();
if (sub_ao) {
CHKERR AOApplicationToPetscIS(sub_ao, is_before);
PetscBool is_the_same;
CHKERR ISEqual(is_before, is_after, &is_the_same);
if (is_the_same == PETSC_FALSE) {
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"IS should be the same if map is correctly implemented");
} else {
MOFEM_LOG("WORLD", Sev::inform) << "Sub data map is correct";
}
} else {
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"AO map should exist");
}
} else {
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Sub DM should exist");
}
};
CHKERR test_is("P1", is_before_remove, is_after_remove);
->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("P1", -2, -2,
0);
CHKERR m_field.add_problem("SUB");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("SUB", bit_level);
CHKERR prb_mng_ptr->buildSubProblem("SUB", {"F1"}, {"F1"}, "P1", PETSC_TRUE,
nullptr, nullptr);
CHKERR prb_mng_ptr->partitionFiniteElements("SUB", true, 0,
m_field.get_comm_size());
CHKERR prb_mng_ptr->removeDofsOnEntities("SUB", "F1", tets_skin, 0, 1, 0,
100, NOISY, true);
SmartPetscObj<IS> is_sub_after_remove;
CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
"SUB", ROW, "F1", 0, MAX_DOFS_ON_ENTITY, is_sub_after_remove, &tets);
CHKERR test_is("SUB", is_before_remove, is_after_remove);
->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("SUB", -2,
-2, 0);
}
// finish work cleaning memory, getting statistics, etc.
return 0;
}
static char help[]
int main()
constexpr int SPACE_DIM
@ NOISY
@ ROW
#define CATCH_ERRORS
Catch errors.
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MYPCOMM_INDEX
default communicator number PCOMM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
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
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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 MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
#define MOFEM_LOG(channel, severity)
Log.
MoFEMErrorCode buildProblemOnDistributedMesh(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures, assuming that mesh is distributed (collective)
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode buildSubProblem(const std::string out_name, const std::vector< std::string > &fields_row, const std::vector< std::string > &fields_col, const std::string main_problem, const bool square_matrix=true, const map< std::string, boost::shared_ptr< Range > > *entityMapRow=nullptr, const map< std::string, boost::shared_ptr< Range > > *entityMapCol=nullptr, int verb=VERBOSE)
build sub problem
MoFEMErrorCode partitionGhostDofsOnDistributedMesh(const std::string name, int verb=VERBOSE)
determine ghost nodes on distributed meshes
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
MoFEMErrorCode removeDofsOnEntities(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 add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
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
char mesh_file_name[255]
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr double t
plate stiffness
Definition plate.cpp:59
Managing BitRefLevels.
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
Matrix manager is used to build and partition problems.
keeps basic data about problem
boost::shared_ptr< SubProblemData > & getSubData() const
Get main problem of sub-problem is.
Problem manager is used to build and partition problems.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.