v0.14.0
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 
7 namespace MoFEM {
8 
10 ISManager::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 
20 MoFEMErrorCode 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 
140 ISManager::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 
149 MoFEMErrorCode 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 
187 MoFEMErrorCode 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 
197 SmartPetscObj<IS> ISManager::isCreateProblem(const std::string problem_name,
198  RowColData rc) const {
200  CHK_THROW_MESSAGE(isCreateProblem(problem_name, rc, is), "IS not created");
201  return is;
202 }
203 
204 MoFEMErrorCode ISManager::isCreateProblemOrder(const std::string problem_name,
205  RowColData rc, int min_order,
206  int max_order, IS *is) const {
207  const MoFEM::Interface &m_field = cOre;
208  const Problem *problem_ptr;
210  CHKERR m_field.get_problem(problem_name, &problem_ptr);
211 
212  typedef multi_index_container<
213  boost::shared_ptr<NumeredDofEntity>,
214 
215  indexed_by<
216 
217  sequenced<>,
218 
219  ordered_non_unique<
220  tag<Order_mi_tag>,
223 
224  >>
225  NumeredDofEntity_order_view_multiIndex;
226 
227  const int rank = m_field.get_comm_rank();
228 
229  NumeredDofEntity_order_view_multiIndex dofs_by_order;
230  auto insert_part_range = [&dofs_by_order, rank](auto &dofs) {
231  dofs_by_order.insert(dofs_by_order.end(), dofs.lower_bound(rank),
232  dofs.upper_bound(rank));
233  };
234 
235  switch (rc) {
236  case ROW:
237  insert_part_range(problem_ptr->numeredRowDofsPtr->get<Part_mi_tag>());
238  break;
239  case COL:
240  insert_part_range(problem_ptr->numeredColDofsPtr->get<Part_mi_tag>());
241  break;
242  default:
243  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
244  }
245 
246  auto lo = dofs_by_order.get<Order_mi_tag>().lower_bound(min_order);
247  auto hi = dofs_by_order.get<Order_mi_tag>().upper_bound(max_order);
248  const int size = std::distance(lo, hi);
249  int *id;
250  CHKERR PetscMalloc(size * sizeof(int), &id);
251  int *id_it = id;
252  for (; lo != hi; ++lo, ++id_it)
253  *id_it = (*lo)->getPetscGlobalDofIdx();
254  sort(id, &id[size]);
255 
256  CHKERR ISCreateGeneral(PETSC_COMM_WORLD, size, id, PETSC_OWN_POINTER, is);
257 
259 }
260 
261 MoFEMErrorCode ISManager::isCreateProblemOrder(const std::string problem_name,
262  RowColData rc, int min_order,
263  int max_order,
264  SmartPetscObj<IS> &is) const {
266  IS raw_is;
267  CHKERR isCreateProblemOrder(problem_name, rc, min_order, max_order, &raw_is);
268  is = SmartPetscObj<IS>(raw_is);
270 }
271 
273  const std::string problem_name, RowColData rc, const std::string field,
274  int min_coeff_idx, int max_coeff_idx, IS *is, Range *ents_ptr) const {
275  const MoFEM::Interface &m_field = cOre;
276  const Problem *problem_ptr;
278  CHKERR m_field.get_problem(problem_name, &problem_ptr);
279  const auto bit_number = m_field.get_field_bit_number(field);
280 
282  DofsByUId::iterator it, hi_it;
283  int nb_loc_dofs;
284  switch (rc) {
285  case ROW:
286  nb_loc_dofs = problem_ptr->getNbLocalDofsRow();
287  it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().lower_bound(
288  FieldEntity::getLoBitNumberUId(bit_number));
289  hi_it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().upper_bound(
290  FieldEntity::getHiBitNumberUId(bit_number));
291  break;
292  case COL:
293  nb_loc_dofs = problem_ptr->getNbLocalDofsCol();
294  it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().lower_bound(
295  FieldEntity::getLoBitNumberUId(bit_number));
296  hi_it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().upper_bound(
297  FieldEntity::getHiBitNumberUId(bit_number));
298  break;
299  default:
300  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
301  }
302 
303  std::vector<int> idx_vec;
304  idx_vec.reserve(std::distance(it, hi_it));
305  for (; it != hi_it; ++it) {
306 
307  auto true_if_dof_on_entity = [&]() {
308  if (ents_ptr) {
309  return ents_ptr->find((*it)->getEnt()) != ents_ptr->end();
310  } else {
311  return true;
312  }
313  };
314 
315  auto check = [&]() {
316  const auto coeff_idx = (*it)->getDofCoeffIdx();
317  if (
318 
319  (*it)->getPetscLocalDofIdx() >= nb_loc_dofs ||
320 
321  coeff_idx < min_coeff_idx || coeff_idx > max_coeff_idx
322 
323  )
324  return false;
325  else
326  return true;
327  };
328 
329  if (check()) {
330  if (true_if_dof_on_entity()) {
331  idx_vec.emplace_back((*it)->getPetscGlobalDofIdx());
332  }
333  }
334  }
335 
336  int *id;
337  CHKERR PetscMalloc(idx_vec.size() * sizeof(int), &id);
338  std::copy(idx_vec.begin(), idx_vec.end(), id);
339  CHKERR ISCreateGeneral(m_field.get_comm(), idx_vec.size(), id,
340  PETSC_OWN_POINTER, is);
341 
343 }
344 
346  const std::string problem_name, RowColData rc, const std::string field,
347  int min_coeff_idx, int max_coeff_idx, SmartPetscObj<IS> &smart_is,
348  Range *ents_ptr) const {
350  IS is;
351  CHKERR isCreateProblemFieldAndRank(problem_name, rc, field, min_coeff_idx,
352  max_coeff_idx, &is, ents_ptr);
353  smart_is = SmartPetscObj<IS>(is);
355 }
356 
358  const std::vector<boost::weak_ptr<NumeredDofEntity>> &dofs_vec,
359  SmartPetscObj<IS> &smart_is, MPI_Comm comm) const {
361 
362  std::vector<int> idx_vec;
363  idx_vec.reserve(dofs_vec.size());
364  for (auto &dof : dofs_vec) {
365  if (auto d = dof.lock()) {
366  idx_vec.emplace_back(d->getPetscGlobalDofIdx());
367  }
368  }
369 
370  IS is_raw;
371  int *id;
372  CHKERR PetscMalloc(idx_vec.size() * sizeof(int), &id);
373  std::copy(idx_vec.begin(), idx_vec.end(), id);
374  CHKERR ISCreateGeneral(comm, idx_vec.size(), id, PETSC_OWN_POINTER, &is_raw);
375 
376  smart_is = SmartPetscObj<IS>(is_raw);
377 
379 }
380 
382  const std::vector<boost::weak_ptr<NumeredDofEntity>> &dofs_vec,
383  SmartPetscObj<IS> &smart_is) const {
385 
386  std::vector<int> idx_vec;
387  idx_vec.reserve(dofs_vec.size());
388  for (auto &dof : dofs_vec) {
389  if (auto d = dof.lock()) {
390  idx_vec.emplace_back(d->getPetscLocalDofIdx());
391  }
392  }
393 
394  IS is_raw;
395  int *id;
396  CHKERR PetscMalloc(idx_vec.size() * sizeof(int), &id);
397  std::copy(idx_vec.begin(), idx_vec.end(), id);
398  CHKERR ISCreateGeneral(PETSC_COMM_SELF, idx_vec.size(), id, PETSC_OWN_POINTER,
399  &is_raw);
400 
401  smart_is = SmartPetscObj<IS>(is_raw);
402 
404 }
405 
407  const std::string problem_name, RowColData rc, const std::string field,
408  int min_coeff_idx, int max_coeff_idx, IS *is, Range *ents_ptr) const {
409  const MoFEM::Interface &m_field = cOre;
410  const Problem *problem_ptr;
412  CHKERR m_field.get_problem(problem_name, &problem_ptr);
413  const auto bit_number = m_field.get_field_bit_number(field);
414 
415  auto get_low_hi_uid = [&]() {
416  return std::make_pair(FieldEntity::getLoBitNumberUId(bit_number),
417  FieldEntity::getHiBitNumberUId(bit_number));
418  };
419 
420  auto get_low_hi_uid_by_entities = [&](auto f, auto s) {
421  return std::make_pair(DofEntity::getLoFieldEntityUId(bit_number, f),
422  DofEntity::getHiFieldEntityUId(bit_number, s));
423  };
424 
425  auto get_low_hi = [&](auto lo_uid, auto hi_uid) {
427  DofsByUId::iterator it, hi_it;
428  switch (rc) {
429  case ROW:
430  it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().lower_bound(
431  lo_uid);
432  hi_it = problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().upper_bound(
433  hi_uid);
434  break;
435  case COL:
436  it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().lower_bound(
437  lo_uid);
438  hi_it = problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>().upper_bound(
439  hi_uid);
440  break;
441  default:
442  THROW_MESSAGE("not implemented");
443  }
444  return std::make_pair(it, hi_it);
445  };
446 
447  auto check = [&](auto it) {
448  const auto coeff_idx = (*it)->getDofCoeffIdx();
449  if (
450 
451  coeff_idx < min_coeff_idx || coeff_idx > max_coeff_idx
452 
453  )
454  return false;
455  else
456  return true;
457  };
458 
459  auto emplace_indices = [&](auto it, auto hi_it, auto &idx_vec) {
460  for (; it != hi_it; ++it) {
461  if (check(it))
462  idx_vec.emplace_back((*it)->getPetscLocalDofIdx());
463  }
464  };
465 
466  auto [lo_uid, hi_uid] = get_low_hi_uid();
467  auto [lo, hi] = get_low_hi(lo_uid, hi_uid);
468  std::vector<int> idx_vec;
469  idx_vec.reserve(std::distance(lo, hi));
470 
471  if (ents_ptr) {
472  for (auto pit = ents_ptr->const_pair_begin();
473  pit != ents_ptr->const_pair_end(); ++pit) {
474  auto [lo_uid, hi_uid] =
475  get_low_hi_uid_by_entities(pit->first, pit->second);
476  auto [lo, hi] = get_low_hi(lo_uid, hi_uid);
477  emplace_indices(lo, hi, idx_vec);
478  }
479  } else {
480  emplace_indices(lo, hi, idx_vec);
481  }
482 
483  int *id;
484  CHKERR PetscMalloc(idx_vec.size() * sizeof(int), &id);
485  std::copy(idx_vec.begin(), idx_vec.end(), id);
486  CHKERR ISCreateGeneral(PETSC_COMM_SELF, idx_vec.size(), id, PETSC_OWN_POINTER,
487  is);
488 
490 }
491 
493  const std::string problem_name, RowColData rc, const std::string field,
494  int min_coeff_idx, int max_coeff_idx, SmartPetscObj<IS> &smart_is,
495  Range *ents_ptr) const {
497  IS is;
498  CHKERR isCreateProblemFieldAndRankLocal(problem_name, rc, field,
499  min_coeff_idx, max_coeff_idx, &is,
500  ents_ptr);
501  smart_is = SmartPetscObj<IS>(is);
503 }
504 
506  const std::string problem_name, RowColData rc, const std::string field,
507  EntityType low_type, EntityType hi_type, int min_coeff_idx,
508  int max_coeff_idx, IS *is, Range *ents_ptr) const {
509  const MoFEM::Interface &m_field = cOre;
511  EntityHandle field_meshset = m_field.get_field_meshset(field);
512  Range ents;
513  for (; low_type <= hi_type; ++low_type)
514  CHKERR m_field.get_moab().get_entities_by_type(field_meshset, low_type,
515  ents, true);
516  if (ents_ptr)
517  ents = intersect(ents, *ents_ptr);
518  CHKERR isCreateProblemFieldAndRank(problem_name, rc, field, min_coeff_idx,
519  max_coeff_idx, is, &ents);
521 }
522 
524  const std::string x_problem, const std::string x_field_name,
525  RowColData x_rc, const std::string y_problem,
526  const std::string y_field_name, RowColData y_rc, std::vector<int> &idx,
527  std::vector<int> &idy) const {
528  const MoFEM::Interface &m_field = cOre;
529  const Problem *px_ptr;
530  const Problem *py_ptr;
532 
533  CHKERR m_field.get_problem(x_problem, &px_ptr);
534  CHKERR m_field.get_problem(y_problem, &py_ptr);
535 
536  typedef multi_index_container<
537  boost::shared_ptr<NumeredDofEntity>,
538 
539  indexed_by<
540 
541  sequenced<>,
542 
543  ordered_non_unique<
544  tag<Composite_Ent_And_EntDofIdx_mi_tag>,
545  composite_key<
551 
552  >>
553  NumeredDofEntity_view_multiIndex;
554 
555  NumeredDofEntity_view_multiIndex dofs_view;
556 
557  auto x_bit_number = m_field.get_field_bit_number(x_field_name);
558 
559  switch (x_rc) {
560  case ROW:
561  dofs_view.insert(
562  dofs_view.end(),
563  px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().lower_bound(
564  FieldEntity::getLoBitNumberUId(x_bit_number)),
565  px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>().upper_bound(
566  FieldEntity::getHiBitNumberUId(x_bit_number)));
567  break;
568  case COL:
569  dofs_view.insert(
570  dofs_view.end(),
571  px_ptr->numeredColDofsPtr->get<Unique_mi_tag>().lower_bound(
572  FieldEntity::getLoBitNumberUId(x_bit_number)),
573  px_ptr->numeredColDofsPtr->get<Unique_mi_tag>().upper_bound(
574  FieldEntity::getHiBitNumberUId(x_bit_number)));
575  break;
576  default:
577  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
578  "only makes sense for ROWS and COLS");
579  }
580 
581  decltype(py_ptr->numeredRowDofsPtr) dofs_ptr;
582  switch (y_rc) {
583  case ROW:
584  dofs_ptr = py_ptr->numeredRowDofsPtr;
585  break;
586  case COL:
587  dofs_ptr = py_ptr->numeredColDofsPtr;
588  break;
589  default:
590  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
591  "only makes sense for ROWS and COLS");
592  }
593 
594  std::map<int, int> global_dofs_map;
595  const auto y_bit_number = m_field.get_field_bit_number(y_field_name);
596  auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(
597  FieldEntity::getLoBitNumberUId(y_bit_number));
598  auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(
599  FieldEntity::getHiBitNumberUId(y_bit_number));
600  const auto rank = m_field.get_comm_rank();
601  for (; dit != hi_dit; ++dit) {
602  if ((*dit)->getPart() == rank) {
603  auto x_dit = dofs_view.get<Composite_Ent_And_EntDofIdx_mi_tag>().find(
604  boost::make_tuple((*dit)->getEnt(), (*dit)->getEntDofIdx()));
605  if (x_dit != dofs_view.get<Composite_Ent_And_EntDofIdx_mi_tag>().end()) {
606  global_dofs_map[(*x_dit)->getPetscGlobalDofIdx()] =
607  (*dit)->getPetscGlobalDofIdx();
608  }
609  }
610  }
611 
612  idx.resize(global_dofs_map.size());
613  idy.resize(global_dofs_map.size());
614  {
615  auto ix = idx.begin();
616  auto iy = idy.begin();
617  for (auto mit = global_dofs_map.begin(); mit != global_dofs_map.end();
618  mit++, ix++, iy++) {
619  *ix = mit->first;
620  *iy = mit->second;
621  }
622  }
624 }
625 
627  const std::string x_problem, const std::string x_field_name,
628  RowColData x_rc, const std::string y_problem,
629  const std::string y_field_name, RowColData y_rc, IS *ix, IS *iy) const {
631  const MoFEM::Interface &m_field = cOre;
632  std::vector<int> idx(0), idy(0);
634  x_problem, x_field_name, x_rc, y_problem, y_field_name, y_rc, idx, idy);
635  if (ix != PETSC_NULL) {
636  CHKERR ISCreateGeneral(m_field.get_comm(), idx.size(), &idx[0],
637  PETSC_COPY_VALUES, ix);
638  }
639  CHKERR ISCreateGeneral(m_field.get_comm(), idy.size(), &idy[0],
640  PETSC_COPY_VALUES, iy);
641  if (dEbug) {
642  ISView(*ix, PETSC_VIEWER_STDOUT_WORLD);
643  ISView(*iy, PETSC_VIEWER_STDOUT_WORLD);
644  }
646 }
647 
649  const std::string x_problem, RowColData x_rc, const std::string y_problem,
650  RowColData y_rc, std::vector<int> &idx, std::vector<int> &idy) const {
651  const MoFEM::Interface &m_field = cOre;
652  const Problem *px_ptr;
653  const Problem *py_ptr;
655  CHKERR m_field.get_problem(x_problem, &px_ptr);
656  CHKERR m_field.get_problem(y_problem, &py_ptr);
657  NumeredDofEntityByLocalIdx::iterator y_dit, hi_y_dit;
658  switch (y_rc) {
659  case ROW:
660  y_dit =
661  py_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(0);
662  hi_y_dit =
663  py_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(
664  py_ptr->getNbLocalDofsRow()); // should be lower
665  break;
666  case COL:
667  y_dit =
668  py_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(0);
669  hi_y_dit =
670  py_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>().lower_bound(
671  py_ptr->getNbLocalDofsCol()); // should be lower
672  break;
673  default:
674  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
675  }
676  const NumeredDofEntityByUId *x_numered_dofs_by_uid;
677  switch (x_rc) {
678  case ROW:
679  x_numered_dofs_by_uid = &(px_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
680  break;
681  case COL:
682  x_numered_dofs_by_uid = &(px_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
683  break;
684  default:
685  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
686  }
687  for (; y_dit != hi_y_dit; y_dit++) {
688  if ((*y_dit)->getPart() != (unsigned int)m_field.get_comm_rank()) {
689  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
690  }
691  NumeredDofEntityByUId::iterator x_dit;
692  x_dit = x_numered_dofs_by_uid->find((*y_dit)->getLocalUniqueId());
693  if (x_dit == x_numered_dofs_by_uid->end())
694  continue;
695  idx.push_back((*x_dit)->getPetscGlobalDofIdx());
696  idy.push_back((*y_dit)->getPetscGlobalDofIdx());
697  }
699 }
700 
702  const std::string x_problem, RowColData x_rc, const std::string y_problem,
703  RowColData y_rc, IS *ix, IS *iy) const {
704  const MoFEM::Interface &m_field = cOre;
706  std::vector<int> idx(0), idy(0);
707  CHKERR isCreateFromProblemToOtherProblem(x_problem, x_rc, y_problem, y_rc,
708  idx, idy);
709  CHKERR ISCreateGeneral(m_field.get_comm(), idx.size(), &idx[0],
710  PETSC_COPY_VALUES, ix);
711  CHKERR ISCreateGeneral(m_field.get_comm(), idy.size(), &idy[0],
712  PETSC_COPY_VALUES, iy);
713  if (dEbug) {
714  ISView(*ix, PETSC_VIEWER_STDOUT_WORLD);
715  ISView(*iy, PETSC_VIEWER_STDOUT_WORLD);
716  }
718 }
719 } // namespace MoFEM
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::ISManager::isCreateProblemFieldAndRank
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:272
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::ISManager::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: ISManager.cpp:10
EntityHandle
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
MoFEM::ISManager::isCreateFromProblemFieldToOtherProblemField
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:523
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::CoreInterface::get_field_bit_number
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
MoFEM::PetscLocalIdx_mi_tag
Definition: TagMultiIndices.hpp:45
MoFEM::NumeredDofEntity::interface_type_DofEntity
interface_DofEntity< DofEntity > interface_type_DofEntity
Definition: DofsMultiIndices.hpp:232
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
MoFEM::Types::BitFieldId
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:42
MoFEM::CoreInterface::get_problem
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
MoFEM::Problem::getBitFEId
BitFEId getBitFEId() const
Get the BitFEIDs in problem
Definition: ProblemsMultiIndices.cpp:45
MoFEM::CoreInterface::get_finite_elements
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
MoFEM::CoreInterface::get_field_meshset
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
MoFEM::Order_mi_tag
MultiIndex Tag for field order.
Definition: TagMultiIndices.hpp:64
THROW_MESSAGE
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:574
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::ISManager::isCreateProblem
MoFEMErrorCode isCreateProblem(const std::string problem_name, RowColData rc, IS *is) const
Create IS for problem.
Definition: ISManager.cpp:149
ROW
@ ROW
Definition: definitions.h:136
MoFEM::ISManager::isCreateProblemBrokenFieldAndRankLocal
MoFEMErrorCode isCreateProblemBrokenFieldAndRankLocal(const std::vector< boost::weak_ptr< NumeredDofEntity >> &dofs_vec, SmartPetscObj< IS > &smart_is) const
IS for given problem, field and rank range (collective)
Definition: ISManager.cpp:381
MoFEM::NumeredDofEntityByUId
NumeredDofEntity_multiIndex::index< Unique_mi_tag >::type NumeredDofEntityByUId
Numbered DoF multi-index by UId.
Definition: DofsMultiIndices.hpp:476
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::NumeredDofEntity
keeps information about indexed dofs for the problem
Definition: DofsMultiIndices.hpp:226
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::Part_mi_tag
Definition: TagMultiIndices.hpp:53
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
MoFEM::ISManager::cOre
const MoFEM::Interface & cOre
Definition: ISManager.hpp:28
MoFEM::Problem::numeredColDofsPtr
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem
Definition: ProblemsMultiIndices.hpp:75
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
RowColData
RowColData
RowColData.
Definition: definitions.h:136
convert.type
type
Definition: convert.py:64
MoFEM::FieldEntity::getLoBitNumberUId
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
Definition: FieldEntsMultiIndices.hpp:222
COL
@ COL
Definition: definitions.h:136
MoFEM::interface_DofEntity< DofEntity >::getEnt
EntityHandle getEnt() const
Definition: DofsMultiIndices.hpp:192
MoFEM::Problem::getNbLocalDofsRow
DofIdx getNbLocalDofsRow() const
Definition: ProblemsMultiIndices.hpp:378
MoFEM::ISManager::dEbug
bool dEbug
Definition: ISManager.hpp:29
MoFEM::ISManager::isCreateProblemFieldAndRankLocal
MoFEMErrorCode isCreateProblemFieldAndRankLocal(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:406
MoFEM::DofEntity::getLoFieldEntityUId
static UId getLoFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
Definition: DofsMultiIndices.hpp:69
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::CoreInterface::get_fields
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
MoFEM::Problem::numeredRowDofsPtr
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
Definition: ProblemsMultiIndices.hpp:73
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
MoFEM::ISManager::ISManager
ISManager(const MoFEM::Core &core)
Definition: ISManager.cpp:17
Range
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
MoFEM::Types::ApproximationOrder
int ApproximationOrder
Approximation on the entity.
Definition: Types.hpp:26
MoFEM::interface_DofEntity< DofEntity >::getDofOrder
ApproximationOrder getDofOrder() const
Definition: DofsMultiIndices.hpp:195
MoFEM::Composite_Ent_And_EntDofIdx_mi_tag
Definition: TagMultiIndices.hpp:81
MoFEM::DofEntity::getHiFieldEntityUId
static UId getHiFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
Definition: DofsMultiIndices.hpp:75
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::ISManager::isCreateProblemBrokenFieldAndRank
MoFEMErrorCode isCreateProblemBrokenFieldAndRank(const std::vector< boost::weak_ptr< NumeredDofEntity >> &dofs_vec, SmartPetscObj< IS > &smart_is, MPI_Comm comm=PETSC_COMM_SELF) const
IS for given problem, field and rank range (collective)
Definition: ISManager.cpp:357
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
MoFEM::FieldEntity::getHiBitNumberUId
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
Definition: FieldEntsMultiIndices.hpp:228
MoFEM::Problem::getNbLocalDofsCol
DofIdx getNbLocalDofsCol() const
Definition: ProblemsMultiIndices.hpp:379
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::ISManager::isCreateFromProblemToOtherProblem
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:648
MoFEM::Unique_mi_tag
Definition: TagMultiIndices.hpp:18
MoFEM::ISManager::sectionCreate
MoFEMErrorCode sectionCreate(const std::string problem_name, PetscSection *s, const RowColData row_col=COL) const
Create global selection.
Definition: ISManager.cpp:20
MoFEM::ISManager::isCreateProblemFieldAndEntityType
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:505
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::SmartPetscObj
intrusive_ptr for managing petsc objects
Definition: PetscSmartObj.hpp:82
MoFEM::ISManager::isCreateProblemOrder
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:204
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::interface_DofEntity< DofEntity >::getEntDofIdx
DofIdx getEntDofIdx() const
Definition: DofsMultiIndices.hpp:186
MoFEM::Types::DofIdx
int DofIdx
Index of DOF.
Definition: Types.hpp:18