v0.13.2
Loading...
Searching...
No Matches
dm_partitioned_no_field.cpp
Go to the documentation of this file.
1/** \file dm_build_partitioned_mesh.cpp
2 \example dm_partitioned_no_field.cpp
3 \brief Testing problem for partitioned mesh with NOFIELD field
4
5*/
6
7#include <MoFEM.hpp>
8
9using namespace MoFEM;
10
11static char help[] = "...\n\n";
12constexpr bool debug = false;
13
14int main(int argc, char *argv[]) {
15 // initialize petsc
16 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
17
18 try {
19
20 PetscBool flg = PETSC_TRUE;
21 char mesh_file_name[255];
22 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
23 255, &flg);
24 if (flg != PETSC_TRUE)
25 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
26
27 // register new dm type, i.e. mofem
28 DMType dm_name = "DMMOFEM";
29 CHKERR DMRegister_MoFEM(dm_name);
30
31 // create dm instance
32 DM dm;
33 CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
34 CHKERR DMSetType(dm, dm_name);
35
36 // read mesh and create moab and mofem data structures
37 moab::Core mb_instance;
38 moab::Interface &moab = mb_instance;
39
40 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
41 auto moab_comm_wrap =
42 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
43 if (pcomm == NULL)
44 pcomm = new ParallelComm(&moab, moab_comm_wrap->get_comm());
45
46 const std::string options = "PARALLEL=READ_PART;"
47 "PARALLEL_RESOLVE_SHARED_ENTS;"
48 "PARTITION=PARALLEL_PARTITION;";
49 CHKERR moab.load_file(mesh_file_name, 0, options.c_str());
50
51 MoFEM::Core core(moab, PETSC_COMM_WORLD);
52 MoFEM::Interface &m_field = core;
53
54 auto *bit_ref_ptr = m_field.getInterface<BitRefManager>();
55 auto *comm_interface_ptr = m_field.getInterface<CommInterface>();
56
57 EntityHandle root_set = moab.get_root_set();
58 // add all entities to database, all of them will be used
59 CHKERR bit_ref_ptr->setBitRefLevelByDim(root_set, 3, BitRefLevel().set(0));
60
61 // add field
62 CHKERR m_field.add_field("FIELD", H1, AINSWORTH_LEGENDRE_BASE, 1);
63 CHKERR m_field.add_field("LAMBDA", NOFIELD, NOBASE, 3);
64
65 // add entities to field
66 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD");
67
68 // set app. order
69 CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD", 1);
70
71 // Create vertices for NOFILE
72 std::array<double, 6> coords = {0, 0, 0, 0, 0, 0};
73 CHKERR m_field.create_vertices_and_add_to_field("LAMBDA", coords.data(), 2);
74 CHKERR comm_interface_ptr->makeFieldEntitiesMultishared("LAMBDA", 0,
75 NOISY);
76 CHKERR bit_ref_ptr->setFieldEntitiesBitRefLevel("LAMBDA",
77 BitRefLevel().set(0));
78
79 // build data structures for fields
80 CHKERR m_field.build_fields();
81
82 auto print_field_ents = [&](const std::string field_name) {
83 auto *field_ents = m_field.get_field_ents();
84 auto field_bit_number = m_field.get_field_bit_number(field_name);
85 auto lo = field_ents->get<Unique_mi_tag>().lower_bound(
86 FieldEntity::getLoBitNumberUId(field_bit_number));
87 auto hi = field_ents->get<Unique_mi_tag>().upper_bound(
88 FieldEntity::getHiBitNumberUId(field_bit_number));
89 for (auto it = lo; it != hi; ++it) {
90 std::ostringstream ss;
91 ss << "Rank " << m_field.get_comm_rank() << " -> " << **it << endl;
92 PetscSynchronizedPrintf(m_field.get_comm(), "%s", ss.str().c_str());
93 PetscSynchronizedFlush(m_field.get_comm(), PETSC_STDOUT);
94 }
95 };
96
97 print_field_ents("LAMBDA");
98
99 if (m_field.get_comm_rank() == 0) {
100 CHKERR m_field.getInterface<FieldBlas>()->setField(1, "LAMBDA");
101 CHKERR m_field.getInterface<FieldBlas>()->setField(1, "FIELD");
102 }
103
104 CHKERR comm_interface_ptr->exchangeFieldData("LAMBDA");
105 CHKERR comm_interface_ptr->exchangeFieldData("FIELD");
106
107 auto check_exchanged_values = [&](const std::string field_name) {
109 if (m_field.get_comm_rank() != 0) {
110 auto *field_ents = m_field.get_field_ents();
111 auto field_bit_number = m_field.get_field_bit_number(field_name);
112 auto lo = field_ents->get<Unique_mi_tag>().lower_bound(
113 FieldEntity::getLoBitNumberUId(field_bit_number));
114 auto hi = field_ents->get<Unique_mi_tag>().upper_bound(
115 FieldEntity::getHiBitNumberUId(field_bit_number));
116 for (auto it = lo; it != hi; ++it) {
117 VectorDouble field_data = (*it)->getEntFieldData();
118 for (auto v : field_data)
119 if (v != 1)
120 SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
121 "Wrong value on field %4.3f", v);
122 }
123 }
125 };
126
127 CHKERR check_exchanged_values("LAMBDA");
128
129 // define & build finite elements
130 CHKERR m_field.add_finite_element("FE");
131 // Define rows/cols and element data
132 CHKERR m_field.modify_finite_element_add_field_row("FE", "FIELD");
133 CHKERR m_field.modify_finite_element_add_field_col("FE", "FIELD");
134 CHKERR m_field.modify_finite_element_add_field_data("FE", "FIELD");
135 CHKERR m_field.modify_finite_element_add_field_row("FE", "LAMBDA");
136 CHKERR m_field.modify_finite_element_add_field_col("FE", "LAMBDA");
137 CHKERR m_field.modify_finite_element_add_field_data("FE", "LAMBDA");
138
139 // add entities to finite element
140 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE");
141 // build finite elemnts
143 // build adjacencies
144 CHKERR m_field.build_adjacencies(BitRefLevel().set());
145
146 // set dm data structure which created mofem data structures
147 CHKERR DMMoFEMCreateMoFEM(dm, &m_field, "TEST_PROBLEM",
148 BitRefLevel().set(0));
150 dm, PETSC_FALSE); // this is for testing (this problem has the same rows
151 // and cols)
152 CHKERR DMSetFromOptions(dm);
153 CHKERR DMMoFEMAddElement(dm, "FE");
154 CHKERR DMSetUp(dm);
155
156 Mat m;
157 CHKERR DMCreateMatrix(dm, &m);
158
160 ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>(
161 "TEST_PROBLEM", -1, -1, 1);
162
163 CHKERR MatDestroy(&m);
164
165 // check if file can be saved
166 if (debug)
167 CHKERR moab.write_file("test_out.h5m", "MOAB", "PARALLEL=WRITE_PART");
168
169 // destry dm
170 CHKERR DMDestroy(&dm);
171 }
173
174 // finish work cleaning memory, getting statistics, ect.
176
177 return 0;
178}
int main()
Definition: adol-c_atom.cpp:46
@ NOISY
Definition: definitions.h:211
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ 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
Definition: definitions.h:215
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
static char help[]
constexpr bool debug
FTensor::Index< 'm', SPACE_DIM > m
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:485
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition: DMMoFEM.cpp:444
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMoFEM.cpp:118
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
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]
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
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
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
Matrix manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.