v0.14.0
remove_entities_from_problem_not_partitioned.cpp
Go to the documentation of this file.
1 /** \file remove_entities_from_problem_not_partitioned.cpp
2  \example remove_entities_from_problem_not_partitioned.cpp
3  \brief Remove field entities from the problem
4 
5 */
6 
7 #include <MoFEM.hpp>
8 using namespace MoFEM;
9 
10 static char help[] = "...\n\n";
11 
12 constexpr int SPACE_DIM = 3;
13 
14 int main(int argc, char *argv[]) {
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 
38  const char *option;
39  option = "";
40  CHKERR moab.load_file(mesh_file_name, 0, option);
41 
42  // Create MoFEM database
43  MoFEM::Core core(moab);
44  MoFEM::Interface &m_field = core;
45 
46  MOFEM_LOG_CHANNEL("WORLD");
47  BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
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 
69  for (EntityType t = CN::TypeDimensionMap[1].first;
70  t <= CN::TypeDimensionMap[SPACE_DIM].second; ++t) {
71  CHKERR m_field.set_field_order(root_set, t, "F1", order);
72  }
73  for (EntityType t = CN::TypeDimensionMap[2].first;
74  t <= CN::TypeDimensionMap[SPACE_DIM].second; ++t) {
75  CHKERR m_field.set_field_order(root_set, t, "F2", order);
76  }
77 
78  CHKERR m_field.build_fields();
79 
80  // add elements
81  CHKERR m_field.add_finite_element("E1");
82 
83  CHKERR m_field.modify_finite_element_add_field_row("E1", "F1");
84  CHKERR m_field.modify_finite_element_add_field_row("E1", "F2");
85  CHKERR m_field.modify_finite_element_add_field_col("E1", "F2");
86  CHKERR m_field.modify_finite_element_add_field_col("E1", "F1");
87  CHKERR m_field.modify_finite_element_add_field_data("E1", "F1");
88  CHKERR m_field.modify_finite_element_add_field_data("E1", "F2");
89 
90  CHKERR m_field.add_ents_to_finite_element_by_dim(root_set, SPACE_DIM, "E1");
91 
92  CHKERR m_field.build_finite_elements();
93  CHKERR m_field.build_adjacencies(bit_level);
94 
95  // Problems
96  CHKERR m_field.add_problem("P1");
97 
98  // set refinement level for problem
99  CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_level);
100  CHKERR m_field.modify_problem_add_finite_element("P1", "E1");
101 
102  // build problems
103  ProblemsManager *prb_mng_ptr;
104  CHKERR m_field.getInterface(prb_mng_ptr);
105  CHKERR prb_mng_ptr->buildProblem("P1", true);
106  CHKERR prb_mng_ptr->partitionProblem("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  MOFEM_LOG("WORLD", Sev::inform) << "Remove dofs";
130 
131  Range tets_skin;
132  CHKERR get_triangles_on_skin(tets_skin);
134  "P1", "F1", tets_skin, 0, 1, 0, 100, VERBOSE, true);
135 
136  MOFEM_LOG("WORLD", Sev::inform) << "Check consistency";
137 
138  SmartPetscObj<IS> is_after_remove;
139  CHKERR m_field.getInterface<ISManager>()->isCreateProblemFieldAndRank(
140  "P1", ROW, "F1", 0, 1, is_after_remove, &tets);
141 
142  auto test_is = [&]() {
144  const Problem *prb_ptr = m_field.get_problem("P1");
145  if (auto sub_data = prb_ptr->getSubData()) {
146  auto sub_ao = sub_data->getSmartRowMap();
147  if (sub_ao) {
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) {
152  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
153  "IS should be the same if map is correctly implemented");
154  } else {
155  MOFEM_LOG("WORLD", Sev::inform) << "Sub data map is correct";
156  }
157  } else {
158  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
159  "AO map should exist");
160  }
161  } else {
162  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
163  "Sub DM should exist");
164  }
166  };
167 
168  CHKERR test_is();
169 
171  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("P1", ROW, v);
172 
174  ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>("P1", -2,
175  -2,
176  0);
177 
178  MOFEM_LOG("WORLD", Sev::inform) << "done";
179 
180  }
181  CATCH_ERRORS;
182 
183  // finish work cleaning memory, getting statistics, etc.
185 
186  return 0;
187 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
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.hpp
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:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
VERBOSE
@ VERBOSE
Definition: definitions.h:222
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
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::ProblemsManager::removeDofsOnEntitiesNotDistributed
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.
Definition: ProblemsManager.cpp:3226
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
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::Problem::getSubData
boost::shared_ptr< SubProblemData > & getSubData() const
Get main problem of sub-problem is.
Definition: ProblemsMultiIndices.hpp:119
help
static char help[]
Definition: remove_entities_from_problem_not_partitioned.cpp:10
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:58
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_field)=0
set finite element field data
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
MoFEM::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
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:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
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
SPACE_DIM
constexpr int SPACE_DIM
Definition: remove_entities_from_problem_not_partitioned.cpp:12
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::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:429
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
main
int main(int argc, char *argv[])
Definition: remove_entities_from_problem_not_partitioned.cpp:14
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.
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::ProblemsManager::partitionProblem
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
Definition: ProblemsManager.cpp:1683