v0.15.5
Loading...
Searching...
No Matches
Files | Functions
Vectors (Vec)

Creating and scattering vectors on the mesh for given problem. More...

Collaboration diagram for Vectors (Vec):

Files

file  VecManager.hpp
 Interface managing vectors.
 

Functions

MoFEMErrorCode MoFEM::VecManager::vecCreateSeq (const std::string name, RowColData rc, Vec *V) const
 create local vector for problem
 
MoFEMErrorCode MoFEM::VecManager::vecCreateGhost (const std::string name, RowColData rc, Vec *V) const
 create ghost vector for problem (collective)
 
MoFEMErrorCode MoFEM::VecManager::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 from one to another problem (collective)
 
MoFEMErrorCode MoFEM::VecManager::vecScatterCreate (Vec xin, const std::string x_problem, RowColData x_rc, Vec yin, const std::string y_problem, RowColData y_rc, VecScatter *newctx) const
 create scatter for vectors from one to another problem (collective)
 
MoFEMErrorCode MoFEM::VecManager::setLocalGhostVector (const Problem *problem_ptr, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
 set values of vector from/to mesh database
 
MoFEMErrorCode MoFEM::VecManager::setLocalGhostVector (const std::string name, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
 set values of vector from/to mesh database
 
MoFEMErrorCode MoFEM::VecManager::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.
 
MoFEMErrorCode MoFEM::VecManager::setOtherLocalGhostVector (const std::string name, 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.
 
MoFEMErrorCode MoFEM::VecManager::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)
 
MoFEMErrorCode MoFEM::VecManager::setOtherGlobalGhostVector (const std::string name, 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)
 

Detailed Description

Creating and scattering vectors on the mesh for given problem.

Function Documentation

◆ setLocalGhostVector() [1/2]

MoFEMErrorCode MoFEM::VecManager::setLocalGhostVector ( const Problem problem_ptr,
RowColData  rc,
Vec  V,
InsertMode  mode,
ScatterMode  scatter_mode 
) const

#include <src/interfaces/VecManager.hpp>

set values of vector from/to mesh database

Parameters
pointerto problem structure
RowColDatafor row or column (ROW or COL)
Vvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

SCATTER_FORWARD set vector V from data field entities

Definition at line 233 of file VecManager.cpp.

236 {
239 DofIdx nb_local_dofs, nb_ghost_dofs;
240 switch (rc) {
241 case ROW:
242 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
243 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
244 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
245 &problem_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>());
246 break;
247 case COL:
248 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
249 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
250 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
251 &problem_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>());
252 break;
253 default:
254 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
255 }
256 Vec Vlocal;
257 CHKERR VecGhostGetLocalForm(V, &Vlocal);
258 int size;
259 CHKERR VecGetLocalSize(V, &size);
260 if (size != nb_local_dofs) {
261 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
262 "data inconsistency: check ghost vector, problem <%s> with nb. of "
263 "local nodes %d != %d",
264 problem_ptr->getName().c_str(), size, nb_local_dofs);
265 }
266 CHKERR VecGetLocalSize(Vlocal, &size);
267 if (size != nb_local_dofs + nb_ghost_dofs) {
268 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
269 "data inconsistency: check ghost vector, problem with nb. of "
270 "ghost nodes %d != %d ",
271 size, nb_local_dofs + nb_ghost_dofs);
272 }
273 auto miit = dofs->lower_bound(0);
274 auto hi_miit = dofs->upper_bound(nb_local_dofs + nb_ghost_dofs);
275 DofIdx ii = 0;
276 switch (scatter_mode) {
277 case SCATTER_FORWARD: {
278 PetscScalar *array;
279 CHKERR VecGetArray(Vlocal, &array);
280 switch (mode) {
281 case INSERT_VALUES:
282 for (; miit != hi_miit; ++miit, ++ii)
283 array[ii] = (*miit)->getFieldData();
284 break;
285 case ADD_VALUES:
286 for (; miit != hi_miit; ++miit, ++ii)
287 array[ii] += (*miit)->getFieldData();
288 break;
289 default:
290 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
291 }
292 CHKERR VecRestoreArray(Vlocal, &array);
293 }; break;
294 case SCATTER_REVERSE: {
295 const PetscScalar *array;
296 CHKERR VecGetArrayRead(Vlocal, &array);
297 switch (mode) {
298 case INSERT_VALUES:
299 for (; miit != hi_miit; ++miit, ++ii) {
300 (*miit)->getFieldData() = array[ii];
301 }
302 break;
303 case ADD_VALUES:
304 for (; miit != hi_miit; ++miit, ++ii)
305 (*miit)->getFieldData() += array[ii];
306 break;
307 default:
308 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
309 }
310 CHKERR VecRestoreArrayRead(Vlocal, &array);
311 }; break;
312 default:
313 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
314 }
315 CHKERR VecGhostRestoreLocalForm(V, &Vlocal);
317}
@ COL
@ ROW
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ 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()
#define CHKERR
Inline error check.
NumeredDofEntity_multiIndex::index< PetscLocalIdx_mi_tag >::type NumeredDofEntityByLocalIdx
Numbered DoF multi-index by local index.
const FTensor::Tensor2< T, Dim, Dim > Vec
int DofIdx
Index of DOF.
Definition Types.hpp:18

