v0.14.0
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 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...
 

Detailed Description

Creating and scattering vectors on the mesh for given problem.

Function Documentation

◆ 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
pointerto porblem struture
RowColDatafor row or column (i.e. Row,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 form V vector.

Definition at line 301 of file VecManager.cpp.

303 {
305 typedef NumeredDofEntity_multiIndex::index<PetscGlobalIdx_mi_tag>::type
306 DofsByGlobalIdx;
307 DofsByGlobalIdx *dofs;
308 DofIdx nb_dofs;
309 switch (rc) {
310 case ROW:
311 nb_dofs = problem_ptr->getNbDofsRow();
312 dofs = const_cast<DofsByGlobalIdx *>(
313 &problem_ptr->numeredRowDofsPtr->get<PetscGlobalIdx_mi_tag>());
314 break;
315 case COL:
316 nb_dofs = problem_ptr->getNbDofsCol();
317 dofs = const_cast<DofsByGlobalIdx *>(
318 &problem_ptr->numeredColDofsPtr->get<PetscGlobalIdx_mi_tag>());
319 break;
320 default:
321 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
322 "Function works only for ROWs and COLs");
323 }
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: {
329 VecScatter ctx;
330 Vec V_glob;
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);
334 int size;
335 CHKERR VecGetSize(V_glob, &size);
336 if (size != nb_dofs) {
337 SETERRQ(comm, MOFEM_DATA_INCONSISTENCY,
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.");
341 }
342 PetscScalar *array;
343 CHKERR VecGetArray(V_glob, &array);
344 switch (mode) {
345 case INSERT_VALUES:
346 for (; miit != hi_miit; miit++)
347 (*miit)->getFieldData() = array[(*miit)->getPetscGlobalDofIdx()];
348 break;
349 case ADD_VALUES:
350 for (; miit != hi_miit; miit++)
351 (*miit)->getFieldData() += array[(*miit)->getPetscGlobalDofIdx()];
352 break;
353 default:
354 SETERRQ(comm, MOFEM_NOT_IMPLEMENTED, "not implemented");
355 }
356 CHKERR VecRestoreArray(V_glob, &array);
357 CHKERR VecScatterDestroy(&ctx);
358 CHKERR VecDestroy(&V_glob);
359 break;
360 }
361 default:
362 SETERRQ(comm, MOFEM_NOT_IMPLEMENTED, "not implemented");
363 }
365}
@ COL
Definition: definitions.h:123
@ ROW
Definition: definitions.h:123
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
const FTensor::Tensor2< T, Dim, Dim > Vec
int DofIdx
Index of DOF.
Definition: Types.hpp:18

