v0.14.0
boundary_marker.cpp

Mark skin entities, i.e. boundary, next check DOFs if are properly marked when accessed form user data operator.

/**
* \file boundary_marker.cpp
* \example boundary_marker.cpp
*
* Mark skin entities, i.e. boundary, next check DOFs if are properly marked
* when accessed form user data operator.
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
struct OpFace : public FaceEleOp {
const Range &skinEnts;
const std::vector<unsigned char> &mArker;
OpFace(const Range &skin_ents, const std::vector<unsigned char> &marker)
: FaceEleOp("FIELD1", OPROW), skinEnts(skin_ents), mArker(marker) {}
MoFEMErrorCode doWork(int side, EntityType type,
const int nb_dofs = data.getIndices().size();
if (nb_dofs == 0)
// This function takes entities on DOFs, and check if those entities are
// contained in the range. Note it is local element vector, which is used
// to validate of global local vector.
auto get_boundary_marker_directly_from_range = [&]() {
std::vector<unsigned char> ents_marker_used_for_testing;
ents_marker_used_for_testing.resize(data.getLocalIndices().size(),0);
switch (type) {
case MBVERTEX: {
for (size_t side = 0; side != getNumberOfNodesOnElement(); ++side) {
EntityHandle ent = getSideEntity(side, MBVERTEX);
ents_marker_used_for_testing[3 * side + 0] =
skinEnts.find(ent) != skinEnts.end() ? 1 : 0;
ents_marker_used_for_testing[3 * side + 2] =
ents_marker_used_for_testing[3 * side + 0];
}
break;
}
default: {
EntityHandle ent = getSideEntity(side, type);
bool is_on_boundary = skinEnts.find(ent) != skinEnts.end();
for (size_t dof = 0; dof != nb_dofs; ++dof)
if ((dof % 3) != 1)
ents_marker_used_for_testing[dof] = is_on_boundary ? 1 : 0;
}
};
return ents_marker_used_for_testing;
};
auto test_marker = get_boundary_marker_directly_from_range();
for (size_t n = 0; n != nb_dofs; ++n) {
const size_t local_index = data.getLocalIndices()[n];
if (test_marker[n] != mArker[local_index])
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong boundary marker");
}
}
};
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
// Create MoFEM instance
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
int order = 4;
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
DMType dm_name = "DMMOFEM";
auto *simple_interface = m_field.getInterface<Simple>();
CHKERR simple_interface->getOptions();
CHKERR simple_interface->loadFile("");
CHKERR simple_interface->addDomainField("FIELD1", H1,
CHKERR simple_interface->setFieldOrder("FIELD1", order);
CHKERR simple_interface->setUp();
auto get_ents_on_mesh_skin = [&]() {
Range faces;
CHKERR m_field.get_moab().get_entities_by_type(0, MBTRI, faces);
Skinner skin(&m_field.get_moab());
Range skin_edges;
CHKERR skin.find_skin(0, faces, false, skin_edges);
Range skin_verts;
CHKERR moab.get_connectivity(skin_edges, skin_verts, true);
skin_edges.merge(skin_verts);
return skin_edges;
};
auto mark_boundary_dofs = [&](Range &skin_edges) {
auto problem_manager = m_field.getInterface<ProblemsManager>();
std::vector<unsigned char> marker;
problem_manager->markDofs(simple_interface->getProblemName(), ROW,
ProblemsManager::OR, skin_edges, marker);
// Unset coefficient 1, e.g. u_y
problem_manager->modifyMarkDofs(simple_interface->getProblemName(), ROW,
"FIELD1", 1, 1, ProblemsManager::AND, 0,
return marker;
};
auto skin_ents = get_ents_on_mesh_skin();
// Get global local vector of marked DOFs. Is global, since is set for all
// DOFs on processor. Is local since only DOFs on processor are in the
// vector. To access DOFs use local indices.
auto marker = mark_boundary_dofs(skin_ents);
boost::shared_ptr<FaceEle> fe(new FaceEle(m_field));
fe->getOpPtrVector().push_back(new OpFace(skin_ents, marker));
auto dm = simple_interface->getDM();
CHKERR DMoFEMLoopFiniteElements(dm, simple_interface->getDomainFEName(),
fe);
}
}
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::EntitiesFieldData::EntData::getLocalIndices
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
Definition: EntitiesFieldData.hpp:1216
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::ProblemsManager::OR
@ OR
Definition: ProblemsManager.hpp:365
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::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
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::FaceElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: FaceElementForcesAndSourcesCore.hpp:86
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
convert.type
type
Definition: convert.py:64
help
static char help[]
Definition: boundary_marker.cpp:14
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1201
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
OpFace
Definition: boundary_marker.cpp:19
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::ProblemsManager::AND
@ AND
Definition: ProblemsManager.hpp:365
convert.n
n
Definition: convert.py:82
FaceEle
MoFEM::FaceElementForcesAndSourcesCore FaceEle
Definition: boundary_marker.cpp:16
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
main
int main(int argc, char *argv[])
Definition: boundary_marker.cpp:74
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
marker
auto marker
set bit to marker
Definition: hanging_node_approx.cpp:82
FaceEleOp
FaceEle::UserDataOperator FaceEleOp
Definition: boundary_marker.cpp:17
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:590
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346