◆ setLocalGhostVector() [2/2]

MoFEMErrorCode MoFEM::VecManager::setLocalGhostVector ( const std::string  name,
RowColData  rc,
Vec  V,
InsertMode  mode,
ScatterMode  scatter_mode 
) const

#include <src/interfaces/VecManager.hpp>

set values of vector from/to mesh database

Parameters
nameof the problem
RowColDatafor row or column (ROW or COL)
Vvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

SCATTER_FORWARD set vector V from data field entities

Definition at line 319 of file VecManager.cpp.

322 {
323 const MoFEM::Interface &m_field = cOre;
324 const Problem *problem_ptr;
326 CHKERR m_field.get_problem(name, &problem_ptr);
327 CHKERR setLocalGhostVector(problem_ptr, rc, V, mode, scatter_mode);
329}
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
MoFEMErrorCode setLocalGhostVector(const Problem *problem_ptr, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
set values of vector from/to mesh database
Deprecated interface functions.
const MoFEM::Interface & cOre

◆ setOtherGlobalGhostVector() [1/2]

MoFEMErrorCode MoFEM::VecManager::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

#include <src/interfaces/VecManager.hpp>

Copy vector to field which is not part of the problem (collective)

collective - need to be run on all processors in communicator

\param problem_ptr pointer to problem
\param field_name field name used for indexing petsc vectors used in the

problem

Parameters
cpy_fieldfield name where data from vector are stored
RowColDatafor row or column
Vvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

Definition at line 575 of file VecManager.cpp.

578 {
579 const MoFEM::Interface &m_field = cOre;
580 auto fields_ptr = m_field.get_fields();
581 auto dofs_ptr = m_field.get_dofs();
584 DofIdx nb_dofs;
585 switch (rc) {
586 case ROW:
587 nb_dofs = problem_ptr->getNbDofsRow();
588 dofs = const_cast<NumeredDofEntityByUId *>(
589 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
590 break;
591 case COL:
592 nb_dofs = problem_ptr->getNbDofsCol();
593 dofs = const_cast<NumeredDofEntityByUId *>(
594 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
595 break;
596 default:
597 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "not implemented");
598 }
599 auto cpy_fit = fields_ptr->get<FieldName_mi_tag>().find(cpy_field_name);
600 if (cpy_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
601 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
602 "cpy field < %s > not found, (top tip: check spelling)",
603 cpy_field_name.c_str());
604 }
605 const auto cpy_bit_number = (*cpy_fit)->getBitNumber();
606
607 auto miit = dofs->lower_bound(
609 if (miit == dofs->end()) {
610 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
611 "problem field < %s > not found, (top tip: check spelling)",
612 field_name.c_str());
613 }
614 auto hi_miit = dofs->upper_bound(
616 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
617 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
618 "fields have to have same space");
619 }
620 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
621 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
622 "fields have to have same rank");
623 }
624 switch (scatter_mode) {
625 case SCATTER_REVERSE: {
626 Vec V_glob;
627 VecScatter ctx;
628 CHKERR VecScatterCreateToAll(V, &ctx, &V_glob);
629 CHKERR VecScatterBegin(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
630 CHKERR VecScatterEnd(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
631 int size;
632 CHKERR VecGetSize(V_glob, &size);
633 if (size != nb_dofs)
634 SETERRQ(
635 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
636 "data inconsistency: nb. of dofs and declared nb. dofs in database");
637 PetscScalar *array;
638 CHKERR VecGetArray(V_glob, &array);
639
640 if (mode == INSERT_VALUES)
641 CHKERR SetOtherGlobalGhostVector<INSERT_VALUES>()(
642 dofs_ptr, array, miit, hi_miit, cpy_bit_number);
643 else if (mode == ADD_VALUES)
644 CHKERR SetOtherGlobalGhostVector<ADD_VALUES>()(dofs_ptr, array, miit,
645 hi_miit, cpy_bit_number);
646 else
647 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
648 "Operation mode not implemented");
649
650 CHKERR VecRestoreArray(V_glob, &array);
651 CHKERR VecDestroy(&V_glob);
652 CHKERR VecScatterDestroy(&ctx);
653 } break;
654 case SCATTER_FORWARD: {
655 for (; miit != hi_miit; miit++) {
656 const auto uid = DofEntity::getUniqueIdCalculate(
657 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
658 cpy_bit_number, (*miit)->getEnt()));
659 auto diiiit = dofs_ptr->get<Unique_mi_tag>().find(uid);
660 if (diiiit == dofs_ptr->get<Unique_mi_tag>().end()) {
661 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
662 "No data to fill the vector (top tip: you want scatter forward "
663 "or scatter reverse?)");
664 }
665 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
666 (*diiiit)->getFieldData(), mode);
667 }
668 CHKERR VecAssemblyBegin(V);
669 CHKERR VecAssemblyEnd(V);
670 } break;
671 default:
672 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
673 }
675}
@ MOFEM_NOT_FOUND
Definition definitions.h:33
NumeredDofEntity_multiIndex::index< Unique_mi_tag >::type NumeredDofEntityByUId
Numbered DoF multi-index by UId.
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
constexpr auto field_name
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
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)