◆ 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
nameof the problem
RowColDatafor row or column (i.e. Row,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 form V vector.

Definition at line 368 of file VecManager.cpp.

370 {
371 const MoFEM::Interface &m_field = cOre;
372 const Problem *problem_ptr;
374 CHKERR m_field.get_problem(name, &problem_ptr);
375 CHKERR setGlobalGhostVector(problem_ptr, rc, V, mode, scatter_mode);
377}
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
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:301
Deprecated interface functions.
const MoFEM::Interface & cOre
Definition: VecManager.hpp:28

◆ 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
pointerto problem struture
RowColDatafor row or column:e (i.e. Row,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 202 of file VecManager.cpp.

205 {
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 auto miit = dofs->lower_bound(0);
243 auto hi_miit = dofs->upper_bound(nb_local_dofs + nb_ghost_dofs);
244 DofIdx ii = 0;
245 switch (scatter_mode) {
246 case SCATTER_FORWARD: {
247 PetscScalar *array;
248 CHKERR VecGetArray(Vlocal, &array);
249 switch (mode) {
250 case INSERT_VALUES:
251 for (; miit != hi_miit; ++miit, ++ii)
252 array[ii] = (*miit)->getFieldData();
253 break;
254 case ADD_VALUES:
255 for (; miit != hi_miit; ++miit, ++ii)
256 array[ii] += (*miit)->getFieldData();
257 break;
258 default:
259 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
260 }
261 CHKERR VecRestoreArray(Vlocal, &array);
262 }; break;
263 case SCATTER_REVERSE: {
264 const PetscScalar *array;
265 CHKERR VecGetArrayRead(Vlocal, &array);
266 switch (mode) {
267 case INSERT_VALUES:
268 for (; miit != hi_miit; ++miit, ++ii) {
269 (*miit)->getFieldData() = array[ii];
270 }
271 break;
272 case ADD_VALUES:
273 for (; miit != hi_miit; ++miit, ++ii)
274 (*miit)->getFieldData() += array[ii];
275 break;
276 default:
277 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
278 }
279 CHKERR VecRestoreArrayRead(Vlocal, &array);
280 }; break;
281 default:
282 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
283 }
284 CHKERR VecGhostRestoreLocalForm(V, &Vlocal);
286}
NumeredDofEntity_multiIndex::index< PetscLocalIdx_mi_tag >::type NumeredDofEntityByLocalIdx
Numbered DoF multi-index by local index.

◆ 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
nameof the problem
RowColDatafor row or column:e (i.e. Row,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 288 of file VecManager.cpp.

291 {
292 const MoFEM::Interface &m_field = cOre;
293 const Problem *problem_ptr;
295 CHKERR m_field.get_problem(name, &problem_ptr);
296 CHKERR setLocalGhostVector(problem_ptr, rc, V, mode, scatter_mode);
298}
MoFEMErrorCode setLocalGhostVector(const Problem *problem_ptr, RowColData rc, Vec V, InsertMode mode, ScatterMode scatter_mode) const
set values of vector from/to mesh database
Definition: VecManager.cpp:202

◆ 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_ptrpointer to problem
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 form V vector.

Definition at line 544 of file VecManager.cpp.

547 {
548 const MoFEM::Interface &m_field = cOre;
549 auto fields_ptr = m_field.get_fields();
550 auto dofs_ptr = m_field.get_dofs();
553 DofIdx nb_dofs;
554 switch (rc) {
555 case ROW:
556 nb_dofs = problem_ptr->getNbDofsRow();
557 dofs = const_cast<NumeredDofEntityByUId *>(
558 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
559 break;
560 case COL:
561 nb_dofs = problem_ptr->getNbDofsCol();
562 dofs = const_cast<NumeredDofEntityByUId *>(
563 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
564 break;
565 default:
566 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "not implemented");
567 }
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()) {
570 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
571 "cpy field < %s > not found, (top tip: check spelling)",
572 cpy_field_name.c_str());
573 }
574 const auto cpy_bit_number = (*cpy_fit)->getBitNumber();
575
576 auto miit = dofs->lower_bound(
578 if (miit == dofs->end()) {
579 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
580 "problem field < %s > not found, (top tip: check spelling)",
581 field_name.c_str());
582 }
583 auto hi_miit = dofs->upper_bound(
585 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
586 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
587 "fields have to have same space");
588 }
589 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
590 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
591 "fields have to have same rank");
592 }
593 switch (scatter_mode) {
594 case SCATTER_REVERSE: {
595 Vec V_glob;
596 VecScatter ctx;
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);
600 int size;
601 CHKERR VecGetSize(V_glob, &size);
602 if (size != nb_dofs)
603 SETERRQ(
604 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
605 "data inconsistency: nb. of dofs and declared nb. dofs in database");
606 PetscScalar *array;
607 CHKERR VecGetArray(V_glob, &array);
608
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);
615 else
616 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
617 "Operation mode not implemented");
618
619 CHKERR VecRestoreArray(V_glob, &array);
620 CHKERR VecDestroy(&V_glob);
621 CHKERR VecScatterDestroy(&ctx);
622 } break;
623 case SCATTER_FORWARD: {
624 for (; miit != hi_miit; miit++) {
625 const auto uid = DofEntity::getUniqueIdCalculate(
626 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
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()) {
630 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
631 "No data to fill the vector (top tip: you want scatter forward "
632 "or scatter reverse?)");
633 }
634 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
635 (*diiiit)->getFieldData(), mode);
636 }
637 CHKERR VecAssemblyBegin(V);
638 CHKERR VecAssemblyEnd(V);
639 } break;
640 default:
641 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
642 }
644}
@ 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

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
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 form V vector.

Definition at line 646 of file VecManager.cpp.

