v0.14.0
Functions | Variables
remove_entities_from_problem.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 
constexpr int SPACE_DIM = 3
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
remove_entities_from_problem.cpp.

Definition at line 14 of file remove_entities_from_problem.cpp.

14  {
15 
16  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
17 
18  try {
19 
20  moab::Core mb_instance;
21  moab::Interface &moab = mb_instance;
22  int rank;
23  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
24 
25  PetscBool flg = PETSC_TRUE;
26  char mesh_file_name[255];
27 #if PETSC_VERSION_GE(3, 6, 4)
28  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
29  255, &flg);
30 #else
31  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
32  mesh_file_name, 255, &flg);
33 #endif
34  if (!flg)
35  SETERRQ(PETSC_COMM_WORLD, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
36 
37  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
38  if (pcomm == NULL)
39  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
40  const char *option;
41  option = "PARALLEL=BCAST_DELETE;PARALLEL_RESOLVE_SHARED_ENTS;PARTITION="
42  "PARALLEL_PARTITION;";
43  CHKERR moab.load_file(mesh_file_name, 0, option);
44 
45  // Create MoFEM database
46  MoFEM::Core core(moab);
47  MoFEM::Interface &m_field = core;
48 
49  // set entities bit level
50  const auto bit_level = BitRefLevel().set(0);
51  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
52  0, 3, bit_level);
53 
54  // Fields
55  CHKERR m_field.add_field("F1", HCURL, DEMKOWICZ_JACOBI_BASE, 3);
56  CHKERR m_field.add_field("F2", HDIV, DEMKOWICZ_JACOBI_BASE, 1);
57 
58  // meshset consisting all entities in mesh
59  EntityHandle root_set = moab.get_root_set();
60  // add entities to field
61  CHKERR m_field.add_ents_to_field_by_dim(root_set, SPACE_DIM, "F1");
62  CHKERR m_field.add_ents_to_field_by_dim(root_set, SPACE_DIM, "F2");
63 
64  // set app. order
65  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
66  // (Mark Ainsworth & Joe Coyle)
67  int order = 2;
68  for (EntityType t = CN::TypeDimensionMap[1].first;
69  t <= CN::TypeDimensionMap[SPACE_DIM].second; ++t) {
70  CHKERR m_field.set_field_order(root_set, t, "F1", order);
71  }
72  for (EntityType t = CN::TypeDimensionMap[2].first;
73  t <= CN::TypeDimensionMap[SPACE_DIM].second; ++t) {
74  CHKERR m_field.set_field_order(root_set, t, "F2", order);
75  }
76 
77  CHKERR m_field.build_fields();
78 
79  // add elements
80  CHKERR m_field.add_finite_element("E1");
81 
82  CHKERR m_field.modify_finite_element_add_field_row("E1", "F1");
83  CHKERR m_field.modify_finite_element_add_field_row("E1", "F2");
84  CHKERR m_field.modify_finite_element_add_field_col("E1", "F2");
85  CHKERR m_field.modify_finite_element_add_field_col("E1", "F1");
86  CHKERR m_field.modify_finite_element_add_field_data("E1", "F1");
87  CHKERR m_field.modify_finite_element_add_field_data("E1", "F2");
88 
89  CHKERR m_field.add_ents_to_finite_element_by_dim(root_set, SPACE_DIM, "E1");
90 
91  CHKERR m_field.build_finite_elements();
92  CHKERR m_field.build_adjacencies(bit_level);
93 
94  // Problems
95  CHKERR m_field.add_problem("P1");
96 
97  // set refinement level for problem
98  CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_level);
99  CHKERR m_field.modify_problem_add_finite_element("P1", "E1");
100 
101  // build problems
102  ProblemsManager *prb_mng_ptr;
103  CHKERR m_field.getInterface(prb_mng_ptr);
104  CHKERR prb_mng_ptr->buildProblemOnDistributedMesh("P1", true);
105  CHKERR prb_mng_ptr->partitionFiniteElements("P1");
106  CHKERR prb_mng_ptr->partitionGhostDofs("P1");
107 
108  auto get_triangles_on_skin = [&](Range &tets_skin) {
110  Range tets;
111  CHKERR m_field.get_moab().get_entities_by_dimension(root_set, SPACE_DIM,
112  tets);
113  Skinner skin(&m_field.get_moab());
114  CHKERR skin.find_skin(0, tets, false, tets_skin);
115  Range adj;
116  CHKERR m_field.get_moab().get_adjacencies(tets_skin, SPACE_DIM - 2, false,
117  adj, moab::Interface::UNION);
118  tets_skin.merge(adj);
120  };
121 
122  Range tets;
123  CHKERR m_field.get_moab().get_entities_by_dimension(root_set, SPACE_DIM,
124  tets);
125  SmartPetscObj<IS> is_before_remove;
126  CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
127  "P1", ROW, "F1", 0, 1, is_before_remove, &tets);
128 
129  Range tets_skin;
130  CHKERR get_triangles_on_skin(tets_skin);
131  CHKERR prb_mng_ptr->removeDofsOnEntities("P1", "F1", tets_skin, 0, 1, 0,
132  100, NOISY, true);
133 
134  SmartPetscObj<IS> is_after_remove;
135  CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
136  "P1", ROW, "F1", 0, 1, is_after_remove, &tets);
137 
138  auto test_is = [&](auto prb_name, auto is_before, auto is_after) {
140  const Problem *prb_ptr = m_field.get_problem(prb_name);
141  if (auto sub_data = prb_ptr->getSubData()) {
142  auto sub_ao = sub_data->getSmartRowMap();
143  if (sub_ao) {
144  CHKERR AOApplicationToPetscIS(sub_ao, is_before);
145  PetscBool is_the_same;
146  CHKERR ISEqual(is_before, is_after, &is_the_same);
147  if (is_the_same == PETSC_FALSE) {
148  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
149  "IS should be the same if map is correctly implemented");
150  } else {
151  MOFEM_LOG("WORLD", Sev::inform) << "Sub data map is correct";
152  }
153  } else {
154  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
155  "AO map should exist");
156  }
157  } else {
158  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
159  "Sub DM should exist");
160  }
162  };
163 
164  CHKERR test_is("P1", is_before_remove, is_after_remove);
165 
167  ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("P1", -2, -2,
168  0);
169 
170  CHKERR m_field.add_problem("SUB");
171  // set refinement level for problem
172  CHKERR m_field.modify_problem_ref_level_add_bit("SUB", bit_level);
173  CHKERR m_field.modify_problem_add_finite_element("SUB", "E1");
174  CHKERR prb_mng_ptr->buildSubProblem("SUB", {"F1"}, {"F1"}, "P1", PETSC_TRUE,
175  nullptr, nullptr);
176  CHKERR prb_mng_ptr->partitionFiniteElements("SUB", true, 0,
177  m_field.get_comm_size());
178  CHKERR prb_mng_ptr->partitionGhostDofsOnDistributedMesh("SUB");
179 
180  CHKERR prb_mng_ptr->removeDofsOnEntities("SUB", "F1", tets_skin, 0, 1, 0,
181  100, NOISY, true);
182  SmartPetscObj<IS> is_sub_after_remove;
183  CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
184  "SUB", ROW, "F1", 0, MAX_DOFS_ON_ENTITY, is_sub_after_remove, &tets);
185  CHKERR test_is("SUB", is_before_remove, is_after_remove);
186 
188  ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("SUB", -2,
189  -2, 0);
190  }
191  CATCH_ERRORS;
192 
193  // finish work cleaning memory, getting statistics, etc.
195 
196  return 0;
197 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static