◆ setOtherGlobalGhostVector() [2/2]

MoFEMErrorCode MoFEM::VecManager::setOtherGlobalGhostVector ( const std::string  name,
const std::string  field_name,
const std::string  cpy_field_name,
RowColData  rc,
Vec  V,
InsertMode  mode,
ScatterMode  scatter_mode 
) const

#include <src/interfaces/VecManager.hpp>

Copy vector to field which is not part of the problem (collective)

Deprecated:
VecManager

collective - need to be run on all processors in communicator

\param name problem name
\param field_name field name used for indexing petsc vectors used in the

problem

Parameters
cpy_fieldfield name where data from vector are stored
RowColDatafor row or column
Vvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

Definition at line 677 of file VecManager.cpp.

680 {
681 const MoFEM::Interface &m_field = cOre;
682 const Problem *problem_ptr;
684 CHKERR m_field.get_problem(name, &problem_ptr);
685 CHKERR setOtherGlobalGhostVector(problem_ptr, field_name, cpy_field_name, rc,
686 V, mode, scatter_mode);
688}
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)

◆ setOtherLocalGhostVector() [1/2]

MoFEMErrorCode MoFEM::VecManager::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

#include <src/interfaces/VecManager.hpp>

Copy vector to field which is not part of the problem.

Parameters
pointerto problem multi_index
field_namefield name used for indexing petsc vectors used in the problem
cpy_fieldfield name where data from vector are stored
RowColDatafor row or column
Vvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

