OGS
|
Wrapper class for PETSc vector.
It can be used to create a global vector for either parallel or serial computing.
Caution: Using it to create a local vector is not allowed, as the created vector will be partitioned and distributed across all ranks in an MPI environment.
Definition at line 39 of file PETScVector.h.
#include <PETScVector.h>
Public Types | |
using | IndexType = PetscInt |
using | PETSc_Vec = Vec |
Public Member Functions | |
PETScVector () | |
PETScVector (const PetscInt vec_size, const bool is_global_size=true) | |
Constructor. | |
PETScVector (const PetscInt vec_size, const std::vector< PetscInt > &ghost_ids, const bool is_global_size=true) | |
Constructor. | |
PETScVector (const PETScVector &existing_vec, const bool deep_copy=true) | |
Copy constructor. | |
PETScVector (PETScVector &&other) | |
~PETScVector () | |
void | finalizeAssembly () |
Perform MPI collection of assembled entries in buffer. | |
PetscInt | size () const |
Get the global size of the vector. | |
PetscInt | getLocalSize () const |
Get the number of entries in the same rank. | |
PetscInt | getGhostSize () const |
Get the number of ghost entries in the same rank. | |
PetscInt | getRangeBegin () const |
Get the start index of the local vector. | |
PetscInt | getRangeEnd () const |
Get the end index of the local vector. | |
void | set (const PetscInt i, const PetscScalar value) |
void | add (const PetscInt i, const PetscScalar value) |
template<class T_SUBVEC > | |
void | add (const std::vector< PetscInt > &e_idxs, const T_SUBVEC &sub_vec) |
template<class T_SUBVEC > | |
void | set (const std::vector< PetscInt > &e_idxs, const T_SUBVEC &sub_vec) |
void | setZero () |
PETScVector & | operator= (PETScVector &&)=delete |
void | setLocalAccessibleVector () const |
std::vector< PetscScalar > | get (std::vector< IndexType > const &indices) const |
PetscScalar | operator[] (PetscInt idx) const |
void | getGlobalVector (std::vector< PetscScalar > &u) const |
PetscScalar | get (const PetscInt idx) const |
PETSc_Vec & | getRawVector () |
Exposes the underlying PETSc vector. | |
PETSc_Vec const & | getRawVector () const |
void | copyValues (std::vector< PetscScalar > &u) const |
void | viewer (const std::string &file_name, const PetscViewerFormat vw_format=PETSC_VIEWER_ASCII_MATLAB) const |
void | shallowCopy (const PETScVector &v) |
Private Member Functions | |
void | destroy () |
void | gatherLocalVectors (PetscScalar local_array[], PetscScalar global_array[]) const |
Collect local vectors. | |
PetscScalar * | getLocalVector () const |
PetscInt | getLocalIndex (const PetscInt global_index) const |
Get local index by a global index. | |
void | restoreArray (PetscScalar *array) const |
void | config () |
A function called by constructors to configure members. | |
Private Attributes | |
PETSc_Vec | v_ = nullptr |
PETSc_Vec | v_loc_ = nullptr |
PetscInt | start_rank_ |
Starting index in a rank. | |
PetscInt | end_rank_ |
Ending index in a rank. | |
PetscInt | size_ |
Size of the vector. | |
PetscInt | size_loc_ |
Size of local entries. | |
PetscInt | size_ghosts_ = 0 |
Size of local ghost entries. | |
bool | created_with_ghost_id_ = false |
Flag to indicate whether the vector is created with ghost entry indices. | |
std::vector< PetscScalar > | entry_array_ |
Array containing the entries of the vector. If the vector is created without given ghost IDs, the array contains all entries of the global vector, v_. Otherwise it only contains the entries of the global vector owned by the current rank. | |
std::map< PetscInt, PetscInt > | global_ids2local_ids_ghost_ |
Map global indices of ghost entries to local indices. | |
Friends | |
void | finalizeVectorAssembly (PETScVector &vec) |
Function to finalize the vector assembly. | |
using MathLib::PETScVector::IndexType = PetscInt |
Definition at line 42 of file PETScVector.h.
using MathLib::PETScVector::PETSc_Vec = Vec |
Definition at line 45 of file PETScVector.h.
|
inline |
Definition at line 48 of file PETScVector.h.
MathLib::PETScVector::PETScVector | ( | const PetscInt | vec_size, |
const bool | is_global_size = true ) |
Constructor.
vec_size | The size of the vector, either global or local |
is_global_size | The flag of the type of vec_size, i.e. whether it is a global size or local size. The default is true. If is_global_size is true, the vector is created by the global size, the local size of the vector is determined by PETSc, and vice versa is the same. |
Definition at line 30 of file PETScVector.cpp.
MathLib::PETScVector::PETScVector | ( | const PetscInt | vec_size, |
const std::vector< PetscInt > & | ghost_ids, | ||
const bool | is_global_size = true ) |
Constructor.
vec_size | The size of the vector, either global or local |
ghost_ids | The global indices of ghost entries |
is_global_size | The flag of the type of vec_size, i.e. whether it is a global size or local size. The default is true. |
Definition at line 47 of file PETScVector.cpp.
References config(), global_ids2local_ids_ghost_, size_loc_, and v_.
|
explicit |
Copy constructor.
existing_vec | The vector to be copied |
deep_copy | The flag for a deep copy, which means to copy the values as well, the default is true |
Definition at line 75 of file PETScVector.cpp.
References shallowCopy(), and v_.
MathLib::PETScVector::PETScVector | ( | PETScVector && | other | ) |
Definition at line 86 of file PETScVector.cpp.
|
inline |
|
inline |
Add a value to an entry.
i | Number of the entry |
value | Value. |
Definition at line 114 of file PETScVector.h.
References v_.
|
inline |
Add values to several entries.
e_idxs | Indices of entries to be added Note: std::size_t cannot be the type of e_idxs template argument. |
sub_vec | Entries to be added. |
Definition at line 127 of file PETScVector.h.
References v_.
|
private |
A function called by constructors to configure members.
Definition at line 99 of file PETScVector.cpp.
References end_rank_, size_, size_loc_, start_rank_, and v_.
Referenced by PETScVector(), PETScVector(), and shallowCopy().
void MathLib::PETScVector::copyValues | ( | std::vector< PetscScalar > & | u | ) | const |
Copy local entries including ghost ones to a vector.
u | a vector for the values of local entries. It will be resized to hold the current vector data. |
Definition at line 192 of file PETScVector.cpp.
References getGhostSize(), getLocalSize(), getLocalVector(), restoreArray(), and MathLib::u.
Referenced by setLocalAccessibleVector().
|
inlineprivate |
Definition at line 229 of file PETScVector.h.
References v_.
Referenced by ~PETScVector(), and shallowCopy().
void MathLib::PETScVector::finalizeAssembly | ( | ) |
Perform MPI collection of assembled entries in buffer.
Definition at line 110 of file PETScVector.cpp.
References v_.
Referenced by MathLib::applyKnownSolution(), and MathLib::LinAlg::finalizeAssembly().
|
private |
Collect local vectors.
local_array | Local array |
global_array | Global array |
Definition at line 116 of file PETScVector.cpp.
References size_loc_.
Referenced by getGlobalVector().
PetscScalar MathLib::PETScVector::get | ( | const PetscInt | idx | ) | const |
Definition at line 201 of file PETScVector.cpp.
References created_with_ghost_id_, entry_array_, and getLocalIndex().
std::vector< PetscScalar > MathLib::PETScVector::get | ( | std::vector< IndexType > const & | indices | ) | const |
Get several entries. setLocalAccessibleVector() must be called beforehand.
Definition at line 211 of file PETScVector.cpp.
Referenced by operator[]().
|
inline |
Get the number of ghost entries in the same rank.
Definition at line 93 of file PETScVector.h.
References size_ghosts_.
Referenced by copyValues().
void MathLib::PETScVector::getGlobalVector | ( | std::vector< PetscScalar > & | u | ) | const |
Get global vector
u | Array to store the global vector. Memory allocation is needed in advance |
Definition at line 143 of file PETScVector.cpp.
References gatherLocalVectors(), OGS_FATAL, size_, MathLib::u, and v_.
Referenced by setLocalAccessibleVector().
|
private |
Get local index by a global index.
Definition at line 238 of file PETScVector.cpp.
References end_rank_, global_ids2local_ids_ghost_, OGS_FATAL, size_, and start_rank_.
Referenced by get().
|
inline |
Get the number of entries in the same rank.
Definition at line 91 of file PETScVector.h.
References size_loc_.
Referenced by copyValues().
|
private |
Get local vector, i.e. entries in the same rank
Definition at line 220 of file PETScVector.cpp.
References created_with_ghost_id_, global_ids2local_ids_ghost_, v_, and v_loc_.
Referenced by copyValues().
|
inline |
Get the start index of the local vector.
Definition at line 95 of file PETScVector.h.
References start_rank_.
|
inline |
Get the end index of the local vector.
Definition at line 97 of file PETScVector.h.
References end_rank_.
|
inline |
Exposes the underlying PETSc vector.
Definition at line 178 of file PETScVector.h.
References v_.
Referenced by MathLib::LinAlg::axpby(), MathLib::LinAlg::axpy(), MathLib::LinAlg::aypx(), MathLib::LinAlg::componentwiseDivide(), MathLib::LinAlg::copy(), MathLib::LinAlg::matMult(), MathLib::LinAlg::matMultAdd(), MathLib::LinAlg::norm1(), MathLib::LinAlg::norm2(), MathLib::LinAlg::normMax(), MathLib::LinAlg::scale(), MathLib::LinAlg::set(), and MathLib::PETScLinearSolver::solve().
|
inline |
Exposes the underlying PETSc vector.
Definition at line 185 of file PETScVector.h.
References v_.
|
delete |
Disallow moving.
|
inline |
Get the value of an entry by [] operator. setLocalAccessibleVector() must be called beforehand.
Definition at line 162 of file PETScVector.h.
References get().
|
inlineprivate |
Restore array after finish access local array
array | Pointer to the local array fetched by VecGetArray |
Definition at line 271 of file PETScVector.cpp.
References created_with_ghost_id_, global_ids2local_ids_ghost_, v_, and v_loc_.
Referenced by copyValues().
|
inline |
Insert a single entry with value.
i | Entry index |
value | Entry value |
Definition at line 104 of file PETScVector.h.
References v_.
Referenced by MathLib::applyKnownSolution().
|
inline |
Add values to several entries
e_idxs | Indices of entries to be added. Note: std::size_t cannot be the type of e_idxs template argument |
sub_vec | Entries to be added |
Definition at line 140 of file PETScVector.h.
References v_.
void MathLib::PETScVector::setLocalAccessibleVector | ( | ) | const |
Set local accessible vector in order to get entries. Call this before call operator[] or get(...).
Definition at line 180 of file PETScVector.cpp.
References copyValues(), created_with_ghost_id_, entry_array_, getGlobalVector(), and size_.
Referenced by MathLib::LinAlg::setLocalAccessibleVector().
|
inline |
void MathLib::PETScVector::shallowCopy | ( | const PETScVector & | v | ) |
Definition at line 302 of file PETScVector.cpp.
References config(), created_with_ghost_id_, destroy(), end_rank_, global_ids2local_ids_ghost_, size_, size_ghosts_, size_loc_, start_rank_, MathLib::v, and v_.
Referenced by PETScVector(), MathLib::LinAlg::copy(), MathLib::LinAlg::matMult(), and MathLib::LinAlg::matMultAdd().
|
inline |
void MathLib::PETScVector::viewer | ( | const std::string & | file_name, |
const PetscViewerFormat | vw_format = PETSC_VIEWER_ASCII_MATLAB ) const |
View the global vector for test purpose. Do not use it for output a big vector.
file_name | File name for output |
vw_format | File format listed as: PETSC_VIEWER_DEFAULT Default format PETSC_VIEWER_ASCII_MATLAB MATLAB format PETSC_VIEWER_ASCII_DENSE Print matrix as dense PETSC_VIEWER_ASCII_IMPL Implementation-specific format (which is in many cases the same as the default) PETSC_VIEWER_ASCII_INFO Basic information about object PETSC_VIEWER_ASCII_INFO_DETAIL More detailed info about object PETSC_VIEWER_ASCII_COMMON Identical output format for all objects of a particular type PETSC_VIEWER_ASCII_INDEX (for vectors) Prints the vector element number next to each vector entry PETSC_VIEWER_ASCII_SYMMODU Print parallel vectors without indicating the processor ranges PETSC_VIEWER_ASCII_VTK Outputs the object to a VTK file PETSC_VIEWER_NATIVE Store the object to the binary file in its native format (for example, dense matrices are stored as dense), DMDA vectors are dumped directly to the file instead of being first put in the natural ordering PETSC_VIEWER_DRAW_BASIC Views the vector with a simple 1d plot PETSC_VIEWER_DRAW_LG Views the vector with a line graph PETSC_VIEWER_DRAW_CONTOUR Views the vector with a contour plot |
Definition at line 284 of file PETScVector.cpp.
Referenced by viewer().
|
friend |
Function to finalize the vector assembly.
Definition at line 319 of file PETScVector.cpp.
|
private |
Flag to indicate whether the vector is created with ghost entry indices.
Definition at line 256 of file PETScVector.h.
Referenced by get(), getLocalVector(), restoreArray(), setLocalAccessibleVector(), and shallowCopy().
|
private |
Ending index in a rank.
Definition at line 246 of file PETScVector.h.
Referenced by config(), getLocalIndex(), getRangeEnd(), and shallowCopy().
|
mutableprivate |
Array containing the entries of the vector. If the vector is created without given ghost IDs, the array contains all entries of the global vector, v_. Otherwise it only contains the entries of the global vector owned by the current rank.
Definition at line 265 of file PETScVector.h.
Referenced by get(), and setLocalAccessibleVector().
|
mutableprivate |
Map global indices of ghost entries to local indices.
Definition at line 268 of file PETScVector.h.
Referenced by PETScVector(), getLocalIndex(), getLocalVector(), restoreArray(), and shallowCopy().
|
private |
Size of the vector.
Definition at line 249 of file PETScVector.h.
Referenced by config(), getGlobalVector(), getLocalIndex(), setLocalAccessibleVector(), shallowCopy(), and size().
|
private |
Size of local ghost entries.
Definition at line 253 of file PETScVector.h.
Referenced by getGhostSize(), and shallowCopy().
|
private |
Size of local entries.
Definition at line 251 of file PETScVector.h.
Referenced by PETScVector(), config(), gatherLocalVectors(), getLocalSize(), and shallowCopy().
|
private |
Starting index in a rank.
Definition at line 244 of file PETScVector.h.
Referenced by config(), getLocalIndex(), getRangeBegin(), and shallowCopy().
|
private |
Definition at line 238 of file PETScVector.h.
Referenced by PETScVector(), PETScVector(), PETScVector(), add(), add(), config(), destroy(), finalizeAssembly(), getGlobalVector(), getLocalVector(), getRawVector(), getRawVector(), restoreArray(), set(), set(), setZero(), shallowCopy(), and viewer().
|
mutableprivate |
Local vector, which is only for the case that v_ is created with ghost entries.
Definition at line 241 of file PETScVector.h.
Referenced by getLocalVector(), and restoreArray().