v0.14.0
partition_mesh.cpp

Atom testing mesh partitioning

/** \file partition_mesh.cpp
\example partition_mesh.cpp
\brief Atom testing mesh partitioning
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
// initialize petsc
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
int nb_vertices;
try {
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 != PETSC_TRUE)
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
// read mesh and create moab and mofem datastrutures
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
{
const char *option;
option = "";
CHKERR moab.load_file(mesh_file_name, 0, option);
}
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
CHKERR moab.get_number_entities_by_dimension(0, 0, nb_vertices, true);
EntityHandle root_set = moab.get_root_set();
Range tets;
moab.get_entities_by_type(root_set, MBTET, tets, false);
Tag th_vertex_weight;
int def_val = 1;
CHKERR moab.tag_get_handle("VERTEX_WEIGHT", 1, MB_TYPE_INTEGER,
th_vertex_weight, MB_TAG_CREAT | MB_TAG_DENSE,
&def_val);
CommInterface *comm_interafce_ptr = m_field.getInterface<CommInterface>();
CHKERR comm_interafce_ptr->partitionMesh(
tets, 3, 2, m_field.get_comm_size(), &th_vertex_weight, NULL, NULL,
VERBOSE, false);
if (!m_field.get_comm_rank()) {
CHKERR moab.write_file("partitioned_mesh.h5m");
}
}
PetscBarrier(PETSC_NULL);
try {
moab::Core mb_instance2;
moab::Interface &moab2 = mb_instance2;
MoFEM::CoreTmp<1> core(moab2);
MoFEM::Interface &m_field = core;
// Register DM Manager
DMType dm_name = "DMMOFEM";
// Test build simple problem
const char *option = "DEBUG_IO;"
"PARALLEL=READ_PART;"
"PARALLEL_RESOLVE_SHARED_ENTS;"
"PARTITION=PARALLEL_PARTITION;";
CHKERR m_field.getInterface<Simple>()->getOptions();
CHKERR m_field.getInterface<Simple>()->loadFile(option,
"partitioned_mesh.h5m");
CHKERR m_field.getInterface<Simple>()->addDomainField(
CHKERR m_field.getInterface<Simple>()->setFieldOrder("U", 1);
CHKERR m_field.getInterface<Simple>()->setUp();
auto dm = m_field.getInterface<Simple>()->getDM();
const MoFEM::Problem *problem_ptr;
CHKERR DMMoFEMGetProblemPtr(dm, &problem_ptr);
if(problem_ptr->nbDofsRow != nb_vertices)
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Number of vertices and DOFs is inconstent");
MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "Atom test")
<< "All is good in this test";
}
return 0;
}
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
EntityHandle
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
MoFEM.hpp
MoFEM::CommInterface::partitionMesh
MoFEMErrorCode partitionMesh(const Range &ents, const int dim, const int adj_dim, const int n_parts, Tag *th_vertex_weights=nullptr, Tag *th_edge_weights=nullptr, Tag *th_part_weights=nullptr, int verb=VERBOSE, const bool debug=false)
Set partition tag to each finite element in the problem.
Definition: CommInterface.cpp:868
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
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
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
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreTmp
Definition: Core.hpp:36
help
static char help[]
Definition: partition_mesh.cpp:13
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::CoreInterface::get_comm_size
virtual int get_comm_size() const =0
main
int main(int argc, char *argv[])
Definition: partition_mesh.cpp:15
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
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_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::Problem::nbDofsRow
DofIdx nbDofsRow
Global number of DOFs in row.
Definition: ProblemsMultiIndices.hpp:65
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::DMMoFEMGetProblemPtr
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:426
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54