Definition at line 444 of file VecManager.cpp.

447 {
448 const MoFEM::Interface &m_field = cOre;
449 auto fields_ptr = m_field.get_fields();
450 auto dofs_ptr = m_field.get_dofs();
452 using DofsByUId = NumeredDofEntity_multiIndex::index<Unique_mi_tag>::type;
453 DofsByUId *dofs;
454
455 switch (rc) {
456 case ROW:
457 dofs = const_cast<DofsByUId *>(
458 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
459 break;
460 case COL:
461 dofs = const_cast<DofsByUId *>(
462 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
463 break;
464 default:
465 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
466 }
467
468 auto cpy_fit = fields_ptr->get<FieldName_mi_tag>().find(cpy_field_name);
469 if (cpy_fit == fields_ptr->get<FieldName_mi_tag>().end())
470 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
471 "cpy field < %s > not found, (top tip: check spelling)",
472 cpy_field_name.c_str());
473 const auto cpy_bit_number = (*cpy_fit)->getBitNumber();
474
475 auto miit = dofs->lower_bound(
477 if (miit == dofs->end()) {
478 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
479 "cpy field < %s > not found, (top tip: check spelling)",
480 field_name.c_str());
481 }
482 auto hi_miit = dofs->upper_bound(
484
485 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
486 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
487 "fields have to have same space (%s) %s != (%s) %s",
488 (*miit)->getName().c_str(), FieldSpaceNames[(*miit)->getSpace()],
489 cpy_field_name.c_str(), FieldSpaceNames[(*cpy_fit)->getSpace()]);
490 }
491 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
492 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
493 "fields have to have same rank");
494 }
495
496 switch (scatter_mode) {
497 case SCATTER_REVERSE: {
498
499 PetscScalar *array;
500 CHKERR VecGetArray(V, &array);
501 if (mode == INSERT_VALUES)
502 CHKERR SetOtherLocalGhostVector<INSERT_VALUES>()(
503 dofs_ptr, array, miit, hi_miit, cpy_bit_number, cpy_field_name);
504 else if (mode == ADD_VALUES)
505 CHKERR SetOtherLocalGhostVector<ADD_VALUES>()(
506 dofs_ptr, array, miit, hi_miit, cpy_bit_number, cpy_field_name);
507 else
508 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
509 "Operation mode not implemented");
510 CHKERR VecRestoreArray(V, &array);
511 } break;
512 case SCATTER_FORWARD: {
513 for (; miit != hi_miit; miit++) {
514 if (!miit->get()->getHasLocalIndex())
515 continue;
516 const auto uid = DofEntity::getUniqueIdCalculate(
517 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
518 cpy_bit_number, (*miit)->getEnt()));
519 auto diiiit = dofs_ptr->get<Unique_mi_tag>().find(uid);
520 if (diiiit == dofs_ptr->get<Unique_mi_tag>().end()) {
521 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
522 "no data to fill the vector (top tip: you want scatter forward "
523 "or scatter reverse?)");
524 }
525 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
526 (*diiiit)->getFieldData(), mode);
527 }
528 CHKERR VecAssemblyBegin(V);
529 CHKERR VecAssemblyEnd(V);
530 } break;
531 default:
532 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "not implemented");
533 }
535}
static const char *const FieldSpaceNames[]
Definition definitions.h:92

◆ setOtherLocalGhostVector() [2/2]

MoFEMErrorCode MoFEM::VecManager::setOtherLocalGhostVector ( const std::string  name,
const std::string  field_name,
const std::string  cpy_field_name,
RowColData  rc,
Vec  V,
InsertMode  mode,
ScatterMode  scatter_mode 
) const

#include <src/interfaces/VecManager.hpp>

Copy vector to field which is not part of the problem.

