v0.13.2
Loading...
Searching...
No Matches
VecManager.cpp
Go to the documentation of this file.
1/** \file VecManager.cpp
2 * \brief Implementation of VecManager
3 */
4
5namespace MoFEM {
6
8VecManager::query_interface(boost::typeindex::type_index type_index,
9 UnknownInterface **iface) const {
10 *iface = const_cast<VecManager *>(this);
11 return 0;
12}
13
15 : cOre(const_cast<MoFEM::Core &>(core)), dEbug(false) {
16 if (!LogManager::checkIfChannelExist("VECWORLD")) {
17 auto core_log = logging::core::get();
18
19 core_log->add_sink(
21 core_log->add_sink(
23 core_log->add_sink(
25
26 LogManager::setLog("VECWORLD");
27 LogManager::setLog("VECSYNC");
28 LogManager::setLog("VECSELF");
29
30 MOFEM_LOG_TAG("VECWORLD", "VecMng");
31 MOFEM_LOG_TAG("VECSYNC", "VecMng");
32 MOFEM_LOG_TAG("VECSELF", "VecMng");
33 }
34
35 MOFEM_LOG("VECWORLD", Sev::noisy) << "Vector manager created";
36}
37
39 Vec *V) const {
40 const MoFEM::Interface &m_field = cOre;
41 const Problem *problem_ptr;
43 CHKERR m_field.get_problem(name, &problem_ptr);
44 DofIdx nb_local_dofs, nb_ghost_dofs;
45 switch (rc) {
46 case ROW:
47 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
48 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
49 break;
50 case COL:
51 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
52 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
53 break;
54 default:
55 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
56 }
57 CHKERR ::VecCreateSeq(PETSC_COMM_SELF, nb_local_dofs + nb_ghost_dofs, V);
58 CHKERR ::VecSetOption(
59 *V, VEC_IGNORE_NEGATIVE_INDICES,
60 PETSC_TRUE); // As default in MOFEM we will assume that this is always on
62}
63
65 Vec *V) const {
66 const MoFEM::Interface &m_field = cOre;
67 const Problem *problem_ptr;
69 CHKERR m_field.get_problem(name, &problem_ptr);
70 DofIdx nb_dofs, nb_local_dofs, nb_ghost_dofs;
72 switch (rc) {
73 case ROW:
74 nb_dofs = problem_ptr->getNbDofsRow();
75 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
76 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
77 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
78 &problem_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>());
79 break;
80 case COL:
81 nb_dofs = problem_ptr->getNbDofsCol();
82 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
83 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
84 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
85 &problem_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>());
86 break;
87 default:
88 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
89 }
90 // get ghost dofs
91 auto miit = dofs->lower_bound(nb_local_dofs);
92 auto hi_miit = dofs->upper_bound(nb_local_dofs + nb_ghost_dofs);
93 int count = std::distance(miit, hi_miit);
94 if (count != nb_ghost_dofs) {
95 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
96 }
97 std::vector<DofIdx> ghost_idx;
98 ghost_idx.reserve(count);
99 for (; miit != hi_miit; ++miit) {
100 if ((*miit)->getActive()) {
101 ghost_idx.emplace_back((*miit)->petscGloablDofIdx);
102
103#ifndef NDEBUG
104 if (PetscUnlikely(ghost_idx.back() > nb_dofs || ghost_idx.back() < 0)) {
106 MOFEM_LOG_ATTRIBUTES("VECSELF",
108 MOFEM_LOG("VECSELF", Sev::error) << "Problem with DOF " << **miit;
109 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong ghost index");
110 }
111#endif
112 }
113 }
114 CHKERR ::VecCreateGhost(m_field.get_comm(), nb_local_dofs, nb_dofs,
115 ghost_idx.size(), &*ghost_idx.begin(), V);
116 CHKERR ::VecSetOption(
117 *V, VEC_IGNORE_NEGATIVE_INDICES,
118 PETSC_TRUE); // As default in MOFEM we will assume that this is always on
120}
121
123 SmartPetscObj<Vec> &v_ptr) const {
125 Vec vec;
126 CHKERR vecCreateGhost(name, rc, &vec);
127 v_ptr.reset(vec, false);
129}
130
132 Vec xin, const std::string x_problem, const std::string x_field_name,
133 RowColData x_rc, Vec yin, const std::string y_problem,
134 const std::string y_field_name, RowColData y_rc, VecScatter *newctx) const {
135 const MoFEM::Interface &m_field = cOre;
137 std::vector<int> idx(0), idy(0);
138 CHKERR m_field.getInterface<ISManager>()
139 ->isCreateFromProblemFieldToOtherProblemField(
140 x_problem, x_field_name, x_rc, y_problem, y_field_name, y_rc, idx,
141 idy);
142 IS ix, iy;
143 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
144 PETSC_USE_POINTER, &ix);
145 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
146 PETSC_USE_POINTER, &iy);
147 CHKERR ::VecScatterCreate(xin, ix, yin, iy, newctx);
148 CHKERR ISDestroy(&ix);
149 CHKERR ISDestroy(&iy);
151}
152
154 Vec xin, const std::string x_problem, RowColData x_rc, Vec yin,
155 const std::string y_problem, RowColData y_rc, VecScatter *newctx) const {
156 const MoFEM::Interface &m_field = cOre;
158 std::vector<int> idx(0), idy(0);
159 CHKERR m_field.getInterface<ISManager>()->isCreateFromProblemToOtherProblem(
160 x_problem, x_rc, y_problem, y_rc, idx, idy);
161 IS ix, iy;
162 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
163 PETSC_USE_POINTER, &ix);
164 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
165 PETSC_USE_POINTER, &iy);
166 if (dEbug) {
167 ISView(ix, PETSC_VIEWER_STDOUT_WORLD);
168 ISView(iy, PETSC_VIEWER_STDOUT_WORLD);
169 }
170 CHKERR ::VecScatterCreate(xin, ix, yin, iy, newctx);
171 CHKERR ISDestroy(&ix);
172 CHKERR ISDestroy(&iy);
174}
175
177VecManager::vecScatterCreate(Vec xin, const std::string x_problem,
178 const std::string x_field_name, RowColData x_rc,
179 Vec yin, const std::string y_problem,
180 const std::string y_field_name, RowColData y_rc,
181 SmartPetscObj<VecScatter> &smart_newctx) const {
183 VecScatter newctx;
184 CHKERR vecScatterCreate(xin, x_problem, x_field_name, x_rc, yin, y_problem,
185 y_field_name, y_rc, &newctx);
186 smart_newctx = SmartPetscObj<VecScatter>(newctx);
188}
189
191VecManager::vecScatterCreate(Vec xin, const std::string x_problem,
192 RowColData x_rc, Vec yin,
193 const std::string y_problem, RowColData y_rc,
194 SmartPetscObj<VecScatter> &smart_newctx) const {
196 VecScatter newctx;
197 CHKERR vecScatterCreate(xin, x_problem, x_rc, yin, y_problem, y_rc, &newctx);
198 smart_newctx = SmartPetscObj<VecScatter>(newctx);
200}
201
203 RowColData rc, Vec V,
204 InsertMode mode,
205 ScatterMode scatter_mode) const {
208 DofIdx nb_local_dofs, nb_ghost_dofs;
209 switch (rc) {
210 case ROW:
211 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
212 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
213 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
214 &problem_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>());
215 break;
216 case COL:
217 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
218 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
219 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
220 &problem_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>());
221 break;
222 default:
223 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
224 }
225 Vec Vlocal;
226 CHKERR VecGhostGetLocalForm(V, &Vlocal);
227 int size;
228 CHKERR VecGetLocalSize(V, &size);
229 if (size != nb_local_dofs) {
230 SETERRQ3(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
231 "data inconsistency: check ghost vector, problem <%s> with nb. of "
232 "local nodes %d != %d",
233 problem_ptr->getName().c_str(), size, nb_local_dofs);
234 }
235 CHKERR VecGetLocalSize(Vlocal, &size);
236 if (size != nb_local_dofs + nb_ghost_dofs) {
237 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
238 "data inconsistency: check ghost vector, problem with nb. of "
239 "ghost nodes %d != ",
240 size, nb_local_dofs + nb_ghost_dofs);
241 }
242 NumeredDofEntityByLocalIdx::iterator miit = dofs->lower_bound(0);
243 NumeredDofEntityByLocalIdx::iterator hi_miit =
244 dofs->upper_bound(nb_local_dofs + nb_ghost_dofs);
245 DofIdx ii = 0;
246 switch (scatter_mode) {
247 case SCATTER_FORWARD: {
248 PetscScalar *array;
249 CHKERR VecGetArray(Vlocal, &array);
250 switch (mode) {
251 case INSERT_VALUES:
252 for (; miit != hi_miit; ++miit, ++ii)
253 array[ii] = (*miit)->getFieldData();
254 break;
255 case ADD_VALUES:
256 for (; miit != hi_miit; ++miit, ++ii)
257 array[ii] += (*miit)->getFieldData();
258 break;
259 default:
260 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
261 }
262 CHKERR VecRestoreArray(Vlocal, &array);
263 }; break;
264 case SCATTER_REVERSE: {
265 const PetscScalar *array;
266 CHKERR VecGetArrayRead(Vlocal, &array);
267 switch (mode) {
268 case INSERT_VALUES:
269 for (; miit != hi_miit; ++miit, ++ii) {
270 // std::cerr << *miit << std::endl;
271 // std::cerr << array[ii] << std::endl;
272 (*miit)->getFieldData() = array[ii];
273 }
274 break;
275 case ADD_VALUES:
276 for (; miit != hi_miit; ++miit, ++ii)
277 (*miit)->getFieldData() += array[ii];
278 break;
279 default:
280 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
281 }
282 CHKERR VecRestoreArrayRead(Vlocal, &array);
283 }; break;
284 default:
285 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
286 }
287 CHKERR VecGhostRestoreLocalForm(V, &Vlocal);
289}
290
292 RowColData rc, Vec V,
293 InsertMode mode,
294 ScatterMode scatter_mode) const {
295 const MoFEM::Interface &m_field = cOre;
296 const Problem *problem_ptr;
298 CHKERR m_field.get_problem(name, &problem_ptr);
299 CHKERR setLocalGhostVector(problem_ptr, rc, V, mode, scatter_mode);
301}
302
305 Vec V, InsertMode mode,
306 ScatterMode scatter_mode) const {
308 typedef NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type
309 DofsByGlobalIdx;
310 DofsByGlobalIdx *dofs;
311 DofIdx nb_dofs;
312 switch (rc) {
313 case ROW:
314 nb_dofs = problem_ptr->getNbDofsRow();
315 dofs = const_cast<DofsByGlobalIdx *>(
316 &problem_ptr->numeredRowDofsPtr->get<PetscGlobalIdx_mi_tag>());
317 break;
318 case COL:
319 nb_dofs = problem_ptr->getNbDofsCol();
320 dofs = const_cast<DofsByGlobalIdx *>(
321 &problem_ptr->numeredColDofsPtr->get<PetscGlobalIdx_mi_tag>());
322 break;
323 default:
324 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
325 "Function works only for ROWs and COLs");
326 }
327 DofsByGlobalIdx::iterator miit = dofs->lower_bound(0);
328 DofsByGlobalIdx::iterator hi_miit = dofs->upper_bound(nb_dofs);
329 auto comm = PetscObjectComm((PetscObject)V);
330 switch (scatter_mode) {
331 case SCATTER_REVERSE: {
332 VecScatter ctx;
333 Vec V_glob;
334 CHKERR VecScatterCreateToAll(V, &ctx, &V_glob);
335 CHKERR VecScatterBegin(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
336 CHKERR VecScatterEnd(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
337 int size;
338 CHKERR VecGetSize(V_glob, &size);
339 if (size != nb_dofs) {
340 SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
341 "Size of vector is inconsistent with size of problem. You could "
342 "use wrong vector with wrong problem, or you created vector "
343 "before you remove DOFs from problem.");
344 }
345 PetscScalar *array;
346 CHKERR VecGetArray(V_glob, &array);
347 switch (mode) {
348 case INSERT_VALUES:
349 for (; miit != hi_miit; miit++)
350 (*miit)->getFieldData() = array[(*miit)->getPetscGlobalDofIdx()];
351 break;
352 case ADD_VALUES:
353 for (; miit != hi_miit; miit++)
354 (*miit)->getFieldData() += array[(*miit)->getPetscGlobalDofIdx()];
355 break;
356 default:
357 SETERRQ(comm, MOFEM_NOT_IMPLEMENTED, "not implemented");
358 }
359 CHKERR VecRestoreArray(V_glob, &array);
360 CHKERR VecScatterDestroy(&ctx);
361 CHKERR VecDestroy(&V_glob);
362 break;
363 }
364 default:
365 SETERRQ(comm, MOFEM_NOT_IMPLEMENTED, "not implemented");
366 }
368}
369
371VecManager::setGlobalGhostVector(const std::string name, RowColData rc, Vec V,
372 InsertMode mode,
373 ScatterMode scatter_mode) const {
374 const MoFEM::Interface &m_field = cOre;
375 const Problem *problem_ptr;
377 CHKERR m_field.get_problem(name, &problem_ptr);
378 CHKERR setGlobalGhostVector(problem_ptr, rc, V, mode, scatter_mode);
380}
381
382template <int MODE> struct SetOtherLocalGhostVector {
383 template <typename A0, typename A1, typename A2, typename A3, typename A4>
384 inline MoFEMErrorCode operator()(A0 dofs_ptr, A1 array, A2 miit, A3 hi_miit,
385 A4 &cpy_bit_number,
386 const std::string cpy_field_name) {
388 for (; miit != hi_miit; miit++) {
389 if (miit->get()->getHasLocalIndex()) {
390 const auto uid = DofEntity::getUniqueIdCalculate(
391 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
392 cpy_bit_number, (*miit)->getEnt()));
393 auto diiiit = dofs_ptr->template get<Unique_mi_tag>().find(uid);
394 if (diiiit == dofs_ptr->template get<Unique_mi_tag>().end()) {
396 MOFEM_LOG_ATTRIBUTES("VECSELF",
398 MOFEM_LOG("VECSELF", Sev::error)
399 << "Problem finding DOFs in the copy field";
400 MOFEM_LOG("VECSELF", Sev::error) << "Field DOF: " << (**miit);
401 MOFEM_LOG("VECSELF", Sev::error)
402 << "Copy field name: " << cpy_field_name;
403 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
404 "Automatic creation of entity and dof not implemented");
405 }
406 if constexpr (MODE == INSERT_VALUES)
407 (*diiiit)->getFieldData() = array[(*miit)->getPetscLocalDofIdx()];
408 else if constexpr (MODE == ADD_VALUES)
409 (*diiiit)->getFieldData() += array[(*miit)->getPetscLocalDofIdx()];
410 }
411 }
413 }
414};
415
417 const Problem *problem_ptr, const std::string field_name,
418 const std::string cpy_field_name, RowColData rc, Vec V, InsertMode mode,
419 ScatterMode scatter_mode) const {
420 const MoFEM::Interface &m_field = cOre;
421 auto fields_ptr = m_field.get_fields();
422 auto dofs_ptr = m_field.get_dofs();
424 using DofsByUId = NumeredDofEntity_multiIndex::index<Unique_mi_tag>::type;
425 DofsByUId *dofs;
426
427 switch (rc) {
428 case ROW:
429 dofs = const_cast<DofsByUId *>(
430 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
431 break;
432 case COL:
433 dofs = const_cast<DofsByUId *>(
434 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
435 break;
436 default:
437 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
438 }
439
440 auto cpy_fit = fields_ptr->get<FieldName_mi_tag>().find(cpy_field_name);
441 if (cpy_fit == fields_ptr->get<FieldName_mi_tag>().end())
442 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
443 "cpy field < %s > not found, (top tip: check spelling)",
444 cpy_field_name.c_str());
445 const auto cpy_bit_number = (*cpy_fit)->getBitNumber();
446
447 auto miit = dofs->lower_bound(
449 if (miit == dofs->end()) {
450 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
451 "cpy field < %s > not found, (top tip: check spelling)",
452 field_name.c_str());
453 }
454 auto hi_miit = dofs->upper_bound(
456
457 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
458 SETERRQ4(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
459 "fields have to have same space (%s) %s != (%s) %s",
460 (*miit)->getName().c_str(), FieldSpaceNames[(*miit)->getSpace()],
461 cpy_field_name.c_str(), FieldSpaceNames[(*cpy_fit)->getSpace()]);
462 }
463 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
464 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
465 "fields have to have same rank");
466 }
467
468 switch (scatter_mode) {
469 case SCATTER_REVERSE: {
470
471 PetscScalar *array;
472 CHKERR VecGetArray(V, &array);
473 if (mode == INSERT_VALUES)
475 dofs_ptr, array, miit, hi_miit, cpy_bit_number, cpy_field_name);
476 else if (mode == ADD_VALUES)
478 dofs_ptr, array, miit, hi_miit, cpy_bit_number, cpy_field_name);
479 else
480 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
481 "Operation mode not implemented");
482 CHKERR VecRestoreArray(V, &array);
483 } break;
484 case SCATTER_FORWARD: {
485 for (; miit != hi_miit; miit++) {
486 if (!miit->get()->getHasLocalIndex())
487 continue;
488 const auto uid = DofEntity::getUniqueIdCalculate(
489 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
490 cpy_bit_number, (*miit)->getEnt()));
491 auto diiiit = dofs_ptr->get<Unique_mi_tag>().find(uid);
492 if (diiiit == dofs_ptr->get<Unique_mi_tag>().end()) {
493 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
494 "no data to fill the vector (top tip: you want scatter forward "
495 "or scatter reverse?)");
496 }
497 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
498 (*diiiit)->getFieldData(), mode);
499 }
500 CHKERR VecAssemblyBegin(V);
501 CHKERR VecAssemblyEnd(V);
502 } break;
503 default:
504 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "not implemented");
505 }
507}
508
510 const std::string name, const std::string field_name,
511 const std::string cpy_field_name, RowColData rc, Vec V, InsertMode mode,
512 ScatterMode scatter_mode) const {
513 const MoFEM::Interface &m_field = cOre;
514 const Problem *problem_ptr;
516 CHKERR m_field.get_problem(name, &problem_ptr);
517 CHKERR setOtherLocalGhostVector(problem_ptr, field_name, cpy_field_name, rc,
518 V, mode, scatter_mode);
520}
521
522template <int MODE> struct SetOtherGlobalGhostVector {
523 template <typename A0, typename A1, typename A2, typename A3, typename A4>
524 inline MoFEMErrorCode operator()(A0 dofs_ptr, A1 array, A2 miit, A3 hi_miit,
525 A4 &cpy_bit_number) {
527 for (; miit != hi_miit; miit++) {
528 const auto uid = DofEntity::getUniqueIdCalculate(
529 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
530 cpy_bit_number, (*miit)->getEnt()));
531 auto diiiit = dofs_ptr->template get<Unique_mi_tag>().find(uid);
532#ifndef NDEBUG
533 if (diiiit == dofs_ptr->template get<Unique_mi_tag>().end()) {
534 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
535 "Automatic creation of entity and dof not implemented");
536 }
537#endif
538 if constexpr (MODE == INSERT_VALUES)
539 (*diiiit)->getFieldData() = array[(*miit)->getPetscGlobalDofIdx()];
540 else if constexpr (MODE == ADD_VALUES)
541 (*diiiit)->getFieldData() += array[(*miit)->getPetscGlobalDofIdx()];
542 }
544 }
545};
546
548 const Problem *problem_ptr, const std::string field_name,
549 const std::string cpy_field_name, RowColData rc, Vec V, InsertMode mode,
550 ScatterMode scatter_mode) const {
551 const MoFEM::Interface &m_field = cOre;
552 auto fields_ptr = m_field.get_fields();
553 auto dofs_ptr = m_field.get_dofs();
554 auto *field_ents = m_field.get_field_ents();
557 DofIdx nb_dofs;
558 switch (rc) {
559 case ROW:
560 nb_dofs = problem_ptr->getNbDofsRow();
561 dofs = const_cast<NumeredDofEntityByUId *>(
562 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
563 break;
564 case COL:
565 nb_dofs = problem_ptr->getNbDofsCol();
566 dofs = const_cast<NumeredDofEntityByUId *>(
567 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
568 break;
569 default:
570 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "not implemented");
571 }
572 auto cpy_fit = fields_ptr->get<FieldName_mi_tag>().find(cpy_field_name);
573 if (cpy_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
574 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
575 "cpy field < %s > not found, (top tip: check spelling)",
576 cpy_field_name.c_str());
577 }
578 const auto cpy_bit_number = (*cpy_fit)->getBitNumber();
579
580 auto miit = dofs->lower_bound(
582 if (miit == dofs->end()) {
583 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
584 "problem field < %s > not found, (top tip: check spelling)",
585 field_name.c_str());
586 }
587 auto hi_miit = dofs->upper_bound(
589 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
590 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
591 "fields have to have same space");
592 }
593 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
594 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
595 "fields have to have same rank");
596 }
597 switch (scatter_mode) {
598 case SCATTER_REVERSE: {
599 Vec V_glob;
600 VecScatter ctx;
601 CHKERR VecScatterCreateToAll(V, &ctx, &V_glob);
602 CHKERR VecScatterBegin(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
603 CHKERR VecScatterEnd(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
604 int size;
605 CHKERR VecGetSize(V_glob, &size);
606 if (size != nb_dofs)
607 SETERRQ(
608 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
609 "data inconsistency: nb. of dofs and declared nb. dofs in database");
610 PetscScalar *array;
611 CHKERR VecGetArray(V_glob, &array);
612
613 if (mode == INSERT_VALUES)
615 dofs_ptr, array, miit, hi_miit, cpy_bit_number);
616 else if (mode == ADD_VALUES)
617 CHKERR SetOtherGlobalGhostVector<ADD_VALUES>()(dofs_ptr, array, miit,
618 hi_miit, cpy_bit_number);
619 else
620 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
621 "Operation mode not implemented");
622
623 CHKERR VecRestoreArray(V_glob, &array);
624 CHKERR VecDestroy(&V_glob);
625 CHKERR VecScatterDestroy(&ctx);
626 } break;
627 case SCATTER_FORWARD: {
628 for (; miit != hi_miit; miit++) {
629 const auto uid = DofEntity::getUniqueIdCalculate(
630 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
631 cpy_bit_number, (*miit)->getEnt()));
632 auto diiiit = dofs_ptr->get<Unique_mi_tag>().find(uid);
633 if (diiiit == dofs_ptr->get<Unique_mi_tag>().end()) {
634 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
635 "No data to fill the vector (top tip: you want scatter forward "
636 "or scatter reverse?)");
637 }
638 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
639 (*diiiit)->getFieldData(), mode);
640 }
641 CHKERR VecAssemblyBegin(V);
642 CHKERR VecAssemblyEnd(V);
643 } break;
644 default:
645 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
646 }
648}
649
651 const std::string name, const std::string field_name,
652 const std::string cpy_field_name, RowColData rc, Vec V, InsertMode mode,
653 ScatterMode scatter_mode) const {
654 const MoFEM::Interface &m_field = cOre;
655 const Problem *problem_ptr;
657 CHKERR m_field.get_problem(name, &problem_ptr);
658 CHKERR setOtherGlobalGhostVector(problem_ptr, field_name, cpy_field_name, rc,
659 V, mode, scatter_mode);
661}
662
663} // namespace MoFEM
RowColData
RowColData.
Definition: definitions.h:123
@ COL
Definition: definitions.h:123
@ ROW
Definition: definitions.h:123
static const char *const FieldSpaceNames[]
Definition: definitions.h:92
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
@ 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
NumeredDofEntity_multiIndex::index< Unique_mi_tag >::type NumeredDofEntityByUId
Numbered DoF multi-index by UId.
NumeredDofEntity_multiIndex::index< PetscLocalIdx_mi_tag >::type NumeredDofEntityByLocalIdx
Numbered DoF multi-index by local index.
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:391
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
Definition: LogManager.hpp:296
#define MOFEM_LOG_FUNCTION()
Set scope.
Definition: LogManager.hpp:325
MoFEMErrorCode vecCreateGhost(const std::string name, RowColData rc, Vec *V) const
create ghost vector for problem (collective)
Definition: VecManager.cpp:64
MoFEMErrorCode setOtherLocalGhostVector(const Problem *problem_ptr, const std::string field_name, const std::string cpy_field_name, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
Copy vector to field which is not part of the problem.
Definition: VecManager.cpp:416
MoFEMErrorCode setGlobalGhostVector(const Problem *problem_ptr, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
set values of vector from/to mesh database (collective)
Definition: VecManager.cpp:304
MoFEMErrorCode vecScatterCreate(Vec xin, const std::string x_problem, const std::string x_field_name, RowColData x_rc, Vec yin, const std::string y_problem, const std::string y_field_name, RowColData y_rc, VecScatter *newctx) const
create scatter for vectors form one to another problem (collective)
Definition: VecManager.cpp:131
MoFEMErrorCode vecCreateSeq(const std::string name, RowColData rc, Vec *V) const
create local vector for problem
Definition: VecManager.cpp:38
MoFEMErrorCode setLocalGhostVector(const Problem *problem_ptr, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
set values of vector from/to meshdatabase
Definition: VecManager.cpp:202
MoFEMErrorCode setOtherGlobalGhostVector(const Problem *problem_ptr, const std::string field_name, const std::string cpy_field_name, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
Copy vector to field which is not part of the problem (collective)
Definition: VecManager.cpp:547
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
int DofIdx
Index of DOF.
Definition: Types.hpp:18
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 MPI_Comm & get_comm() const =0
Core (interface) class.
Definition: Core.hpp:82
Deprecated interface functions.
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
MultiIndex Tag for field name.
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:300
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:346
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:350
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
Definition: LogManager.cpp:401
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:342
keeps basic data about problem
DofIdx getNbLocalDofsRow() const
DofIdx getNbDofsRow() const
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
DofIdx getNbGhostDofsCol() const
DofIdx getNbDofsCol() const
DofIdx getNbLocalDofsCol() const
DofIdx getNbGhostDofsRow() const
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem
MoFEMErrorCode operator()(A0 dofs_ptr, A1 array, A2 miit, A3 hi_miit, A4 &cpy_bit_number)
Definition: VecManager.cpp:524
MoFEMErrorCode operator()(A0 dofs_ptr, A1 array, A2 miit, A3 hi_miit, A4 &cpy_bit_number, const std::string cpy_field_name)
Definition: VecManager.cpp:384
intrusive_ptr for managing petsc objects
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: VecManager.cpp:8
const MoFEM::Interface & cOre
Definition: VecManager.hpp:28
VecManager(const MoFEM::Core &core)
Definition: VecManager.cpp:14