v0.13.0
boundary_marker.cpp
Go to the documentation of this file.
1 /**
2  * \file boundary_marker.cpp
3  * \example boundary_marker.cpp
4  *
5  * Mark skin entities, i.e. boundary, next check DOFs if are properly marked
6  * when accessed form user data operator.
7  *
8  */
9 
10 /* This file is part of MoFEM.
11  * MoFEM is free software: you can redistribute it and/or modify it under
12  * the terms of the GNU Lesser General Public License as published by the
13  * Free Software Foundation, either version 3 of the License, or (at your
14  * option) any later version.
15  *
16  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
17  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
19  * License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
23 
24 #include <MoFEM.hpp>
25 
26 using namespace MoFEM;
27 
28 static char help[] = "...\n\n";
29 
32 
33 struct OpFace : public FaceEleOp {
34  const Range &skinEnts;
35  const std::vector<unsigned char> &mArker;
36 
37  OpFace(const Range &skin_ents, const std::vector<unsigned char> &marker)
38  : FaceEleOp("FIELD1", OPROW), skinEnts(skin_ents), mArker(marker) {}
39 
40  MoFEMErrorCode doWork(int side, EntityType type,
43 
44  const int nb_dofs = data.getIndices().size();
45  if (nb_dofs == 0)
47 
48  // This function takes entities on DOFs, and check if those entities are
49  // contained in the range. Note it is local element vector, which is used
50  // to validate of global local vector.
51  auto get_boundary_marker_directly_from_range = [&]() {
52  std::vector<unsigned char> ents_marker_used_for_testing;
53  ents_marker_used_for_testing.resize(data.getLocalIndices().size(),0);
54  switch (type) {
55  case MBVERTEX: {
56  for (size_t side = 0; side != getNumberOfNodesOnElement(); ++side) {
57  EntityHandle ent = getSideEntity(side, MBVERTEX);
58  ents_marker_used_for_testing[3 * side + 0] =
59  skinEnts.find(ent) != skinEnts.end() ? 1 : 0;
60  ents_marker_used_for_testing[3 * side + 2] =
61  ents_marker_used_for_testing[3 * side + 0];
62  }
63  break;
64  }
65  default: {
66  EntityHandle ent = getSideEntity(side, type);
67  bool is_on_boundary = skinEnts.find(ent) != skinEnts.end();
68  for (size_t dof = 0; dof != nb_dofs; ++dof)
69  if ((dof % 3) != 1)
70  ents_marker_used_for_testing[dof] = is_on_boundary ? 1 : 0;
71  }
72  };
73  return ents_marker_used_for_testing;
74  };
75 
76  auto test_marker = get_boundary_marker_directly_from_range();
77  for (size_t n = 0; n != nb_dofs; ++n) {
78  const size_t local_index = data.getLocalIndices()[n];
79  if (test_marker[n] != mArker[local_index])
80  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
81  "Wrong boundary marker");
82  }
83 
85  }
86 };
87 
88 int main(int argc, char *argv[]) {
89 
90  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
91 
92  try {
93 
94  moab::Core mb_instance;
95  moab::Interface &moab = mb_instance;
96 
97  // Create MoFEM instance
98  MoFEM::Core core(moab);
99  MoFEM::Interface &m_field = core;
100 
101  int order = 4;
102  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
103 
104  DMType dm_name = "DMMOFEM";
105  CHKERR DMRegister_MoFEM(dm_name);
106 
107  auto *simple_interface = m_field.getInterface<Simple>();
108  CHKERR simple_interface->getOptions();
109  CHKERR simple_interface->loadFile("");
110  CHKERR simple_interface->addDomainField("FIELD1", H1,
112  CHKERR simple_interface->setFieldOrder("FIELD1", order);
113  CHKERR simple_interface->setUp();
114 
115  auto get_ents_on_mesh_skin = [&]() {
116  Range faces;
117  CHKERR m_field.get_moab().get_entities_by_type(0, MBTRI, faces);
118  Skinner skin(&m_field.get_moab());
119  Range skin_edges;
120  CHKERR skin.find_skin(0, faces, false, skin_edges);
121  Range skin_verts;
122  CHKERR moab.get_connectivity(skin_edges, skin_verts, true);
123  skin_edges.merge(skin_verts);
124  return skin_edges;
125  };
126 
127  auto mark_boundary_dofs = [&](Range &skin_edges) {
128  auto problem_manager = m_field.getInterface<ProblemsManager>();
129  std::vector<unsigned char> marker;
130  problem_manager->markDofs(simple_interface->getProblemName(), ROW,
131  ProblemsManager::OR, skin_edges, marker);
132  // Unset coefficient 1, e.g. u_y
133  problem_manager->modifyMarkDofs(simple_interface->getProblemName(), ROW,
134  "FIELD1", 1, 1, ProblemsManager::AND, 0,
135  marker);
136  return marker;
137  };
138 
139  auto skin_ents = get_ents_on_mesh_skin();
140 
141  // Get global local vector of marked DOFs. Is global, since is set for all
142  // DOFs on processor. Is local since only DOFs on processor are in the
143  // vector. To access DOFs use local indices.
144  auto marker = mark_boundary_dofs(skin_ents);
145 
146  boost::shared_ptr<FaceEle> fe(new FaceEle(m_field));
147  fe->getOpPtrVector().push_back(new OpFace(skin_ents, marker));
148 
149  auto dm = simple_interface->getDM();
150  CHKERR DMoFEMLoopFiniteElements(dm, simple_interface->getDomainFEName(),
151  fe);
152  }
153  CATCH_ERRORS;
154 
156 }
static Index< 'n', 3 > n
ForcesAndSourcesCore::UserDataOperator UserDataOperator
int main(int argc, char *argv[])
static char help[]
MoFEM::FaceElementForcesAndSourcesCore FaceEle
@ ROW
Definition: definitions.h:136
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:73
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
@ H1
continuous field
Definition: definitions.h:98
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:53
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:59
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:544
FaceElementForcesAndSourcesCoreSwitch< 0 > FaceElementForcesAndSourcesCore
Face finite element default.
auto marker
set bit to marker
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
CoreTmp< 0 > Core
Definition: Core.hpp:1096
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1965
virtual moab::Interface & get_moab()=0
Core (interface) class.
Definition: Core.hpp:92
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:85
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:125
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition: Simple.hpp:33
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
OpFace(const Range &skin_ents, const std::vector< unsigned char > &marker)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
const std::vector< unsigned char > & mArker
const Range & skinEnts