Parameters
nameproblem name
field_namefield name used for indexing petsc vectors used in the problem
cpy_fieldfield name where data from vector are stored
RowColDatafor row or column
Vvector
modesee petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
scatter_modesee petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE)

SCATTER_REVERSE set data to field entities from V vector.

Definition at line 537 of file VecManager.cpp.

540 {
541 const MoFEM::Interface &m_field = cOre;
542 const Problem *problem_ptr;
544 CHKERR m_field.get_problem(name, &problem_ptr);
545 CHKERR setOtherLocalGhostVector(problem_ptr, field_name, cpy_field_name, rc,
546 V, mode, scatter_mode);
548}
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.

◆ vecCreateGhost()

MoFEMErrorCode MoFEM::VecManager::vecCreateGhost ( const std::string  name,
RowColData  rc,
Vec *  V 
) const

#include <src/interfaces/VecManager.hpp>

create ghost vector for problem (collective)

collective - need to be run on all processors in communicator

Parameters
nameproblem name
RowColDataspecify what data is taken from Row, Col or Data
Vecthe vector where data is stored

Definition at line 75 of file VecManager.cpp.

76 {
77 const MoFEM::Interface &m_field = cOre;
78 const Problem *problem_ptr;
80 CHKERR m_field.get_problem(name, &problem_ptr);
81 DofIdx nb_dofs, nb_local_dofs, nb_ghost_dofs;
83 switch (rc) {
84 case ROW:
85 nb_dofs = problem_ptr->getNbDofsRow();
86 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
87 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
88 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
89 &problem_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>());
90 break;
91 case COL:
92 nb_dofs = problem_ptr->getNbDofsCol();
93 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
94 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
95 dofs = const_cast<NumeredDofEntityByLocalIdx *>(
96 &problem_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>());
97 break;
98 default:
99 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
100 }
101 // get ghost dofs
102 auto miit = dofs->lower_bound(nb_local_dofs);
103 auto hi_miit = dofs->upper_bound(nb_local_dofs + nb_ghost_dofs);
104 int count = std::distance(miit, hi_miit);
105 if (count != nb_ghost_dofs) {
106 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
107 }
108 std::vector<DofIdx> ghost_idx;
109 ghost_idx.reserve(count);
110 for (; miit != hi_miit; ++miit) {
111 if ((*miit)->getActive()) {
112 ghost_idx.emplace_back((*miit)->petscGloablDofIdx);
113
114 if (PetscUnlikely(ghost_idx.back() > nb_dofs || ghost_idx.back() < 0)) {
115
116#ifndef NDEBUG
117 auto save_no_global_id_entities = [&]() {
119 auto part = m_field.get_comm_rank();
120 Range ghost_range;
121 for (auto &it : (*dofs)) {
122 int global_idx = it->getPetscGlobalDofIdx();
123 if (global_idx > nb_dofs || global_idx < 0) {
124 ghost_range.insert(it->getEnt());
125 }
126 }
127 CHKERR save_range(const_cast<MoFEM::Interface &>(m_field).get_moab(),
128 "error_ghost_dofs_range_" + std::to_string(part) +
129 ".vtk",
130 ghost_range);
132 };
133
134 CHKERR save_no_global_id_entities();
135#endif
136
138 MOFEM_LOG_ATTRIBUTES("VECSELF",
140 MOFEM_LOG("VECSELF", Sev::error) << "Problem with DOF " << **miit;
141 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong ghost index");
142 }
143 }
144 }
145 CHKERR ::VecCreateGhost(m_field.get_comm(), nb_local_dofs, nb_dofs,
146 ghost_idx.size(), &*ghost_idx.begin(), V);
147 CHKERR ::VecSetOption(
148 *V, VEC_IGNORE_NEGATIVE_INDICES,
149 PETSC_TRUE); // As default in MOFEM we will assume that this is always on
151}
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
#define MOFEM_LOG_FUNCTION()
Set scope.
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
DofIdx getNbDofsRow() const
auto save_range

