v0.14.0
Loading...
Searching...
No Matches
dm_partitioned_no_field.cpp

Testing problem for partitioned mesh with NOFIELD field.

Testing problem for partitioned mesh with NOFIELD field

/** \file dm_build_partitioned_mesh.cpp
\example dm_partitioned_no_field.cpp
\brief Testing problem for partitioned mesh with NOFIELD field
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
constexpr bool debug = false;
int main(int argc, char *argv[]) {
// initialize petsc
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
if (flg != PETSC_TRUE)
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
// register new dm type, i.e. mofem
DMType dm_name = "DMMOFEM";
CHKERR DMRegister_MoFEM(dm_name);
// create dm instance
DM dm;
CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
CHKERR DMSetType(dm, dm_name);
// read mesh and create moab and mofem data structures
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
auto moab_comm_wrap =
boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
if (pcomm == NULL)
pcomm = new ParallelComm(&moab, moab_comm_wrap->get_comm());
const std::string options = "PARALLEL=READ_PART;"
"PARALLEL_RESOLVE_SHARED_ENTS;"
"PARTITION=PARALLEL_PARTITION;";
CHKERR moab.load_file(mesh_file_name, 0, options.c_str());
MoFEM::Core core(moab, PETSC_COMM_WORLD);
MoFEM::Interface &m_field = core;
auto *bit_ref_ptr = m_field.getInterface<BitRefManager>();
auto *comm_interface_ptr = m_field.getInterface<CommInterface>();
EntityHandle root_set = moab.get_root_set();
// add all entities to database, all of them will be used
CHKERR bit_ref_ptr->setBitRefLevelByDim(root_set, 3, BitRefLevel().set(0));
// add field
CHKERR m_field.add_field("FIELD", H1, AINSWORTH_LEGENDRE_BASE, 1);
CHKERR m_field.add_field("LAMBDA", NOFIELD, NOBASE, 3);
// add entities to field
CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD");
// set app. order
CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD", 1);
// Create vertices for NOFILE
std::array<double, 6> coords = {0, 0, 0, 0, 0, 0};
CHKERR m_field.create_vertices_and_add_to_field("LAMBDA", coords.data(), 2);
CHKERR comm_interface_ptr->makeFieldEntitiesMultishared("LAMBDA", 0,
CHKERR bit_ref_ptr->setFieldEntitiesBitRefLevel("LAMBDA",
BitRefLevel().set(0));
// build data structures for fields
CHKERR m_field.build_fields();
auto print_field_ents = [&](const std::string field_name) {
auto *field_ents = m_field.get_field_ents();
auto field_bit_number = m_field.get_field_bit_number(field_name);
auto lo = field_ents->get<Unique_mi_tag>().lower_bound(
FieldEntity::getLoBitNumberUId(field_bit_number));
auto hi = field_ents->get<Unique_mi_tag>().upper_bound(
FieldEntity::getHiBitNumberUId(field_bit_number));
for (auto it = lo; it != hi; ++it) {
std::ostringstream ss;
ss << "Rank " << m_field.get_comm_rank() << " -> " << **it << endl;
PetscSynchronizedPrintf(m_field.get_comm(), "%s", ss.str().c_str());
PetscSynchronizedFlush(m_field.get_comm(), PETSC_STDOUT);
}
};
print_field_ents("LAMBDA");
if (m_field.get_comm_rank() == 0) {
CHKERR m_field.getInterface<FieldBlas>()->setField(1, "LAMBDA");
CHKERR m_field.getInterface<FieldBlas>()->setField(1, "FIELD");
}
CHKERR comm_interface_ptr->exchangeFieldData("LAMBDA");
CHKERR comm_interface_ptr->exchangeFieldData("FIELD");
auto check_exchanged_values = [&](const std::string field_name) {
if (m_field.get_comm_rank() != 0) {
auto *field_ents = m_field.get_field_ents();
auto field_bit_number = m_field.get_field_bit_number(field_name);
auto lo = field_ents->get<Unique_mi_tag>().lower_bound(
FieldEntity::getLoBitNumberUId(field_bit_number));
auto hi = field_ents->get<Unique_mi_tag>().upper_bound(
FieldEntity::getHiBitNumberUId(field_bit_number));
for (auto it = lo; it != hi; ++it) {
VectorDouble field_data = (*it)->getEntFieldData();
for (auto v : field_data)
if (v != 1)
SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong value on field %4.3f", v);
}
}
};
CHKERR check_exchanged_values("LAMBDA");
// define & build finite elements
CHKERR m_field.add_finite_element("FE");
// Define rows/cols and element data
// add entities to finite element
CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE");
// build finite elemnts
// build adjacencies
CHKERR m_field.build_adjacencies(BitRefLevel().set());
// set dm data structure which created mofem data structures
CHKERR DMMoFEMCreateMoFEM(dm, &m_field, "TEST_PROBLEM",
BitRefLevel().set(0));
CHKERR DMMoFEMSetSquareProblem(
dm, PETSC_FALSE); // this is for testing (this problem has the same rows
// and cols)
CHKERR DMSetFromOptions(dm);
CHKERR DMMoFEMAddElement(dm, "FE");
CHKERR DMSetUp(dm);
Mat m;
CHKERR DMCreateMatrix(dm, &m);
->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>(
"TEST_PROBLEM", -1, -1, 1);
CHKERR MatDestroy(&m);
// check if file can be saved
if (debug)
CHKERR moab.write_file("test_out.h5m", "MOAB", "PARALLEL=WRITE_PART");
// destry dm
CHKERR DMDestroy(&dm);
}
// finish work cleaning memory, getting statistics, ect.
return 0;
}
static char help[]
int main()
@ NOISY
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ NOBASE
Definition definitions.h:59
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition definitions.h:84
@ H1
continuous field
Definition definitions.h:85
#define MYPCOMM_INDEX
default communicator number PCOMM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
static const bool debug
FTensor::Index< 'm', SPACE_DIM > m
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
MoFEMErrorCode setFieldEntitiesBitRefLevel(const std::string field_name, const BitRefLevel bit=BitRefLevel(), int verb=QUIET) const
Set the bit ref level to entities in the field meshset.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
const double v
phase velocity of light in medium (cm/ns)
char mesh_file_name[255]
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr auto field_name
Managing BitRefLevels.
Managing BitRefLevels.
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MPI_Comm & get_comm() const =0
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.
virtual int get_comm_rank() const =0
virtual MoFEMErrorCode create_vertices_and_add_to_field(const std::string name, const double coords[], int size, int verb=DEFAULT_VERBOSITY)=0
Create a vertices and add to field object.
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
Basic algebra on fields.
Definition FieldBlas.hpp:21
Matrix manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.