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::isCreateProblemOrder(const std::string problem_name,
150 RowColData rc, int min_order,
151 int max_order, IS *is) const {
152 const MoFEM::Interface &m_field = cOre;
153 const Problem *problem_ptr;
155 CHKERR m_field.get_problem(problem_name, &problem_ptr);
156
157 typedef multi_index_container<
158 boost::shared_ptr<NumeredDofEntity>,
159
160 indexed_by<
161
162 sequenced<>,
163
164 ordered_non_unique<
165 tag<Order_mi_tag>,
168
169 >>
170 NumeredDofEntity_order_view_multiIndex;
171
172 const int rank = m_field.get_comm_rank();
173
174 NumeredDofEntity_order_view_multiIndex dofs_by_order;
175 auto insert_part_range = [&dofs_by_order, rank](auto &dofs) {
176 dofs_by_order.insert(dofs_by_order.end(), dofs.lower_bound(rank),
177 dofs.upper_bound(rank));
178 };
179
180 switch (rc) {
181 case ROW:
182 insert_part_range(problem_ptr->numeredRowDofsPtr->get<Part_mi_tag>());
183 break;
184 case COL:
185 insert_part_range(problem_ptr->numeredColDofsPtr->get<Part_mi_tag>());
186 break;
187 default:
188 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
189 }
190
191 auto lo = dofs_by_order.get<Order_mi_tag>().lower_bound(min_order);
192 auto hi = dofs_by_order.get<Order_mi_tag>().upper_bound(max_order);
193 const int size = std::distance(lo, hi);
194 int *id;
195 CHKERR PetscMalloc(size * sizeof(int), &id);
196 int *id_it = id;
197 for (; lo != hi; ++lo, ++id_it)
198 *id_it = (*lo)->getPetscGlobalDofIdx();
199 sort(id, &id[size]);
200
201 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, size, id, PETSC_OWN_POINTER, is);
202
204}
205
206MoFEMErrorCode ISManager::isCreateProblemOrder(const std::string problem_name,
207 RowColData rc, int min_order,
208 int max_order,
209 SmartPetscObj<IS> &is) const {
211 IS raw_is;
212 CHKERR isCreateProblemOrder(problem_name, rc, min_order, max_order, &raw_is);
213 is = SmartPetscObj<IS>(raw_is);
215}
216
218 const std::string problem_name, RowColData rc, const std::string field,
219 int min_coeff_idx, int max_coeff_idx, IS *is, Range *ents_ptr) const {
220 const MoFEM::Interface &m_field = cOre;
221 const Problem *problem_ptr;
223 CHKERR m_field.get_problem(problem_name, &problem_ptr);
224 const auto bit_number = m_field.get_field_bit_number(field);
225
226 typedef NumeredDofEntity_multiIndex::index<Unique_mi_tag>::type DofsByUId;
227 DofsByUId::iterator it, hi_it;
228 int nb_loc_dofs;
229 switch (rc) {
230 case ROW:
231 nb_loc_dofs = problem_ptr->getNbLocalDofsRow();
232 it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().lower_bound(
234 hi_it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().upper_bound(
236 break;
237 case COL:
238 nb_loc_dofs = problem_ptr->getNbLocalDofsCol();
239 it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().lower_bound(
241 hi_it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().upper_bound(
243 break;
244 default:
245 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
246 }
247
248 std::vector<int> idx_vec;
249 idx_vec.reserve(std::distance(it, hi_it));
250 for (; it != hi_it; ++it) {
251
252 auto true_if_dof_on_entity = [&]() {
253 if (ents_ptr) {
254 return ents_ptr->find((*it)->getEnt()) != ents_ptr->end();
255 } else {
256 return true;
257 }
258 };
259
260 auto check = [&]() {
261 const auto ceff_idx = (*it)->getDofCoeffIdx();
262 if (
263
264 (*it)->getPetscLocalDofIdx() >= nb_loc_dofs ||
265
266 ceff_idx < min_coeff_idx || ceff_idx > max_coeff_idx
267
268 )
269 return false;
270 else
271 return true;
272 };
273
274 if (check()) {
275 if (true_if_dof_on_entity()) {
276 idx_vec.emplace_back((*it)->getPetscGlobalDofIdx());
277 }
278 }
279 }
280
281 int *id;
282 CHKERR PetscMalloc(idx_vec.size() * sizeof(int), &id);
283 std::copy(idx_vec.begin(), idx_vec.end(), id);
284 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx_vec.size(), id,
285 PETSC_OWN_POINTER, is);
286
288}
289
291 const std::string problem_name, RowColData rc, const std::string field,
292 int min_coeff_idx, int max_coeff_idx, SmartPetscObj<IS> &smart_is,
293 Range *ents_ptr) const {
295 IS is;
296 CHKERR isCreateProblemFieldAndRank(problem_name, rc, field, min_coeff_idx,
297 max_coeff_idx, &is, ents_ptr);
298 smart_is = SmartPetscObj<IS>(is);
300}
301
303 const std::string problem_name, RowColData rc, const std::string field,
304 EntityType low_type, EntityType hi_type, int min_coeff_idx,
305 int max_coeff_idx, IS *is, Range *ents_ptr) const {
306 const MoFEM::Interface &m_field = cOre;
308 EntityHandle field_meshset = m_field.get_field_meshset(field);
309 Range ents;
310 for (; low_type <= hi_type; ++low_type)
311 CHKERR m_field.get_moab().get_entities_by_type(field_meshset, low_type,
312 ents, true);
313 if (ents_ptr)
314 ents = intersect(ents, *ents_ptr);
315 CHKERR isCreateProblemFieldAndRank(problem_name, rc, field, min_coeff_idx,
316 max_coeff_idx, is, &ents);
318}
319
321 const std::string x_problem, const std::string x_field_name,
322 RowColData x_rc, const std::string y_problem,
323 const std::string y_field_name, RowColData y_rc, std::vector<int> &idx,
324 std::vector<int> &idy) const {
325 const MoFEM::Interface &m_field = cOre;
326 const Problem *px_ptr;
327 const Problem *py_ptr;
329
330 CHKERR m_field.get_problem(x_problem, &px_ptr);
331 CHKERR m_field.get_problem(y_problem, &py_ptr);
332
333 typedef multi_index_container<
334 boost::shared_ptr<NumeredDofEntity>,
335
336 indexed_by<
337
338 sequenced<>,
339
340 ordered_non_unique<
341 tag<Composite_Ent_And_EntDofIdx_mi_tag>,
342 composite_key<
348
349 >>
350 NumeredDofEntity_view_multiIndex;
351
352 NumeredDofEntity_view_multiIndex dofs_view;
353
354 auto x_bit_number = m_field.get_field_bit_number(x_field_name);
355
356 switch (x_rc) {
357 case ROW:
358 dofs_view.insert(
359 dofs_view.end(),
360 px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().lower_bound(
361 FieldEntity::getLoBitNumberUId(x_bit_number)),
362 px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().upper_bound(
363 FieldEntity::getHiBitNumberUId(x_bit_number)));
364 break;
365 case COL:
366 dofs_view.insert(
367 dofs_view.end(),
368 px_ptr->numeredColDofsPtr->get<Unique_mi_tag>().lower_bound(
369 FieldEntity::getLoBitNumberUId(x_bit_number)),
370 px_ptr->numeredColDofsPtr->get<Unique_mi_tag>().upper_bound(
371 FieldEntity::getHiBitNumberUId(x_bit_number)));
372 break;
373 default:
374 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
375 "only makes sense for ROWS and COLS");
376 }
377
378 decltype(py_ptr->numeredRowDofsPtr) dofs_ptr;
379 switch (y_rc) {
380 case ROW:
381 dofs_ptr = py_ptr->numeredRowDofsPtr;
382 break;
383 case COL:
384 dofs_ptr = py_ptr->numeredColDofsPtr;
385 break;
386 default:
387 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
388 "only makes sense for ROWS and COLS");
389 }
390
391 std::map<int, int> global_dofs_map;
392 const auto y_bit_number = m_field.get_field_bit_number(y_field_name);
393 auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(
394 FieldEntity::getLoBitNumberUId(y_bit_number));
395 auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(
396 FieldEntity::getHiBitNumberUId(y_bit_number));
397 const auto rank = m_field.get_comm_rank();
398 for (; dit != hi_dit; ++dit) {
399 if ((*dit)->getPart() == rank) {
400 auto x_dit = dofs_view.get<Composite_Ent_And_EntDofIdx_mi_tag>().find(
401 boost::make_tuple((*dit)->getEnt(), (*dit)->getEntDofIdx()));
402 if (x_dit != dofs_view.get<Composite_Ent_And_EntDofIdx_mi_tag>().end()) {
403 global_dofs_map[(*x_dit)->getPetscGlobalDofIdx()] =
404 (*dit)->getPetscGlobalDofIdx();
405 }
406 }
407 }
408
409 idx.resize(global_dofs_map.size());
410 idy.resize(global_dofs_map.size());
411 {
412 auto ix = idx.begin();
413 auto iy = idy.begin();
414 for (auto mit = global_dofs_map.begin(); mit != global_dofs_map.end();
415 mit++, ix++, iy++) {
416 *ix = mit->first;
417 *iy = mit->second;
418 }
419 }
421}
422
424 const std::string x_problem, const std::string x_field_name,
425 RowColData x_rc, const std::string y_problem,
426 const std::string y_field_name, RowColData y_rc, IS *ix, IS *iy) const {
428 std::vector<int> idx(0), idy(0);
430 x_problem, x_field_name, x_rc, y_problem, y_field_name, y_rc, idx, idy);
431 if (ix != PETSC_NULL) {
432 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
433 PETSC_COPY_VALUES, ix);
434 }
435 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
436 PETSC_COPY_VALUES, iy);
437 if (dEbug) {
438 ISView(*ix, PETSC_VIEWER_STDOUT_WORLD);
439 ISView(*iy, PETSC_VIEWER_STDOUT_WORLD);
440 }
442}
443
445 const std::string x_problem, RowColData x_rc, const std::string y_problem,
446 RowColData y_rc, std::vector<int> &idx, std::vector<int> &idy) const {
447 const MoFEM::Interface &m_field = cOre;
448 const Problem *px_ptr;
449 const Problem *py_ptr;
451 CHKERR m_field.get_problem(x_problem, &px_ptr);
452 CHKERR m_field.get_problem(y_problem, &py_ptr);
453 NumeredDofEntityByLocalIdx::iterator y_dit, hi_y_dit;
454 switch (y_rc) {
455 case ROW:
456 y_dit =
457 py_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(0);
458 hi_y_dit =
459 py_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(
460 py_ptr->getNbLocalDofsRow()); // should be lower
461 break;
462 case COL:
463 y_dit =
464 py_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(0);
465 hi_y_dit =
466 py_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(
467 py_ptr->getNbLocalDofsCol()); // should be lower
468 break;
469 default:
470 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
471 }
472 const NumeredDofEntityByUId *x_numered_dofs_by_uid;
473 switch (x_rc) {
474 case ROW:
475 x_numered_dofs_by_uid = &(px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
476 break;
477 case COL:
478 x_numered_dofs_by_uid = &(px_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
479 break;
480 default:
481 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
482 }
483 for (; y_dit != hi_y_dit; y_dit++) {
484 if ((*y_dit)->getPart() != (unsigned int)m_field.get_comm_rank()) {
485 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
486 }
487 NumeredDofEntityByUId::iterator x_dit;
488 x_dit = x_numered_dofs_by_uid->find((*y_dit)->getLocalUniqueId());
489 if (x_dit == x_numered_dofs_by_uid->end())
490 continue;
491 idx.push_back((*x_dit)->getPetscGlobalDofIdx());
492 idy.push_back((*y_dit)->getPetscGlobalDofIdx());
493 }
495}
496
498 const std::string x_problem, RowColData x_rc, const std::string y_problem,
499 RowColData y_rc, IS *ix, IS *iy) const {
501 std::vector<int> idx(0), idy(0);
502 CHKERR isCreateFromProblemToOtherProblem(x_problem, x_rc, y_problem, y_rc,
503 idx, idy);
504 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
505 PETSC_COPY_VALUES, ix);
506 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
507 PETSC_COPY_VALUES, iy);
508 if (dEbug) {
509 ISView(*ix, PETSC_VIEWER_STDOUT_WORLD);
510 ISView(*iy, PETSC_VIEWER_STDOUT_WORLD);
511 }
513}
514} // namespace MoFEM
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:304
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:320
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:302
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:444
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:149
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:217
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 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