v0.15.0
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_NULLPTR, "", "-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 SETERRQ(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 // add entities to finite element
139 CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE");
140
141 // add meshset finite element
142 CHKERR m_field.add_finite_element("MESHSET_FE");
143 auto create_meshset = [&]() {
144 EntityHandle fe_meshset;
145 auto create = [&]() {
147 CHKERR m_field.get_moab().create_meshset(MESHSET_SET, fe_meshset);
148 {
149 ParallelComm *pcomm =
150 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
151 Tag part_tag = pcomm->part_tag();
152 int rank = pcomm->rank();
153 CHKERR m_field.get_moab().tag_set_data(part_tag, &fe_meshset, 1, &rank);
154 CHKERR m_field.getInterface<BitRefManager>()->setBitLevelToMeshset(
155 fe_meshset, BitRefLevel().set());
156
157 }
159 };
160 CHK_THROW_MESSAGE(create(), "Failed to create meshset");
161 return fe_meshset;
162 };
163
164 auto add_verts_to_meshset = [&](std::string field_name,
165 EntityHandle meshset) {
167 auto field_meshest = m_field.get_field_meshset(field_name);
168 Range verts;
169 CHKERR m_field.get_moab().get_entities_by_dimension(field_meshest, 0,
170 verts, true);
171 CHKERR m_field.get_moab().add_entities(meshset, verts);
173 };
174
175 auto fe_meshset = create_meshset();
176 CHKERR add_verts_to_meshset("FIELD", fe_meshset);
177
179 "MESHSET_FE", false);
180
181 CHKERR m_field.modify_finite_element_add_field_row("MESHSET_FE", "FIELD");
182 CHKERR m_field.modify_finite_element_add_field_col("MESHSET_FE", "FIELD");
183 CHKERR m_field.modify_finite_element_add_field_data("MESHSET_FE", "FIELD");
184 CHKERR m_field.modify_finite_element_add_field_row("MESHSET_FE", "LAMBDA");
185 CHKERR m_field.modify_finite_element_add_field_col("MESHSET_FE", "LAMBDA");
186 CHKERR m_field.modify_finite_element_add_field_data("MESHSET_FE", "LAMBDA");
187
188 // build finite elements
189 CHKERR m_field.build_finite_elements("FE");
190 CHKERR m_field.build_finite_elements("MESHSET_FE");
191 // build adjacencies
192 CHKERR m_field.build_adjacencies(BitRefLevel().set());
193
194 // set dm data structure which created mofem data structures
195 CHKERR DMMoFEMCreateMoFEM(dm, &m_field, "TEST_PROBLEM",
196 BitRefLevel().set(0));
198 dm, PETSC_FALSE); // this is for testing (this problem has the same rows
199 // and cols)
200 CHKERR DMSetFromOptions(dm);
201 CHKERR DMMoFEMAddElement(dm, "FE");
202 CHKERR DMMoFEMAddElement(dm, "MESHSET_FE");
203 CHKERR DMSetUp(dm);
204
205 Mat m;
206 CHKERR DMCreateMatrix(dm, &m);
207
209 ->checkMPIAIJWithArraysMatrixFillIn<PetscGlobalIdx_mi_tag>(
210 "TEST_PROBLEM", -1, -1, 1);
211
212 CHKERR MatDestroy(&m);
213
214 struct PrintField : public FEMethod {
215
216 PrintField(int bit_number) : bitNumber(bit_number) {}
217
218 MoFEMErrorCode preProcess() { return 0; }
221 MOFEM_LOG("SYNC", Sev::inform) << "Element: getDataDofsPtr";
222 auto data_dofs = getDataDofsPtr();
223 auto dit = data_dofs->get<Unique_mi_tag>().lower_bound(
225 auto hi_dit = data_dofs->get<Unique_mi_tag>().upper_bound(
227 for (; dit != hi_dit; dit++) {
228 MOFEM_LOG("SYNC", Sev::inform)
229 << "Ent " << (*dit)->getEnt() << " dof idx "
230 << (*dit)->getEntDofIdx();
231 }
232 MOFEM_LOG("SYNC", Sev::inform) << "Element: getRowDofsPtr";
233 auto row_dofs = getRowDofsPtr();
234 auto r_dit = row_dofs->get<Unique_mi_tag>().lower_bound(
236 auto hi_r_dit = row_dofs->get<Unique_mi_tag>().upper_bound(
238 for (; r_dit != hi_r_dit; r_dit++) {
239 MOFEM_LOG("SYNC", Sev::inform)
240 << "Dof idx " << (*r_dit)->getEntDofIdx() << " global idx "
241 << (*r_dit)->getPetscGlobalDofIdx() << " local idx "
242 << (*r_dit)->getPetscLocalDofIdx();
243 }
244 MOFEM_LOG("SYNC", Sev::inform) << "Element: getColDofsPtr";
245 auto col_dofs = getColDofsPtr();
246 auto c_dit = col_dofs->get<Unique_mi_tag>().lower_bound(
248 auto hi_c_dit = col_dofs->get<Unique_mi_tag>().upper_bound(
250 for (; c_dit != hi_c_dit; c_dit++) {
251 MOFEM_LOG("SYNC", Sev::inform)
252 << "Dof idx " << (*c_dit)->getEntDofIdx() << " global idx "
253 << (*c_dit)->getPetscGlobalDofIdx() << " local idx "
254 << (*c_dit)->getPetscLocalDofIdx();
255 }
257 }
258 MoFEMErrorCode postProcess() { return 0; }
259
260 private:
261 int bitNumber;
262 };
263
264 PrintField print_no_field(m_field.get_field_bit_number("LAMBDA"));
265
266 MOFEM_LOG("SYNC", Sev::inform) << "Loop elements MESHSET_FE (LAMBDA)";
267 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "MESHSET_FE",
268 print_no_field);
269
270 PrintField print_field(m_field.get_field_bit_number("FIELD"));
271 MOFEM_LOG("SYNC", Sev::inform) << "Loop elements MESHSET_FE (FIELD)";
272 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "MESHSET_FE",
273 print_field);
274
275
277
278 // check if file can be saved
279 if (debug)
280 CHKERR moab.write_file("test_out.h5m", "MOAB", "PARALLEL=WRITE_PART");
281
282 // destry dm
283 CHKERR DMDestroy(&dm);
284 }
286
287 // finish work cleaning memory, getting statistics, ect.
289
290 return 0;
291}
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
int main()
@ NOISY
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
@ NOBASE
Definition definitions.h:59
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ 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 char help[]
constexpr bool debug
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition DMMoFEM.cpp:450
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:114
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
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_MESHSET(const EntityHandle meshset, const std::string &name, const bool recursive=false)=0
add MESHSET element to finite element database given by name
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_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
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.
#define MOFEM_LOG(channel, severity)
Log.
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
const double v
phase velocity of light in medium (cm/ns)
char mesh_file_name[255]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
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
FTensor::Index< 'm', 3 > m
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
virtual MoFEMErrorCode operator()()
function is run for every finite element
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
Managing BitRefLevels.
Managing BitRefLevels.
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual moab::Interface & get_moab()=0
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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 MPI_Comm & get_comm() const =0
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:118
Deprecated interface functions.
structure for User Loop Methods on finite elements
auto getDataDofsPtr() const
auto getColDofsPtr() const
auto getRowDofsPtr() const
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 reference to pointer of interface.