|
| v0.14.0
|
Creating and scattering vectors on the mesh for given problem.
More...
|
MoFEMErrorCode | MoFEM::VecManager::vecCreateSeq (const std::string name, RowColData rc, Vec *V) const |
| create local vector for problem More...
|
|
MoFEMErrorCode | MoFEM::VecManager::vecCreateGhost (const std::string name, RowColData rc, Vec *V) const |
| create ghost vector for problem (collective) More...
|
|
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 form one to another problem (collective) More...
|
|
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) More...
|
|
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 More...
|
|
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 More...
|
|
MoFEMErrorCode | MoFEM::VecManager::setGlobalGhostVector (const Problem *problem_ptr, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const |
| set values of vector from/to mesh database (collective) More...
|
|
MoFEMErrorCode | MoFEM::VecManager::setGlobalGhostVector (const std::string name, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const |
| set values of vector from/to mesh database (collective) More...
|
|
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. More...
|
|
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. More...
|
|
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) More...
|
|
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) More...
|
|
Creating and scattering vectors on the mesh for given problem.
◆ setGlobalGhostVector() [1/2]
MoFEMErrorCode MoFEM::VecManager::setGlobalGhostVector |
( |
const Problem * |
problem_ptr, |
|
|
RowColData |
rc, |
|
|
Vec |
V, |
|
|
InsertMode |
mode, |
|
|
ScatterMode |
scatter_mode |
|
) |
| const |
set values of vector from/to mesh database (collective)
collective - need tu be run on all processors in communicator
- Parameters
-
pointer | to porblem struture |
RowColData | for row or column (i.e. Row,Col) |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities form V vector.
Definition at line 301 of file VecManager.cpp.
307 DofsByGlobalIdx *dofs;
311 nb_dofs = problem_ptr->getNbDofsRow();
312 dofs =
const_cast<DofsByGlobalIdx *
>(
313 &problem_ptr->numeredRowDofsPtr->get<PetscGlobalIdx_mi_tag>());
316 nb_dofs = problem_ptr->getNbDofsCol();
317 dofs =
const_cast<DofsByGlobalIdx *
>(
318 &problem_ptr->numeredColDofsPtr->get<PetscGlobalIdx_mi_tag>());
322 "Function works only for ROWs and COLs");
324 DofsByGlobalIdx::iterator miit = dofs->lower_bound(0);
325 DofsByGlobalIdx::iterator hi_miit = dofs->upper_bound(nb_dofs);
326 auto comm = PetscObjectComm((PetscObject)V);
327 switch (scatter_mode) {
328 case SCATTER_REVERSE: {
331 CHKERR VecScatterCreateToAll(V, &ctx, &V_glob);
332 CHKERR VecScatterBegin(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
333 CHKERR VecScatterEnd(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
335 CHKERR VecGetSize(V_glob, &size);
336 if (size != nb_dofs) {
338 "Size of vector is inconsistent with size of problem. You could "
339 "use wrong vector with wrong problem, or you created vector "
340 "before you remove DOFs from problem.");
343 CHKERR VecGetArray(V_glob, &array);
346 for (; miit != hi_miit; miit++)
347 (*miit)->getFieldData() = array[(*miit)->getPetscGlobalDofIdx()];
350 for (; miit != hi_miit; miit++)
351 (*miit)->getFieldData() += array[(*miit)->getPetscGlobalDofIdx()];
356 CHKERR VecRestoreArray(V_glob, &array);
357 CHKERR VecScatterDestroy(&ctx);
358 CHKERR VecDestroy(&V_glob);
◆ setGlobalGhostVector() [2/2]
MoFEMErrorCode MoFEM::VecManager::setGlobalGhostVector |
( |
const std::string |
name, |
|
|
RowColData |
rc, |
|
|
Vec |
V, |
|
|
InsertMode |
mode, |
|
|
ScatterMode |
scatter_mode |
|
) |
| const |
set values of vector from/to mesh database (collective)
collective - need tu be run on all processors in communicator
- Parameters
-
name | of the problem |
RowColData | for row or column (i.e. Row,Col) |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities form V vector.
Definition at line 368 of file VecManager.cpp.
372 const Problem *problem_ptr;
◆ setLocalGhostVector() [1/2]
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
- Parameters
-
pointer | to problem struture |
RowColData | for row or column:e (i.e. Row,Col) |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see 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 202 of file VecManager.cpp.
208 DofIdx nb_local_dofs, nb_ghost_dofs;
211 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
212 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
214 &problem_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>());
217 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
218 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
220 &problem_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>());
226 CHKERR VecGhostGetLocalForm(V, &Vlocal);
228 CHKERR VecGetLocalSize(V, &size);
229 if (size != nb_local_dofs) {
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);
235 CHKERR VecGetLocalSize(Vlocal, &size);
236 if (size != nb_local_dofs + nb_ghost_dofs) {
238 "data inconsistency: check ghost vector, problem with nb. of "
239 "ghost nodes %d != ",
240 size, nb_local_dofs + nb_ghost_dofs);
242 auto miit = dofs->lower_bound(0);
243 auto hi_miit = dofs->upper_bound(nb_local_dofs + nb_ghost_dofs);
245 switch (scatter_mode) {
246 case SCATTER_FORWARD: {
248 CHKERR VecGetArray(Vlocal, &array);
251 for (; miit != hi_miit; ++miit, ++ii)
252 array[ii] = (*miit)->getFieldData();
255 for (; miit != hi_miit; ++miit, ++ii)
256 array[ii] += (*miit)->getFieldData();
261 CHKERR VecRestoreArray(Vlocal, &array);
263 case SCATTER_REVERSE: {
264 const PetscScalar *array;
265 CHKERR VecGetArrayRead(Vlocal, &array);
268 for (; miit != hi_miit; ++miit, ++ii) {
269 (*miit)->getFieldData() = array[ii];
273 for (; miit != hi_miit; ++miit, ++ii)
274 (*miit)->getFieldData() += array[ii];
279 CHKERR VecRestoreArrayRead(Vlocal, &array);
284 CHKERR VecGhostRestoreLocalForm(V, &Vlocal);
◆ setLocalGhostVector() [2/2]
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
- Parameters
-
name | of the problem |
RowColData | for row or column:e (i.e. Row,Col) |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see 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 288 of file VecManager.cpp.
293 const Problem *problem_ptr;
◆ 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 |
Copy vector to field which is not part of the problem (collective)
collective - need tu be run on all processors in communicator
- Parameters
-
problem_ptr | pointer to problem |
field_name | field name used for indexing petsc vectors used in the problem |
cpy_field | field name where data from vector are stored |
RowColData | for row or column |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities form V vector.
Definition at line 544 of file VecManager.cpp.
556 nb_dofs = problem_ptr->getNbDofsRow();
558 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
561 nb_dofs = problem_ptr->getNbDofsCol();
563 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
568 auto cpy_fit = fields_ptr->get<FieldName_mi_tag>().find(cpy_field_name);
569 if (cpy_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
571 "cpy field < %s > not found, (top tip: check spelling)",
572 cpy_field_name.c_str());
574 const auto cpy_bit_number = (*cpy_fit)->getBitNumber();
576 auto miit = dofs->lower_bound(
578 if (miit == dofs->end()) {
580 "problem field < %s > not found, (top tip: check spelling)",
583 auto hi_miit = dofs->upper_bound(
585 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
587 "fields have to have same space");
589 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
591 "fields have to have same rank");
593 switch (scatter_mode) {
594 case SCATTER_REVERSE: {
597 CHKERR VecScatterCreateToAll(V, &ctx, &V_glob);
598 CHKERR VecScatterBegin(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
599 CHKERR VecScatterEnd(ctx, V, V_glob, INSERT_VALUES, SCATTER_FORWARD);
601 CHKERR VecGetSize(V_glob, &size);
605 "data inconsistency: nb. of dofs and declared nb. dofs in database");
607 CHKERR VecGetArray(V_glob, &array);
609 if (mode == INSERT_VALUES)
610 CHKERR SetOtherGlobalGhostVector<INSERT_VALUES>()(
611 dofs_ptr, array, miit, hi_miit, cpy_bit_number);
612 else if (mode == ADD_VALUES)
613 CHKERR SetOtherGlobalGhostVector<ADD_VALUES>()(dofs_ptr, array, miit,
614 hi_miit, cpy_bit_number);
617 "Operation mode not implemented");
619 CHKERR VecRestoreArray(V_glob, &array);
620 CHKERR VecDestroy(&V_glob);
621 CHKERR VecScatterDestroy(&ctx);
623 case SCATTER_FORWARD: {
624 for (; miit != hi_miit; miit++) {
627 cpy_bit_number, (*miit)->getEnt()));
628 auto diiiit = dofs_ptr->get<Unique_mi_tag>().find(uid);
629 if (diiiit == dofs_ptr->get<Unique_mi_tag>().end()) {
631 "No data to fill the vector (top tip: you want scatter forward "
632 "or scatter reverse?)");
634 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
635 (*diiiit)->getFieldData(), mode);
637 CHKERR VecAssemblyBegin(V);
◆ 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 |
Copy vector to field which is not part of the problem (collective)
- Deprecated:
- VecManager
collective - need tu be run on all processors in communicator
- Parameters
-
name | problem name |
field_name | field name used for indexing petsc vectors used in the problem |
cpy_field | field name where data from vector are stored |
RowColData | for row or column |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities form V vector.
Definition at line 646 of file VecManager.cpp.
651 const Problem *problem_ptr;
655 V, mode, scatter_mode);
◆ 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 |
Copy vector to field which is not part of the problem.
- Parameters
-
pointer | to problem multi_index |
field_name | field name used for indexing petsc vectors used in the problem |
cpy_field | field name where data from vector are stored |
RowColData | for row or column |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities form V vector.
Definition at line 413 of file VecManager.cpp.
426 dofs =
const_cast<DofsByUId *
>(
427 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
430 dofs =
const_cast<DofsByUId *
>(
431 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
437 auto cpy_fit = fields_ptr->get<FieldName_mi_tag>().find(cpy_field_name);
438 if (cpy_fit == fields_ptr->get<FieldName_mi_tag>().end())
440 "cpy field < %s > not found, (top tip: check spelling)",
441 cpy_field_name.c_str());
442 const auto cpy_bit_number = (*cpy_fit)->getBitNumber();
444 auto miit = dofs->lower_bound(
446 if (miit == dofs->end()) {
448 "cpy field < %s > not found, (top tip: check spelling)",
451 auto hi_miit = dofs->upper_bound(
454 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
456 "fields have to have same space (%s) %s != (%s) %s",
460 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
462 "fields have to have same rank");
465 switch (scatter_mode) {
466 case SCATTER_REVERSE: {
469 CHKERR VecGetArray(V, &array);
470 if (mode == INSERT_VALUES)
471 CHKERR SetOtherLocalGhostVector<INSERT_VALUES>()(
472 dofs_ptr, array, miit, hi_miit, cpy_bit_number, cpy_field_name);
473 else if (mode == ADD_VALUES)
474 CHKERR SetOtherLocalGhostVector<ADD_VALUES>()(
475 dofs_ptr, array, miit, hi_miit, cpy_bit_number, cpy_field_name);
478 "Operation mode not implemented");
479 CHKERR VecRestoreArray(V, &array);
481 case SCATTER_FORWARD: {
482 for (; miit != hi_miit; miit++) {
483 if (!miit->get()->getHasLocalIndex())
487 cpy_bit_number, (*miit)->getEnt()));
488 auto diiiit = dofs_ptr->get<Unique_mi_tag>().find(uid);
489 if (diiiit == dofs_ptr->get<Unique_mi_tag>().end()) {
491 "no data to fill the vector (top tip: you want scatter forward "
492 "or scatter reverse?)");
494 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
495 (*diiiit)->getFieldData(), mode);
497 CHKERR VecAssemblyBegin(V);
◆ 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 |
Copy vector to field which is not part of the problem.
- Parameters
-
name | problem name |
field_name | field name used for indexing petsc vectors used in the problem |
cpy_field | field name where data from vector are stored |
RowColData | for row or column |
V | vector |
mode | see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES) |
scatter_mode | see petsc manual for ScatterMode (The available modes are: SCATTER_FORWARD or SCATTER_REVERSE) |
SCATTER_REVERSE set data to field entities form V vector.
Definition at line 506 of file VecManager.cpp.
511 const Problem *problem_ptr;
515 V, mode, scatter_mode);
◆ vecCreateGhost()
create ghost vector for problem (collective)
collective - need to be run on all processors in communicator
- Parameters
-
name | problem name |
RowColData | specify what data is taken from Row, Col or Data |
Vec | the vector where data is stored |
Definition at line 64 of file VecManager.cpp.
67 const Problem *problem_ptr;
70 DofIdx nb_dofs, nb_local_dofs, nb_ghost_dofs;
74 nb_dofs = problem_ptr->getNbDofsRow();
75 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
76 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
78 &problem_ptr->numeredRowDofsPtr->get<PetscLocalIdx_mi_tag>());
81 nb_dofs = problem_ptr->getNbDofsCol();
82 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
83 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
85 &problem_ptr->numeredColDofsPtr->get<PetscLocalIdx_mi_tag>());
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) {
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);
104 if (PetscUnlikely(ghost_idx.back() > nb_dofs || ghost_idx.back() < 0)) {
108 MOFEM_LOG(
"VECSELF", Sev::error) <<
"Problem with DOF " << **miit;
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,
◆ vecCreateSeq()
create local vector for problem
- Parameters
-
name | problem name |
RowColData | specify what data is taken from Row, Col or Data |
Vec | the vector where data is stored |
Definition at line 38 of file VecManager.cpp.
41 const Problem *problem_ptr;
44 DofIdx nb_local_dofs, nb_ghost_dofs;
47 nb_local_dofs = problem_ptr->getNbLocalDofsRow();
48 nb_ghost_dofs = problem_ptr->getNbGhostDofsRow();
51 nb_local_dofs = problem_ptr->getNbLocalDofsCol();
52 nb_ghost_dofs = problem_ptr->getNbGhostDofsCol();
57 CHKERR ::VecCreateSeq(PETSC_COMM_SELF, nb_local_dofs + nb_ghost_dofs, V);
58 CHKERR ::VecSetOption(
59 *V, VEC_IGNORE_NEGATIVE_INDICES,
◆ 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 |
create scatter for vectors form one to another problem (collective)
User specify what name of field on one problem is scattered to another.
- Parameters
-
xin | vector |
x_proble | problem name |
x_field | name |
yin | vector |
y_problem | problem name |
y_field_name | |
- Return values
-
Definition at line 131 of file VecManager.cpp.
137 std::vector<int> idx(0), idy(0);
139 ->isCreateFromProblemFieldToOtherProblemField(
140 x_problem, x_field_name, x_rc, y_problem, y_field_name, y_rc, idx,
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);
◆ 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 |
create scatter for vectors from one to another problem (collective)
- Parameters
-
xin | vector |
x_proble | problem name |
yin | vector |
y_problem | problem name |
- Return values
-
Definition at line 153 of file VecManager.cpp.
158 std::vector<int> idx(0), idy(0);
160 x_problem, x_rc, y_problem, y_rc, idx, idy);
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);
167 ISView(ix, PETSC_VIEWER_STDOUT_WORLD);
168 ISView(iy, PETSC_VIEWER_STDOUT_WORLD);
170 CHKERR ::VecScatterCreate(xin, ix, yin, iy, newctx);
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
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.
virtual MPI_Comm & get_comm() const =0
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
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)
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
NumeredDofEntity_multiIndex::index< PetscLocalIdx_mi_tag >::type NumeredDofEntityByLocalIdx
Numbered DoF multi-index by local index.
Deprecated interface functions.
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
#define MOFEM_LOG_FUNCTION()
Set scope.
NumeredDofEntity_multiIndex::index< Unique_mi_tag >::type NumeredDofEntityByUId
Numbered DoF multi-index by UId.
#define CHKERR
Inline error check.
const MoFEM::Interface & cOre
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
constexpr auto field_name
MoFEMErrorCode setLocalGhostVector(const Problem *problem_ptr, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
set values of vector from/to mesh database
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
const static char *const FieldSpaceNames[]
#define MOFEM_LOG(channel, severity)
Log.
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)
const FTensor::Tensor2< T, Dim, Dim > Vec
@ MOFEM_DATA_INCONSISTENCY
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...