◆ SPACE_DIM

constexpr int SPACE_DIM = 3
constexpr
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
SPACE_DIM
constexpr int SPACE_DIM
Definition: remove_entities_from_problem.cpp:12
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
NOISY
@ NOISY
Definition: definitions.h:211
MoFEM::CoreInterface::modify_finite_element_add_field_row
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
MoFEM::CoreInterface::get_problem
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::CoreInterface::add_ents_to_field_by_dim
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.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
ROW
@ ROW
Definition: definitions.h:123
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::add_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
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::CoreInterface::modify_finite_element_add_field_col
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
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::ProblemsManager::partitionGhostDofsOnDistributedMesh
MoFEMErrorCode partitionGhostDofsOnDistributedMesh(const std::string name, int verb=VERBOSE)
determine ghost nodes on distributed meshes
Definition: ProblemsManager.cpp:2473
MoFEM::CoreInterface::add_field
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.
MoFEM::ProblemsManager::removeDofsOnEntities
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.
Definition: ProblemsManager.cpp:2605
MoFEM::Problem::getSubData
boost::shared_ptr< SubProblemData > & getSubData() const
Get main problem of sub-problem is.
Definition: ProblemsMultiIndices.hpp:119
MoFEM::CoreInterface::get_comm_size
virtual int get_comm_size() const =0
MoFEM::CoreInterface::add_ents_to_finite_element_by_dim
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
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
MoFEM::ProblemsManager::partitionFiniteElements
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
Definition: ProblemsManager.cpp:2167
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
MAX_DOFS_ON_ENTITY
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
Definition: definitions.h:236
Range
MoFEM::CoreTmp< 0 >::Initialize
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
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
help
static char help[]
Definition: remove_entities_from_problem.cpp:10
MoFEM::ProblemsManager::buildProblemOnDistributedMesh
MoFEMErrorCode buildProblemOnDistributedMesh(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures, assuming that mesh is distributed (collective)
Definition: ProblemsManager.cpp:290
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
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
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
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::CoreInterface::set_field_order
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.
MoFEM::SmartPetscObj< IS >
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::ProblemsManager::buildSubProblem
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
Definition: ProblemsManager.cpp:981