◆ vecCreateSeq()

MoFEMErrorCode MoFEM::VecManager::vecCreateSeq ( const std::string  name,
RowColData  rc,
Vec *  V 
) const

#include <src/interfaces/VecManager.hpp>

create local vector for problem

Parameters
nameproblem name
RowColDataspecify what data is taken from Row, Col or Data
Vecthe vector where data is stored

Definition at line 49 of file VecManager.cpp.

50 {
51 const MoFEM::Interface &m_field = cOre;
52 const Problem *problem_ptr;
54 CHKERR m_field.get_problem(name, &problem_ptr);
55 DofIdx nb_local_dofs, nb_ghost_dofs;
56 switch (rc) {
57 case ROW:
58 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
59 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
60 break;
61 case COL:
62 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
63 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
64 break;
65 default:
66 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
67 }
68 CHKERR ::VecCreateSeq(PETSC_COMM_SELF, nb_local_dofs + nb_ghost_dofs, V);
69 CHKERR ::VecSetOption(
70 *V, VEC_IGNORE_NEGATIVE_INDICES,
71 PETSC_TRUE); // As default in MOFEM we will assume that this is always on
73}
DofIdx getNbLocalDofsRow() const

◆ vecScatterCreate() [1/2]

MoFEMErrorCode MoFEM::VecManager::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

#include <src/interfaces/VecManager.hpp>

create scatter for vectors from one to another problem (collective)

User specify what name of field on one problem is scattered to another.

Parameters
xinvector
x_problemproblem name
x_fieldname
yinvector
y_problemproblem name
y_field_name
Return values
newctxscatter

Definition at line 162 of file VecManager.cpp.

165 {
166 const MoFEM::Interface &m_field = cOre;
168 std::vector<int> idx(0), idy(0);
169 CHKERR m_field.getInterface<ISManager>()
170 ->isCreateFromProblemFieldToOtherProblemField(
171 x_problem, x_field_name, x_rc, y_problem, y_field_name, y_rc, idx,
172 idy);
173 IS ix, iy;
174 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
175 PETSC_USE_POINTER, &ix);
176 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
177 PETSC_USE_POINTER, &iy);
178 CHKERR ::VecScatterCreate(xin, ix, yin, iy, newctx);
179 CHKERR ISDestroy(&ix);
180 CHKERR ISDestroy(&iy);
182}
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ vecScatterCreate() [2/2]

MoFEMErrorCode MoFEM::VecManager::vecScatterCreate ( Vec  xin,
const std::string  x_problem,
RowColData  x_rc,
Vec  yin,
const std::string  y_problem,
RowColData  y_rc,
VecScatter *  newctx 
) const

#include <src/interfaces/VecManager.hpp>

create scatter for vectors from one to another problem (collective)

Parameters
xinvector
x_problemproblem name
yinvector
y_problemproblem name
Return values
newctxscatter

Definition at line 184 of file VecManager.cpp.

186 {
187 const MoFEM::Interface &m_field = cOre;
189 std::vector<int> idx(0), idy(0);
190 CHKERR m_field.getInterface<ISManager>()->isCreateFromProblemToOtherProblem(
191 x_problem, x_rc, y_problem, y_rc, idx, idy);
192 IS ix, iy;
193 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idx.size(), &idx[0],
194 PETSC_USE_POINTER, &ix);
195 CHKERR ISCreateGeneral(PETSC_COMM_WORLD, idy.size(), &idy[0],
196 PETSC_USE_POINTER, &iy);
197 if (dEbug) {
198 ISView(ix, PETSC_VIEWER_STDOUT_WORLD);
199 ISView(iy, PETSC_VIEWER_STDOUT_WORLD);
200 }
201 CHKERR ::VecScatterCreate(xin, ix, yin, iy, newctx);
202 CHKERR ISDestroy(&ix);
203 CHKERR ISDestroy(&iy);
205}