v0.13.2
Loading...
Searching...
No Matches
ISManager.cpp
Go to the documentation of this file.
1/** \file ISManager.cpp
2 * \brief IS creating
3 * \ingroup mofem_is_managers
4 */
5
6
7namespace MoFEM {
8
10ISManager::query_interface(boost::typeindex::type_index type_index,
11 UnknownInterface **iface) const {
13 *iface = const_cast<ISManager *>(this);
15}
16
18 : cOre(const_cast<MoFEM::Core &>(core)), dEbug(false) {}
19
20MoFEMErrorCode ISManager::sectionCreate(const std::string problem_name,
21 PetscSection *s,
22 const RowColData row_col) const {
23 const MoFEM::Interface &m_field = cOre;
24 const Problem *problem_ptr = m_field.get_problem(problem_name);
25 auto fields_ptr = m_field.get_fields();
26 auto fe_ptr = m_field.get_finite_elements();
28 boost::shared_ptr<NumeredDofEntity_multiIndex> dofs;
29 BitFieldId fields_ids(0);
30 switch (row_col) {
31 case ROW:
32 dofs = problem_ptr->numeredRowDofsPtr;
33 for (auto fit = fe_ptr->begin(); fit != fe_ptr->end(); fit++) {
34 if ((fit->get()->getId() & problem_ptr->getBitFEId()).any()) {
35 fields_ids |= fit->get()->getBitFieldIdRow();
36 }
37 }
38 break;
39 case COL:
40 dofs = problem_ptr->numeredColDofsPtr;
41 for (auto fit = fe_ptr->begin(); fit != fe_ptr->end(); fit++) {
42 if ((fit->get()->getId() & problem_ptr->getBitFEId()).any()) {
43 fields_ids |= fit->get()->getBitFieldIdCol();
44 }
45 }
46 break;
47 default:
48 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
49 "Has to be ROW or COLUMN");
50 }
51 // get fields names on the problem
52 map<std::string, std::pair<int, int>> fields_map;
53 {
54 int field = 0;
55 for (auto fit = fields_ptr->begin(); fit != fields_ptr->end(); fit++) {
56 if ((fit->get()->getId() & fields_ids).any()) {
57 fields_map[fit->get()->getName()].first = field++;
58 fields_map[fit->get()->getName()].second = fit->get()->getNbOfCoeffs();
59 }
60 }
61 }
62 const int proc = m_field.get_comm_rank();
63 CHKERR PetscSectionCreate(PETSC_COMM_WORLD, s);
64 CHKERR PetscSectionSetNumFields(*s, fields_map.size());
65 for (auto mit = fields_map.begin(); mit != fields_map.end(); mit++) {
66 CHKERR PetscSectionSetFieldName(*s, mit->second.first, mit->first.c_str());
67 CHKERR PetscSectionSetFieldComponents(*s, mit->second.first,
68 mit->second.second);
69 }
70 // determine number of points
71 int nb_charts = 0;
72 {
73 auto dit = dofs->begin();
74 auto hi_dit = dofs->end();
75 for (; dit != hi_dit;) {
76 if (static_cast<int>(dit->get()->getPart()) == proc) {
77 const auto &ent_uid = dit->get()->getEntLocalUniqueId();
78 while (dit != hi_dit && dit->get()->getEntLocalUniqueId() == ent_uid) {
79 ++dit;
80 }
81 ++nb_charts;
82 } else {
83 ++dit;
84 }
85 }
86 }
87 // get layout, i.e. chart
88 PetscLayout layout;
89 CHKERR PetscLayoutCreate(PETSC_COMM_WORLD, &layout);
90 CHKERR PetscLayoutSetBlockSize(layout, 1);
91 CHKERR PetscLayoutSetLocalSize(layout, nb_charts);
92 CHKERR PetscLayoutSetUp(layout);
93 int rstart, rend;
94 CHKERR PetscLayoutGetRange(layout, &rstart, &rend);
95 CHKERR PetscLayoutDestroy(&layout);
96 CHKERR PetscSectionSetChart(*s, rstart, rend);
97
98 // loop of all dofs
99 {
100 auto dit = dofs->begin();
101 auto hi_dit = dofs->end();
102 int point = rstart;
103 for (; dit != hi_dit;) {
104 if (static_cast<int>(dit->get()->getPart()) == proc) {
105
106 const auto &field_name = dit->get()->getName();
107
108 int dd = 0;
109 const auto &ent_uid = dit->get()->getEntLocalUniqueId();
110 while (dit != hi_dit && dit->get()->getEntLocalUniqueId() == ent_uid) {
111 const DofIdx loc_idx = dit->get()->getPetscLocalDofIdx();
112 if (loc_idx >= 0)
113 ++dd;
114 ++dit;
115 }
116
117 if (fields_map.find(field_name) == fields_map.end()) {
118 MOFEM_LOG_C("SELF", Sev::warning, "Warning: Field %s not found",
119 dit->get()->getName().c_str());
120 } else {
121 CHKERR PetscSectionAddDof(*s, point, dd);
122 int field = fields_map.at(field_name).first;
123 CHKERR PetscSectionSetFieldDof(*s, point, field, dd);
124 }
125
126 ++point;
127
128 } else {
129 ++dit;
130 }
131 }
132 }
133 // cerr << "done " << proc << endl;
134 CHKERR PetscSectionSetUp(*s);
135 // cerr << "end " << proc << endl;
137}
138
140ISManager::sectionCreate(const std::string problem_name,
141 const RowColData row_col) const {
142
143 PetscSection s;
144 CHK_THROW_MESSAGE(sectionCreate(problem_name, &s, row_col),
145 "Section not created");
146 return SmartPetscObj<PetscSection>(s, false);
147}
148
149MoFEMErrorCode ISManager::isCreateProblem(const std::string problem_name,
150 RowColData rc, IS *is) const {
151 const MoFEM::Interface &m_field = cOre;
152 const Problem *problem_ptr;
154 CHKERR m_field.get_problem(problem_name, &problem_ptr);
155
156 const int rank = m_field.get_comm_rank();
157
158 decltype(problem_ptr->numeredRowDofsPtr) dofs_ptr;
159
160 switch (rc) {
161 case ROW:
162 dofs_ptr = problem_ptr->numeredRowDofsPtr;
163 break;
164 case COL:
165 dofs_ptr = problem_ptr->numeredColDofsPtr;
166 break;
167 default:
168 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
169 }
170
171 auto lo = dofs_ptr->get<Part_mi_tag>().lower_bound(rank);
172 auto hi = dofs_ptr->get<Part_mi_tag>().upper_bound(rank);
173 const int size = std::distance(lo, hi);
174
175 int *id;
176 CHKERR PetscMalloc(size * sizeof(int), &id);
177 int *id_it = id;
178 for (; lo != hi; ++lo, ++id_it)
179 *id_it = (*lo)->getPetscGlobalDofIdx();
180 sort(id, &id[size]);
181
182 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, size, id, PETSC_OWN_POINTER, is);
183
185}
186
187MoFEMErrorCode ISManager::isCreateProblem(const std::string problem_name,
188 RowColData rc,
189 SmartPetscObj<IS> &is) const {
191 IS raw_is;
192 CHKERR isCreateProblem(problem_name, rc, &raw_is);
193 is = SmartPetscObj<IS>(raw_is);
195}
196
197MoFEMErrorCode ISManager::isCreateProblemOrder(const std::string problem_name,
198 RowColData rc, int min_order,
199 int max_order, IS *is) const {
200 const MoFEM::Interface &m_field = cOre;
201 const Problem *problem_ptr;
203 CHKERR m_field.get_problem(problem_name, &problem_ptr);
204
205 typedef multi_index_container<
206 boost::shared_ptr<NumeredDofEntity>,
207
208 indexed_by<
209
210 sequenced<>,
211
212 ordered_non_unique<
213 tag<Order_mi_tag>,
216
217 >>
218 NumeredDofEntity_order_view_multiIndex;
219
220 const int rank = m_field.get_comm_rank();
221
222 NumeredDofEntity_order_view_multiIndex dofs_by_order;
223 auto insert_part_range = [&dofs_by_order, rank](auto &dofs) {
224 dofs_by_order.insert(dofs_by_order.end(), dofs.lower_bound(rank),
225 dofs.upper_bound(rank));
226 };
227
228 switch (rc) {
229 case ROW:
230 insert_part_range(problem_ptr->numeredRowDofsPtr->get<Part_mi_tag>());
231 break;
232 case COL:
233 insert_part_range(problem_ptr->numeredColDofsPtr->get<Part_mi_tag>());
234 break;
235 default:
236 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
237 }
238
239 auto lo = dofs_by_order.get<Order_mi_tag>().lower_bound(min_order);
240 auto hi = dofs_by_order.get<Order_mi_tag>().upper_bound(max_order);
241 const int size = std::distance(lo, hi);
242 int *id;
243 CHKERR PetscMalloc(size * sizeof(int), &id);
244 int *id_it = id;
245 for (; lo != hi; ++lo, ++id_it)
246 *id_it = (*lo)->getPetscGlobalDofIdx();
247 sort(id, &id[size]);
248
249 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, size, id, PETSC_OWN_POINTER, is);
250
252}
253
254MoFEMErrorCode ISManager::isCreateProblemOrder(const std::string problem_name,
255 RowColData rc, int min_order,
256 int max_order,
257 SmartPetscObj<IS> &is) const {
259 IS raw_is;
260 CHKERR isCreateProblemOrder(problem_name, rc, min_order, max_order, &raw_is);
261 is = SmartPetscObj<IS>(raw_is);
263}
264
266 const std::string problem_name, RowColData rc, const std::string field,
267 int min_coeff_idx, int max_coeff_idx, IS *is, Range *ents_ptr) const {
268 const MoFEM::Interface &m_field = cOre;
269 const Problem *problem_ptr;
271 CHKERR m_field.get_problem(problem_name, &problem_ptr);
272 const auto bit_number = m_field.get_field_bit_number(field);
273
274 typedef NumeredDofEntity_multiIndex::index<Unique_mi_tag>::type DofsByUId;
275 DofsByUId::iterator it, hi_it;
276 int nb_loc_dofs;
277 switch (rc) {
278 case ROW:
279 nb_loc_dofs = problem_ptr->getNbLocalDofsRow();
280 it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().lower_bound(
282 hi_it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().upper_bound(
284 break;
285 case COL:
286 nb_loc_dofs = problem_ptr->getNbLocalDofsCol();
287 it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().lower_bound(
289 hi_it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().upper_bound(
291 break;
292 default:
293 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
294 }
295
296 std::vector<int> idx_vec;
297 idx_vec.reserve(std::distance(it, hi_it));
298 for (; it != hi_it; ++it) {
299
300 auto true_if_dof_on_entity = [&]() {
301 if (ents_ptr) {
302 return ents_ptr->find((*it)->getEnt()) != ents_ptr->end();
303 } else {
304 return true;
305 }
306 };
307
308 auto check = [&]() {
309 const auto ceff_idx = (*it)->getDofCoeffIdx();
310 if (
311
312 (*it)->getPetscLocalDofIdx() >= nb_loc_dofs ||
313
314 ceff_idx < min_coeff_idx || ceff_idx > max_coeff_idx
315
316 )
317 return false;
318 else
319 return true;
320 };
321
322 if (check()) {
323 if (true_if_dof_on_entity()) {
324 idx_vec.emplace_back((*it)->getPetscGlobalDofIdx());
325 }
326 }
327 }
328
329 int *id;
330 CHKERR PetscMalloc(idx_vec.size() * sizeof(int), &id);
331 std::copy(idx_vec.begin(), idx_vec.end(), id);
332 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx_vec.size(), id,
333 PETSC_OWN_POINTER, is);
334
336}
337
339 const std::string problem_name, RowColData rc, const std::string field,
340 int min_coeff_idx, int max_coeff_idx, SmartPetscObj<IS> &smart_is,
341 Range *ents_ptr) const {
343 IS is;
344 CHKERR isCreateProblemFieldAndRank(problem_name, rc, field, min_coeff_idx,
345 max_coeff_idx, &is, ents_ptr);
346 smart_is = SmartPetscObj<IS>(is);
348}
349
351 const std::string problem_name, RowColData rc, const std::string field,
352 EntityType low_type, EntityType hi_type, int min_coeff_idx,
353 int max_coeff_idx, IS *is, Range *ents_ptr) const {
354 const MoFEM::Interface &m_field = cOre;
356 EntityHandle field_meshset = m_field.get_field_meshset(field);
357 Range ents;
358 for (; low_type <= hi_type; ++low_type)
359 CHKERR m_field.get_moab().get_entities_by_type(field_meshset, low_type,
360 ents, true);
361 if (ents_ptr)
362 ents = intersect(ents, *ents_ptr);
363 CHKERR isCreateProblemFieldAndRank(problem_name, rc, field, min_coeff_idx,
364 max_coeff_idx, is, &ents);
366}
367
369 const std::string x_problem, const std::string x_field_name,
370 RowColData x_rc, const std::string y_problem,
371 const std::string y_field_name, RowColData y_rc, std::vector<int> &idx,
372 std::vector<int> &idy) const {
373 const MoFEM::Interface &m_field = cOre;
374 const Problem *px_ptr;
375 const Problem *py_ptr;
377
378 CHKERR m_field.get_problem(x_problem, &px_ptr);
379 CHKERR m_field.get_problem(y_problem, &py_ptr);
380
381 typedef multi_index_container<
382 boost::shared_ptr<NumeredDofEntity>,
383
384 indexed_by<
385
386 sequenced<>,
387
388 ordered_non_unique<
389 tag<Composite_Ent_And_EntDofIdx_mi_tag>,
390 composite_key<
396
397 >>
398 NumeredDofEntity_view_multiIndex;
399
400 NumeredDofEntity_view_multiIndex dofs_view;
401
402 auto x_bit_number = m_field.get_field_bit_number(x_field_name);
403
404 switch (x_rc) {
405 case ROW:
406 dofs_view.insert(
407 dofs_view.end(),
408 px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().lower_bound(
409 FieldEntity::getLoBitNumberUId(x_bit_number)),
410 px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().upper_bound(
411 FieldEntity::getHiBitNumberUId(x_bit_number)));
412 break;
413 case COL:
414 dofs_view.insert(
415 dofs_view.end(),
416 px_ptr->numeredColDofsPtr->get<Unique_mi_tag>().lower_bound(
417 FieldEntity::getLoBitNumberUId(x_bit_number)),
418 px_ptr->numeredColDofsPtr->get<Unique_mi_tag>().upper_bound(
419 FieldEntity::getHiBitNumberUId(x_bit_number)));
420 break;
421 default:
422 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
423 "only makes sense for ROWS and COLS");
424 }
425
426 decltype(py_ptr->numeredRowDofsPtr) dofs_ptr;
427 switch (y_rc) {
428 case ROW:
429 dofs_ptr = py_ptr->numeredRowDofsPtr;
430 break;
431 case COL:
432 dofs_ptr = py_ptr->numeredColDofsPtr;
433 break;
434 default:
435 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
436 "only makes sense for ROWS and COLS");
437 }
438
439 std::map<int, int> global_dofs_map;
440 const auto y_bit_number = m_field.get_field_bit_number(y_field_name);
441 auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(
442 FieldEntity::getLoBitNumberUId(y_bit_number));
443 auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(
444 FieldEntity::getHiBitNumberUId(y_bit_number));
445 const auto rank = m_field.get_comm_rank();
446 for (; dit != hi_dit; ++dit) {
447 if ((*dit)->getPart() == rank) {
448 auto x_dit = dofs_view.get<Composite_Ent_And_EntDofIdx_mi_tag>().find(
449 boost::make_tuple((*dit)->getEnt(), (*dit)->getEntDofIdx()));
450 if (x_dit != dofs_view.get<Composite_Ent_And_EntDofIdx_mi_tag>().end()) {
451 global_dofs_map[(*x_dit)->getPetscGlobalDofIdx()] =
452 (*dit)->getPetscGlobalDofIdx();
453 }
454 }
455 }
456
457 idx.resize(global_dofs_map.size());
458 idy.resize(global_dofs_map.size());
459 {
460 auto ix = idx.begin();
461 auto iy = idy.begin();
462 for (auto mit = global_dofs_map.begin(); mit != global_dofs_map.end();
463 mit++, ix++, iy++) {
464 *ix = mit->first;
465 *iy = mit->second;
466 }
467 }
469}
470
472 const std::string x_problem, const std::string x_field_name,
473 RowColData x_rc, const std::string y_problem,
474 const std::string y_field_name, RowColData y_rc, IS *ix, IS *iy) const {
476 std::vector<int> idx(0), idy(0);
478 x_problem, x_field_name, x_rc, y_problem, y_field_name, y_rc, idx, idy);
479 if (ix != PETSC_NULL) {
480 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
481 PETSC_COPY_VALUES, ix);
482 }
483 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
484 PETSC_COPY_VALUES, iy);
485 if (dEbug) {
486 ISView(*ix, PETSC_VIEWER_STDOUT_WORLD);
487 ISView(*iy, PETSC_VIEWER_STDOUT_WORLD);
488 }
490}
491
493 const std::string x_problem, RowColData x_rc, const std::string y_problem,
494 RowColData y_rc, std::vector<int> &idx, std::vector<int> &idy) const {
495 const MoFEM::Interface &m_field = cOre;
496 const Problem *px_ptr;
497 const Problem *py_ptr;
499 CHKERR m_field.get_problem(x_problem, &px_ptr);
500 CHKERR m_field.get_problem(y_problem, &py_ptr);
501 NumeredDofEntityByLocalIdx::iterator y_dit, hi_y_dit;
502 switch (y_rc) {
503 case ROW:
504 y_dit =
505 py_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(0);
506 hi_y_dit =
507 py_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(
508 py_ptr->getNbLocalDofsRow()); // should be lower
509 break;
510 case COL:
511 y_dit =
512 py_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(0);
513 hi_y_dit =
514 py_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(
515 py_ptr->getNbLocalDofsCol()); // should be lower
516 break;
517 default:
518 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
519 }
520 const NumeredDofEntityByUId *x_numered_dofs_by_uid;
521 switch (x_rc) {
522 case ROW:
523 x_numered_dofs_by_uid = &(px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
524 break;
525 case COL:
526 x_numered_dofs_by_uid = &(px_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
527 break;
528 default:
529 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
530 }
531 for (; y_dit != hi_y_dit; y_dit++) {
532 if ((*y_dit)->getPart() != (unsigned int)m_field.get_comm_rank()) {
533 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
534 }
535 NumeredDofEntityByUId::iterator x_dit;
536 x_dit = x_numered_dofs_by_uid->find((*y_dit)->getLocalUniqueId());
537 if (x_dit == x_numered_dofs_by_uid->end())
538 continue;
539 idx.push_back((*x_dit)->getPetscGlobalDofIdx());
540 idy.push_back((*y_dit)->getPetscGlobalDofIdx());
541 }
543}
544
546 const std::string x_problem, RowColData x_rc, const std::string y_problem,
547 RowColData y_rc, IS *ix, IS *iy) const {
549 std::vector<int> idx(0), idy(0);
550 CHKERR isCreateFromProblemToOtherProblem(x_problem, x_rc, y_problem, y_rc,
551 idx, idy);
552 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
553 PETSC_COPY_VALUES, ix);
554 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
555 PETSC_COPY_VALUES, iy);
556 if (dEbug) {
557 ISView(*ix, PETSC_VIEWER_STDOUT_WORLD);
558 ISView(*iy, PETSC_VIEWER_STDOUT_WORLD);
559 }
561}
562} // namespace MoFEM
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
RowColData
RowColData.
Definition: definitions.h:123
@ COL
Definition: definitions.h:123
@ ROW
Definition: definitions.h:123
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:595
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
#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
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
NumeredDofEntity_multiIndex::index< Unique_mi_tag >::type NumeredDofEntityByUId
Numbered DoF multi-index by UId.
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
MoFEMErrorCode isCreateFromProblemFieldToOtherProblemField(const std::string x_problem, const std::string x_field_name, RowColData x_rc, const std::string y_problem, const std::string y_field_name, RowColData y_rc, std::vector< int > &idx, std::vector< int > &idy) const
create IS for give two problems and field
Definition: ISManager.cpp:368
MoFEMErrorCode sectionCreate(const std::string problem_name, PetscSection *s, const RowColData row_col=COL) const
Create global selection.
Definition: ISManager.cpp:20
MoFEMErrorCode isCreateProblemFieldAndEntityType(const std::string problem_name, RowColData rc, const std::string field, EntityType low_type, EntityType hi_type, int min_coeff_idx, int max_coeff_idx, IS *is, Range *ents=nullptr) const
create IS for given problem, field and type range (collective)
Definition: ISManager.cpp:350
MoFEMErrorCode isCreateFromProblemToOtherProblem(const std::string x_problem, RowColData x_rc, const std::string y_problem, RowColData y_rc, std::vector< int > &idx, std::vector< int > &idy) const
Create is from one problem to other problem.
Definition: ISManager.cpp:492
MoFEMErrorCode isCreateProblemOrder(const std::string problem_name, RowColData rc, int min_order, int max_order, IS *is) const
create IS for given order range (collective)
Definition: ISManager.cpp:197
MoFEMErrorCode isCreateProblemFieldAndRank(const std::string problem_name, RowColData rc, const std::string field, int min_coeff_idx, int max_coeff_idx, IS *is, Range *ents=nullptr) const
create IS for given problem, field and rank range (collective)
Definition: ISManager.cpp:265
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
int DofIdx
Index of DOF.
Definition: Types.hpp:18
int ApproximationOrder
Approximation on the entity.
Definition: Types.hpp:26
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:42
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
constexpr auto field_name
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 int get_comm_rank() const =0
Core (interface) class.
Definition: Core.hpp:82
Deprecated interface functions.
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
MoFEMErrorCode isCreateProblem(const std::string problem_name, RowColData rc, IS *is) const
Create IS for problem.
Definition: ISManager.cpp:149
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: ISManager.cpp:10
const MoFEM::Interface & cOre
Definition: ISManager.hpp:28
ISManager(const MoFEM::Core &core)
Definition: ISManager.cpp:17
keeps information about indexed dofs for the problem
interface_DofEntity< DofEntity > interface_type_DofEntity
MultiIndex Tag for field order.
keeps basic data about problem
DofIdx getNbLocalDofsRow() const
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
BitFEId getBitFEId() const
DofIdx getNbLocalDofsCol() const
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem
intrusive_ptr for managing petsc objects
base class for all interface classes