649 {
650 const MoFEM::Interface &m_field = cOre;
651 const Problem *problem_ptr;
653 CHKERR m_field.get_problem(name, &problem_ptr);
654 CHKERR setOtherGlobalGhostVector(problem_ptr, field_name, cpy_field_name, rc,
655 V, mode, scatter_mode);
657}
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:544

◆ 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
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 form V vector.

Definition at line 413 of file VecManager.cpp.

416 {
417 const MoFEM::Interface &m_field = cOre;
418 auto fields_ptr = m_field.get_fields();
419 auto dofs_ptr = m_field.get_dofs();
421 using DofsByUId = NumeredDofEntity_multiIndex::index<Unique_mi_tag>::type;
422 DofsByUId *dofs;
423
424 switch (rc) {
425 case ROW:
426 dofs = const_cast<DofsByUId *>(
427 &problem_ptr->numeredRowDofsPtr->get<Unique_mi_tag>());
428 break;
429 case COL:
430 dofs = const_cast<DofsByUId *>(
431 &problem_ptr->numeredColDofsPtr->get<Unique_mi_tag>());
432 break;
433 default:
434 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
435 }
436
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())
439 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
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();
443
444 auto miit = dofs->lower_bound(
446 if (miit == dofs->end()) {
447 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
448 "cpy field < %s > not found, (top tip: check spelling)",
449 field_name.c_str());
450 }
451 auto hi_miit = dofs->upper_bound(
453
454 if ((*miit)->getSpace() != (*cpy_fit)->getSpace()) {
455 SETERRQ4(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
456 "fields have to have same space (%s) %s != (%s) %s",
457 (*miit)->getName().c_str(), FieldSpaceNames[(*miit)->getSpace()],
458 cpy_field_name.c_str(), FieldSpaceNames[(*cpy_fit)->getSpace()]);
459 }
460 if ((*miit)->getNbOfCoeffs() != (*cpy_fit)->getNbOfCoeffs()) {
461 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
462 "fields have to have same rank");
463 }
464
465 switch (scatter_mode) {
466 case SCATTER_REVERSE: {
467
468 PetscScalar *array;
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);
476 else
477 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
478 "Operation mode not implemented");
479 CHKERR VecRestoreArray(V, &array);
480 } break;
481 case SCATTER_FORWARD: {
482 for (; miit != hi_miit; miit++) {
483 if (!miit->get()->getHasLocalIndex())
484 continue;
485 const auto uid = DofEntity::getUniqueIdCalculate(
486 (*miit)->getEntDofIdx(), FieldEntity::getLocalUniqueIdCalculate(
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()) {
490 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
491 "no data to fill the vector (top tip: you want scatter forward "
492 "or scatter reverse?)");
493 }
494 CHKERR VecSetValue(V, (*miit)->getPetscGlobalDofIdx(),
495 (*diiiit)->getFieldData(), mode);
496 }
497 CHKERR VecAssemblyBegin(V);
498 CHKERR VecAssemblyEnd(V);
499 } break;
500 default:
501 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "not implemented");
502 }
504}
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

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 form V vector.

Definition at line 506 of file VecManager.cpp.

509 {
510 const MoFEM::Interface &m_field = cOre;
511 const Problem *problem_ptr;
513 CHKERR m_field.get_problem(name, &problem_ptr);
514 CHKERR setOtherLocalGhostVector(problem_ptr, field_name, cpy_field_name, rc,
515 V, mode, scatter_mode);
517}
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:413

◆ vecCreateGhost()

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

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 64 of file VecManager.cpp.

65 {
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}
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
Definition: LogManager.hpp:296
#define MOFEM_LOG_FUNCTION()
Set scope.
Definition: LogManager.hpp:325
virtual MPI_Comm & get_comm() const =0

◆ vecCreateSeq()

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

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 38 of file VecManager.cpp.

39 {
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}

◆ 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
xinvector
x_probleproblem name
x_fieldname
yinvector
y_problemproblem name
y_field_name
Return values
newctxscatter

Definition at line 131 of file VecManager.cpp.

134 {
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}
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce 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

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

Parameters
xinvector
x_probleproblem name
yinvector
y_problemproblem name
Return values
newctxscatter

Definition at line 153 of file VecManager.cpp.

155 {
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}