The SplitSystemSolver class allows to solve a split linear system. More...


Public Types | |
typedef SplitSystemMatrixAssembler | Assembler |
typedef SystemSplitStrategy | splitStrategy |
![]() | |
typedef SystemMatrixAssembler | Assembler |
enum | FileFormat { FILE_BINARY , FILE_ASCII } |
Public Member Functions | |
SplitSystemSolver (bool debug=false) | |
SplitSystemSolver (bool flatten, bool transpose, bool debug) | |
SplitSystemSolver (bool transpose, bool debug) | |
SplitSystemSolver (const std::string &prefix, bool debug=false) | |
SplitSystemSolver (const std::string &prefix, bool flatten, bool transpose, bool debug) | |
SplitSystemSolver (const std::string &prefix, bool transpose, bool debug) | |
void | assembly (const Assembler &assembler) |
void | assembly (const Assembler &assembler, const SystemMatrixOrdering &reordering) |
void | assembly (const SparseMatrix &matrix, SystemSplitStrategy splitStrategy, const std::vector< int > &splitSizes) |
void | assembly (const SparseMatrix &matrix, SystemSplitStrategy splitStrategy, const std::vector< int > &splitSizes, const SystemMatrixOrdering &reordering) |
void | exportMatrix (const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const override |
int | getBlockSize () const override |
int | getSplitCount () const |
KSPOptions & | getSplitKSPOptions (int split) |
const KSPOptions & | getSplitKSPOptions (int split) const |
const KSPStatus & | getSplitKSPStatus (int split) const |
std::vector< int > | getSplitOffsets () const |
std::vector< int > | getSplitSizes () const |
SystemSplitStrategy | getSplitStrategy () const |
SystemSolver (bool debug=false) | |
SystemSolver (bool flatten, bool transpose, bool debug) | |
SystemSolver (bool transpose, bool debug) | |
SystemSolver (const std::string &prefix, bool debug=false) | |
SystemSolver (const std::string &prefix, bool flatten, bool transpose, bool debug) | |
SystemSolver (const std::string &prefix, bool transpose, bool debug) | |
void | update (const Assembler &assembler) |
void | update (const SparseMatrix &elements) |
void | update (long nRows, const long *rows, const Assembler &assembler) |
void | update (long nRows, const long *rows, const SparseMatrix &elements) |
![]() | |
SystemSolver (bool debug=false) | |
SystemSolver (bool flatten, bool transpose, bool debug) | |
SystemSolver (bool transpose, bool debug) | |
SystemSolver (const std::string &prefix, bool debug=false) | |
SystemSolver (const std::string &prefix, bool flatten, bool transpose, bool debug) | |
SystemSolver (const std::string &prefix, bool transpose, bool debug) | |
virtual | ~SystemSolver () |
void | assembly (const Assembler &assembler) |
void | assembly (const Assembler &assembler, const SystemMatrixOrdering &reordering) |
void | assembly (const SparseMatrix &matrix) |
void | assembly (const SparseMatrix &matrix, const SystemMatrixOrdering &reordering) |
virtual void | clear () |
virtual void | clearWorkspace () |
void | dumpSystem (const std::string &header, const std::string &directory, const std::string &prefix="") const |
void | enableForceConsistency (bool enable) |
virtual void | exportRHS (const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const |
virtual void | exportSolution (const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const |
long | getColCount () const |
long | getColElementCount () const |
long | getColGlobalCount () const |
long | getColGlobalElementCount () const |
KSPOptions & | getKSPOptions () |
const KSPOptions & | getKSPOptions () const |
const KSPStatus & | getKSPStatus () const |
double * | getRHSRawPtr () |
const double * | getRHSRawPtr () const |
const double * | getRHSRawReadPtr () const |
long | getRowCount () const |
long | getRowElementCount () const |
long | getRowGlobalCount () const |
long | getRowGlobalElementCount () const |
double * | getSolutionRawPtr () |
const double * | getSolutionRawPtr () const |
const double * | getSolutionRawReadPtr () const |
bool | getTranspose () const |
virtual void | importMatrix (const std::string &filePath) |
virtual void | importRHS (const std::string &filePath) |
virtual void | importSolution (const std::string &filePath) |
bool | isAssembled () const |
bool | isForceConsistencyEnabled () const |
bool | isPartitioned () const |
void | restoreRHSRawPtr (double *raw_rhs) |
void | restoreRHSRawReadPtr (const double *raw_rhs) const |
void | restoreSolutionRawPtr (double *raw_solution) |
void | restoreSolutionRawReadPtr (const double *raw_solution) const |
void | restoreSystem (MPI_Comm communicator, bool redistribute, const std::string &directory, const std::string &prefix="") |
virtual void | setNullSpace () |
void | setTranspose (bool transpose) |
void | solve () |
void | solve (const std::vector< double > &rhs, std::vector< double > *solution) |
void | unsetNullSpace () |
void | update (const Assembler &assembler) |
void | update (const SparseMatrix &elements) |
void | update (long nRows, const long *rows, const Assembler &assembler) |
void | update (long nRows, const long *rows, const SparseMatrix &elements) |
Protected Member Functions | |
void | assembly (const Assembler &assembler) |
void | assembly (const Assembler &assembler, const SystemMatrixOrdering &reordering) |
void | assembly (const SparseMatrix &matrix) |
void | assembly (const SparseMatrix &matrix, const SystemMatrixOrdering &reordering) |
template<typename DerivedSystemSolver> | |
void | assembly (const typename DerivedSystemSolver::Assembler &assembler, const SystemMatrixOrdering &reordering) |
void | destroyKSPOptions () override |
void | destroyKSPStatus () override |
virtual void | destroySplitKSPOptions () |
virtual void | destroySplitKSPStatuses () |
void | dumpInfo (std::ostream &systemStream) const override |
void | exportMatrix (Mat matrix, const std::string &filePath, FileFormat fileFormat) const |
void | fillKSPStatus () override |
virtual void | fillKSPStatus (KSP ksp, KSPStatus *status) const |
virtual void | fillSplitKSPStatuses () |
void | generateSplitIndexes (int split, long nItems, std::vector< std::size_t > *indexes) const |
std::string | generateSplitPath (const std::string &path, const std::string &index) const |
std::string | generateSplitPath (const std::string &path, int i) const |
std::string | generateSplitPath (const std::string &path, int i, int j) const |
void | generateSplitPermutation (long nItems, long itemGlobalOffset, IS *splitReordering) const |
void | initializeKSPOptions () override |
void | initializeKSPStatus () override |
virtual void | initializeSplitKSPOptions () |
virtual void | initializeSplitKSPStatuses () |
void | matrixAssembly (const Assembler &assembler) |
void | matrixDestroy () override |
void | matrixDump (std::ostream &systemStream, const std::string &directory, const std::string &prefix) const override |
void | matrixFill (const std::string &filePath) override |
void | matrixRestore (std::istream &systemStream, const std::string &directory, const std::string &prefix, bool redistribute) override |
void | matrixUpdate (long nRows, const long *rows, const Assembler &assembler) |
void | postKSPSolveActions () override |
void | preKrylovSetupActions () override |
virtual void | resetKSPOptions (KSPOptions *options) const |
void | resetKSPStatus () override |
virtual void | resetKSPStatus (KSPStatus *status) const |
virtual void | resetSplitKSPStatuses () |
void | restoreInfo (std::istream &systemStream) override |
void | setupKrylov () override |
virtual void | setupKrylov (KSP ksp, const KSPOptions &options) const |
void | setupPreconditioner () override |
virtual void | setupPreconditioner (PC pc, const KSPOptions &options) const |
virtual void | setupSplitKrylovs () |
virtual void | setupSplitPreconditioners () |
void | update (const Assembler &assembler) |
void | update (long nRows, const long *rows, const Assembler &assembler) |
template<typename DerivedSystemSolver> | |
void | update (long nRows, const long *rows, const typename DerivedSystemSolver::Assembler &assembler) |
void | vectorsCreate () override |
void | vectorsCreateSplitPermutations () |
void | vectorsDestroy () override |
void | vectorsReorder (bool invert) override |
void | vectorsRestore (std::istream &systemStream, const std::string &directory, const std::string &prefix) override |
![]() | |
template<typename DerivedSystemSolver> | |
void | assembly (const typename DerivedSystemSolver::Assembler &assembler, const SystemMatrixOrdering &reordering) |
void | clearReordering () |
virtual void | createKSP () |
void | createMatrix (int rowBlockSize, int colBlockSize, int nNestRows, int nNestCols, Mat *subMatrices, Mat *matrix) const |
void | createMatrix (int rowBlockSize, int colBlockSize, Mat *matrix) const |
void | createVector (int blockSize, int nestSize, Vec *subVectors, Vec *vector) const |
void | createVector (int blockSize, Vec *vector) const |
virtual void | destroyKSP () |
void | destroyMatrix (Mat *matrix) const |
void | destroyVector (Vec *vector) const |
void | dumpMatrix (Mat matrix, const std::string &directory, const std::string &name) const |
void | dumpVector (Vec vector, const std::string &directory, const std::string &name) const |
void | dumpVector (Vec vector, VectorSide side, const std::string &directory, const std::string &name) const |
void | exportMatrix (Mat matrix, const std::string &filePath, FileFormat fileFormat) const |
void | exportVector (Vec vector, const std::string &filePath, FileFormat fileFormat) const |
void | exportVector (Vec vector, std::vector< double > *data) const |
void | fillMatrix (Mat matrix, const std::string &filePath) const |
void | fillVector (Vec vector, const std::string &filePath) const |
void | fillVector (Vec vector, const std::vector< double > &data) const |
virtual void | finalizeKSP () |
const MPI_Comm & | getCommunicator () const |
std::string | getDataFilePath (const std::string &directory, const std::string &name) const |
virtual int | getDumpVersion () const |
std::string | getFilePath (const std::string &directory, const std::string &name) const |
std::string | getInfoFilePath (const std::string &directory, const std::string &name) const |
void | matrixAssembly (const Assembler &assembler) |
void | matrixUpdate (long nRows, const long *rows, const Assembler &assembler) |
virtual void | postKrylovSetupActions () |
virtual void | postPreconditionerSetupActions () |
virtual void | preKSPSolveActions () |
virtual void | prepareKSP () |
virtual void | prePreconditionerSetupActions () |
void | reorderVector (Vec vector, IS permutations, bool invert) const |
void | restoreMatrix (const std::string &directory, const std::string &name, bool redistribute, Mat *matrix) const |
void | restoreVector (const std::string &directory, const std::string &name, bool redistribute, Vec *vector) const |
void | restoreVector (const std::string &directory, const std::string &name, Mat matrix, VectorSide side, Vec *vector) const |
void | restoreVector (const std::string &directory, const std::string &name, std::size_t localSize, std::size_t globalSize, Vec *vector) const |
void | setReordering (long nRows, long nCols, const SystemMatrixOrdering &reordering) |
virtual void | solveKSP () |
template<typename DerivedSystemSolver> | |
void | update (long nRows, const long *rows, const typename DerivedSystemSolver::Assembler &assembler) |
virtual void | vectorsDump (std::ostream &systemStream, const std::string &directory, const std::string &prefix) const |
virtual void | vectorsFill (const std::string &rhsFilePath, const std::string &solutionFilePath) |
virtual void | vectorsFill (const std::vector< double > &rhs, const std::vector< double > &solution) |
Protected Attributes | |
std::vector< Mat > | m_splitAs |
std::vector< KSPOptions > | m_splitKSPOptions |
std::vector< KSPStatus > | m_splitKSPStatuses |
SystemSplitStrategy | m_splitStrategy |
![]() | |
Mat | m_A |
IS | m_colReordering |
bool | m_convergenceMonitorEnabled |
bool | m_flatten |
KSP | m_KSP |
bool | m_KSPDirty |
KSPOptions | m_KSPOptions |
KSPStatus | m_KSPStatus |
Vec | m_rhs |
IS | m_rowReordering |
Vec | m_solution |
bool | m_transpose |
Friends | |
class | SystemSolver |
Additional Inherited Members | |
![]() | |
enum | VectorSide { VECTOR_SIDE_RIGHT , VECTOR_SIDE_LEFT } |
Detailed Description
The SplitSystemSolver class allows to solve a split linear system.
Block matrices represent an important class of problems in numerical linear algebra and offer the possibility of far more efficient iterative solvers than just treating the entire matrix as black box. Following the common linear algebra definition, a block matrix is a matrix divided in a small, problem-size independent (two, three or so) number of very large blocks. These blocks arise naturally from the underlying physics or discretization of the problem, for example, they may be associated with different variables of a physical problem. Under a certain numbering of unknowns the matrix can be written as
[ A00 A01 A02 ] [ A10 A11 A12 ] [ A20 A21 A22 ]
where each A_ij is an entire block (see the paragraph "Solving Block Matrices" in the PETSc manual).
When assembling the matrix, a monolithic matrix should be provided. For example, assuming to group the elements of the matrix in five-by-five groups (here, five is the so-called block size of the matrix and usually rises when coupling variables with different meaning, for example pressure, the three components of the velocity and temperature) the assembler will provide the following system:
[ | | | ] [ ] [ ] [ (5 x 5) | (5 x 5) | ... | (5 x 5) ] [ (5 x 5) ] [ (5 x 5) ] [ | | | ] [ ] [ ] [----------------------------------------------] [-----------] [-----------] [ | | | ] [ ] [ ] [ ... | ... | ... | (5 x 5) ] [ ... ] = [ ... ] [ | | | ] [ ] [ ] [----------------------------------------------] [-----------] [-----------] [ | | | ] [ ] [ ] [ (5 x 5) | (5 x 5) | ... | (5 x 5) ] [ (5 x 5) ] [ (5 x 5) ] [ | | | ] [ ] [ ]
Internally, the monolithic matrix will be split into blocks. For example, considering two splits, the first one that group together the first four variables and the second one that holds the last variable (i.e., split sizes equal to [4, 1]), the internal split system will be:
[ | ] [ ] [ ] [ (4 * N) x (4 * M) | (4 * N) x (1 * M) ] [ (4 x M) ] [ (4 x N) ] [ | ] [ ] [ ] [-------------------------------------------] [-----------] = [-----------] [ | ] [ ] [ ] [ (1 * N) x (4 * M) | (1 * N) x (1 * M) ] [ (1 x M) ] [ (1 x N) ] [ | ] [ ] [ ]
where M and N are the number of rows and columns respectively.
The PETSc PCFIELDSPLIT preconditioner is used to solve the split system. There are different split strategies available:
1) DIAGONAL: this strategy assumes that the only non-zero blocks are the diagonal ones.
Considering a two-by-two block block matrix
[ A00 0 ] [ 0 A11 ],
the preconditioned problem will look like
[ KSPSolve(A00) 0 ] [ A00 0 ] [ 0 KSPSolve(A00) ] [ 0 A11 ],
in other words the preconditioner is:
( [ A00 0 ] ) approximate inverse ( ) ( [ 0 A11 ] ).
The system is solved efficiently by solving each block independently from the others.
Definition at line 66 of file system_solvers_split.hpp.
Member Typedef Documentation
◆ Assembler
Definition at line 69 of file system_solvers_split.hpp.
◆ splitStrategy
typedef SystemSplitStrategy bitpit::SplitSystemSolver::splitStrategy |
Definition at line 71 of file system_solvers_split.hpp.
Constructor & Destructor Documentation
◆ SplitSystemSolver() [1/6]
bitpit::SplitSystemSolver::SplitSystemSolver | ( | bool | debug = false | ) |
Blocks are solved using a flexible GMRES iterative method. If the system is partitioned each block is preconditioned using the (restricted) additive Schwarz method (ASM). On each block of the ASM preconditioner an incomplete LU factorization (ILU) is used. There is one ASM block per process. If the system is not partitioned it is preconditioned using the incomplete LU factorization (ILU).
2) LOWER: this strategy assumes that the only non-zero blocks are the ones on an below the diagonal.
Considering a two-by-two block block matrix
[ A00 0 ] [ A01 A11 ],
the preconditioner is
( [ A00 0 ] ) approximate inverse ( ) ( [ A01 A11 ] ).
The system is solved efficiently by first solving with A00, then applying A01 to that solution, removing it from the right hand side of the second block and then solving with A11, in other words
[ I 0 ] [ I 0 ] [ A00^-1 0 ] [ 0 A11^-1 ] [ - A10 I ] [ 0 I ].
This strategy can be seen as a "block" Gauss-Seidel with the blocks being the splits.
Blocks are solved using a flexible GMRES iterative method. If the system is partitioned each block is preconditioned using the (restricted) additive Schwarz method (ASM). On each block of the ASM preconditioner an incomplete LU factorization (ILU) is used. There is one ASM block per process. If the system is not partitioned it is preconditioned using the incomplete LU factorization (ILU).
3) FULL: this strategy doesn't make assumptions on the structure of the blocks.
Considering a two-by-two block block matrix
[ A00 A01 ] [ A11 A11 ],
the preconditioned problem will look like
[ KSPSolve(A00) 0 ] [ A00 A01 ] [ 0 KSPSolve(A00) ] [ A11 A11 ],
in other words the preconditioner is:
( [ A00 0 ] ) approximate inverse ( ) ( [ 0 A11 ] )
The preconditioner is evaluated considering only the diagonal blocks and then the full system is solved.
The system is solved using a flexible GMRES iterative method. If the system is partitioned each diagonal block is preconditioned using the (restricted) additive Schwarz method (ASM). On each block of the ASM preconditioner an incomplete LU factorization (ILU) is used. There is one ASM block per process. If the system is not partitioned it is preconditioned using the incomplete LU factorization (ILU).
Constructor.
- Parameters
-
debug if set to true, debug information will be printed
Definition at line 306 of file system_solvers_split.cpp.
◆ SplitSystemSolver() [2/6]
bitpit::SplitSystemSolver::SplitSystemSolver | ( | bool | transpose, |
bool | debug ) |
Constuctor
- Parameters
-
transpose if set to true, transposed system will be solved debug if set to true, debug information will be printed
Definition at line 317 of file system_solvers_split.cpp.
◆ SplitSystemSolver() [3/6]
bitpit::SplitSystemSolver::SplitSystemSolver | ( | bool | flatten, |
bool | transpose, | ||
bool | debug ) |
Constuctor
- Parameters
-
flatten if set to true, block size will not be taken into account when allocating the internal storage of the system matrix. Even when the internal storage is flat, block size information are available and can be used by the solver to speed up the solution of the system. However, since the internal storage doesn't take blocks into account, some low level operations (e.g., matrix-matrix multiplications) cannot use block information. Some algorithms for the solution of the system (e.g., multigrid) requires a flat storage. transpose if set to true, transposed system will be solved debug if set to true, debug information will be printed
Definition at line 334 of file system_solvers_split.cpp.
◆ SplitSystemSolver() [4/6]
bitpit::SplitSystemSolver::SplitSystemSolver | ( | const std::string & | prefix, |
bool | debug = false ) |
Constructor.
- Parameters
-
prefix is the prefix string to prepend to all option requests debug if set to true, debug information will be printed
Definition at line 345 of file system_solvers_split.cpp.
◆ SplitSystemSolver() [5/6]
bitpit::SplitSystemSolver::SplitSystemSolver | ( | const std::string & | prefix, |
bool | transpose, | ||
bool | debug ) |
Constructor.
- Parameters
-
prefix is the prefix string to prepend to all option requests debug if set to true, debug information will be printed
Definition at line 356 of file system_solvers_split.cpp.
◆ SplitSystemSolver() [6/6]
bitpit::SplitSystemSolver::SplitSystemSolver | ( | const std::string & | prefix, |
bool | flatten, | ||
bool | transpose, | ||
bool | debug ) |
Constuctor
- Parameters
-
prefix is the prefix string to prepend to all option requests flatten if set to true, block size will not be taken into account when allocating the internal storage of the system matrix. Even when the internal storage is flat, block size information are available and can be used by the solver to speed up the solution of the system. However, since the internal storage doesn't take blocks into account, some low level operations (e.g., matrix-matrix multiplications) cannot use block information. Some algorithms for the solution of the system (e.g., multigrid) requires a flat storage. transpose if set to true, transposed system will be solved debug if set to true, debug information will be printed
Definition at line 374 of file system_solvers_split.cpp.
Member Function Documentation
◆ assembly() [1/9]
void bitpit::SplitSystemSolver::assembly | ( | const Assembler & | assembler | ) |
Assembly the system.
After assembying th system solver, its options will be reset.
- Parameters
-
assembler is the matrix assembler
Definition at line 734 of file system_solvers_split.cpp.
◆ assembly() [2/9]
|
protected |
Assembly the system.
After assembying th system solver, its options will be reset.
- Parameters
-
assembler is the matrix assembler
Definition at line 285 of file system_solvers_large.cpp.
◆ assembly() [3/9]
void bitpit::SplitSystemSolver::assembly | ( | const Assembler & | assembler, |
const SystemMatrixOrdering & | reordering ) |
Assembly the system.
- Parameters
-
assembler is the matrix assembler reordering is the reordering that will be applied when assembling the system
Definition at line 745 of file system_solvers_split.cpp.
◆ assembly() [4/9]
|
protected |
Assembly the system.
After assembying th system solver, its options will be reset.
- Parameters
-
assembler is the matrix assembler reordering is the reordering that will be applied when assembling the system
Definition at line 286 of file system_solvers_large.cpp.
◆ assembly() [5/9]
|
protected |
Assembly the system.
After assembying th system solver, its options will be reset.
- Parameters
-
matrix is the matrix
Definition at line 283 of file system_solvers_large.cpp.
◆ assembly() [6/9]
|
protected |
Assembly the system.
After assembying th system solver, its options will be reset.
- Parameters
-
matrix is the matrix reordering is the reordering that will be applied when assembling the system
Definition at line 284 of file system_solvers_large.cpp.
◆ assembly() [7/9]
void bitpit::SplitSystemSolver::assembly | ( | const SparseMatrix & | matrix, |
SystemSplitStrategy | splitStrategy, | ||
const std::vector< int > & | splitSizes ) |
Assembly the system.
After assembying th system solver, its options will be reset.
- Parameters
-
matrix is the matrix splitStrategy the type of split that will be applied to the system splitSizes are the sizes of the splits
Definition at line 699 of file system_solvers_split.cpp.
◆ assembly() [8/9]
void bitpit::SplitSystemSolver::assembly | ( | const SparseMatrix & | matrix, |
SystemSplitStrategy | splitStrategy, | ||
const std::vector< int > & | splitSizes, | ||
const SystemMatrixOrdering & | reordering ) |
Assembly the system.
- Parameters
-
matrix is the matrix splitStrategy the type of split that will be applied to the system splitSizes are the sizes of the splits reordering is the reordering that will be applied when assemblying the system
Definition at line 713 of file system_solvers_split.cpp.
◆ assembly() [9/9]
|
protected |
Assembly the system.
After assembying th system solver, its options will be reset.
- Parameters
-
assembler is the matrix assembler reordering is the reordering that will be applied when assemblying the system
Definition at line 400 of file system_solvers_large.tpp.
◆ destroyKSPOptions()
|
overrideprotectedvirtual |
Destroy the options associated with the KSP.
Reimplemented from bitpit::SystemSolver.
Definition at line 548 of file system_solvers_split.cpp.
◆ destroyKSPStatus()
|
overrideprotectedvirtual |
Destroy the status of the KSP.
Reimplemented from bitpit::SystemSolver.
Definition at line 649 of file system_solvers_split.cpp.
◆ destroySplitKSPOptions()
|
protectedvirtual |
Destroy the options associated with the split KSPs.
Definition at line 557 of file system_solvers_split.cpp.
◆ destroySplitKSPStatuses()
|
protectedvirtual |
Destroy the status of the split KSPs.
Definition at line 658 of file system_solvers_split.cpp.
◆ dumpInfo()
|
overrideprotectedvirtual |
Dump system information.
- Parameters
-
systemStream is the stream in which system information is written
Reimplemented from bitpit::SystemSolver.
Definition at line 669 of file system_solvers_split.cpp.
◆ exportMatrix() [1/2]
|
overridevirtual |
Export the matrix to the specified file.
- Parameters
-
filePath is the path of the file fileFormat is the file format that will be used for exporting the matrix, note that the ASCII format may not be able to handle large matrices
Reimplemented from bitpit::SystemSolver.
Definition at line 1458 of file system_solvers_split.cpp.
◆ exportMatrix() [2/2]
|
protected |
Export the given matrix to the specified file.
- Parameters
-
matrix is the matrix that will be filled filePath is the path of the file that will contain matrix data fileFormat is the file format that will be used for exporting the matrix, note that the ASCII format may not be able to handle large matrices
Definition at line 447 of file system_solvers_large.cpp.
◆ fillKSPStatus() [1/2]
|
overrideprotectedvirtual |
Fill the status of the KSP.
Reimplemented from bitpit::SystemSolver.
Definition at line 604 of file system_solvers_split.cpp.
◆ fillKSPStatus() [2/2]
|
protectedvirtual |
Fill the status of the specified KSP.
- Parameters
-
ksp is the KSP the status will be fetched from status on output will contain the status of the KSP
Reimplemented from bitpit::SystemSolver.
Definition at line 432 of file system_solvers_large.cpp.
◆ fillSplitKSPStatuses()
|
protectedvirtual |
Fill the status of the split KSPs.
Definition at line 613 of file system_solvers_split.cpp.
◆ generateSplitIndexes()
|
protected |
Generate split indexes.
Split indexes define the indexes of the matrix that are assigned to the specified split. The k-th element of a split matrix row/column corresponds to the indexes[k]-th element of the corresponding row/column of the full matrix.
- Parameters
-
split is the split nItems is the number of items in the row/column indexes on output will contain the split indexes
Definition at line 1680 of file system_solvers_split.cpp.
◆ generateSplitPath() [1/3]
|
protected |
Generate the split path associated with the specified indexes.
The split file path is generated adding "_split_ij" before the extension of the specified path. For example, given the path "data/matrix.dat", the split path associated with the index INDEX is "data/matrix_split_INDEX.dat".
- Parameters
-
path is the path index is the index associated with the split
Definition at line 1741 of file system_solvers_split.cpp.
◆ generateSplitPath() [2/3]
|
protected |
Generate the split path associated with the specified indexes.
The split file path is generated adding "_split_ij" before the extension of the specified path. For example, given the path "data/matrix.dat", the split path associated with the index 1 is "data/matrix_split_1.dat".
- Parameters
-
path is the path i is the index associated with the split
Definition at line 1710 of file system_solvers_split.cpp.
◆ generateSplitPath() [3/3]
|
protected |
Generate the split path associated with the specified indexes.
The split file path is generated adding "_split_ij" before the extension of the specified path. For example, given the path "data/matrix.dat", the split path associated with the indexes (1,0) is "data/matrix_split_10.dat".
- Parameters
-
path is the path i is the first index associated with the split j is the second index associated with the split
Definition at line 1726 of file system_solvers_split.cpp.
◆ generateSplitPermutation()
|
protected |
Generate the split permutation.
Given a vector whose elements are ordered according to the full matrix, the split permutations allow to obtain a vector whose elements are ordered according to the split matrices.
- Parameters
-
nItems is the number of items in the row/column itemsGlobalOffset is the global offset of the items that belong to the current process split is the split [out] permutations on output will contain the permutations
Definition at line 1620 of file system_solvers_split.cpp.
◆ getBlockSize()
|
overridevirtual |
Get the block size of the system.
- Returns
- The block size of the system.
Reimplemented from bitpit::SystemSolver.
Definition at line 385 of file system_solvers_split.cpp.
◆ getSplitCount()
int bitpit::SplitSystemSolver::getSplitCount | ( | ) | const |
Get the number of splits that will be applied.
- Returns
- The number of splits that will be applied.
Definition at line 412 of file system_solvers_split.cpp.
◆ getSplitKSPOptions() [1/2]
KSPOptions & bitpit::SplitSystemSolver::getSplitKSPOptions | ( | int | split | ) |
Get a reference to the options associated to the Krylov solver of the specified split.
The options associated with the Krylov solver can only be accessed after assembling the system.
- Parameters
-
split is the split
- Returns
- A reference to the options associated to the Krylov solver of the specified split.
Definition at line 497 of file system_solvers_split.cpp.
◆ getSplitKSPOptions() [2/2]
const KSPOptions & bitpit::SplitSystemSolver::getSplitKSPOptions | ( | int | split | ) | const |
Get a constant reference to the options associated to the Krylov solver of the specified split.
The options associated with the Krylov solver can only be accessed after assembling the system.
- Parameters
-
split is the split
- Returns
- A constant reference to the options associated to the Krylov solver of the specified split.
Definition at line 515 of file system_solvers_split.cpp.
◆ getSplitKSPStatus()
const KSPStatus & bitpit::SplitSystemSolver::getSplitKSPStatus | ( | int | split | ) | const |
Get a constant reference to the status of the Krylov solver of the specified split.
The status of the Krylov solver can only be accessed after assembling the system.
- Parameters
-
split is the split
- Returns
- A constant reference to the status of the Krylov solver of the specified split.
Definition at line 570 of file system_solvers_split.cpp.
◆ getSplitOffsets()
std::vector< int > bitpit::SplitSystemSolver::getSplitOffsets | ( | ) | const |
Get the block offtest of the splits.
- Returns
- The block offsets of the splits.
Definition at line 467 of file system_solvers_split.cpp.
◆ getSplitSizes()
std::vector< int > bitpit::SplitSystemSolver::getSplitSizes | ( | ) | const |
Get the size of each split.
Each block of elements is split into pieces whose size is defined by the split sizes. The sum of the split sizes must be equal to the block size. Given a matrix with a block size equal to 5 and split sizes equal to [4, 1], the split matrix will look like this:
[ | ] [ (4 * N) * (4 * M) | (4 * N) * (1 * M) ] [ | ] [-------------------------------------------] [ | ] [ (1 * N) * (4 * M) | (1 * N) * (1 * M) ] [ | ]
where M and N are the number of rows and columns respectively.
- Returns
- The size of each split.
Definition at line 443 of file system_solvers_split.cpp.
◆ getSplitStrategy()
SystemSplitStrategy bitpit::SplitSystemSolver::getSplitStrategy | ( | ) | const |
Get the split strategy will be used when solving the system.
- Returns
- The split strategy will be used when solving the system.
Definition at line 402 of file system_solvers_split.cpp.
◆ initializeKSPOptions()
|
overrideprotectedvirtual |
Initialize the options associated with the KSP.
Reimplemented from bitpit::SystemSolver.
Definition at line 527 of file system_solvers_split.cpp.
◆ initializeKSPStatus()
|
overrideprotectedvirtual |
Initialize the status of the split KSPs.
Reimplemented from bitpit::SystemSolver.
Definition at line 581 of file system_solvers_split.cpp.
◆ initializeSplitKSPOptions()
|
protectedvirtual |
Initialize the options associated with the sub KSP.
Definition at line 536 of file system_solvers_split.cpp.
◆ initializeSplitKSPStatuses()
|
protectedvirtual |
Initialize the options associated with the split KSPs.
Definition at line 590 of file system_solvers_split.cpp.
◆ matrixAssembly()
|
protected |
Assemble the matrix.
- Parameters
-
assembler is the matrix assembler
Definition at line 818 of file system_solvers_split.cpp.
◆ matrixDestroy()
|
overrideprotectedvirtual |
Destroy the matrix.
Reimplemented from bitpit::SystemSolver.
Definition at line 1336 of file system_solvers_split.cpp.
◆ matrixDump()
|
overrideprotectedvirtual |
Dump the matrix.
- Parameters
-
systemStream is the stream in which system information is written directory is the directory in which the matrix data file will be written. prefix is the prefix added to the name of the file containing matrix data
Reimplemented from bitpit::SystemSolver.
Definition at line 1260 of file system_solvers_split.cpp.
◆ matrixFill()
|
overrideprotectedvirtual |
Fill the matrix reading its contents form the specified file.
The input file should contain a compatible vector stored in PETSc binary format. It's up to the caller of this routine to make sure the loaded vector is compatible with the matrix. If the matrix file cannot be read an exception is thrown.
- Parameters
-
filePath is the path of the file
Reimplemented from bitpit::SystemSolver.
Definition at line 1005 of file system_solvers_split.cpp.
◆ matrixRestore()
|
overrideprotectedvirtual |
Restore the matrix.
- Parameters
-
systemStream is the stream from which system information is read directory is the directory from which the matrix data file will be read prefix is the prefix that will be was added to the files during redistribute if set to true, the matrix will be redistributed among the available processes, allowing to restore the matrix with a different number of processes than those used to dump it
Reimplemented from bitpit::SystemSolver.
Definition at line 1300 of file system_solvers_split.cpp.
◆ matrixUpdate()
|
protected |
Update the specified rows of the matrix.
The contents of the specified rows will be replaced by the data provided by the given assembler. If the matrix has not been assembled yet, both the pattern and the values of the matrix will be updated. After the matrix has been assembled only the values will be updated.
The block size of the assembler should be equal to the block size of the matrix.
- Parameters
-
nRows is the number of rows that will be updated rows are the local indices of the rows that will be updated, if a null pointer is passed, the rows that will be updated are the rows from 0 to (nRows - 1). assembler is the matrix assembler for the rows that will be updated
Definition at line 1039 of file system_solvers_split.cpp.
◆ postKSPSolveActions()
|
overrideprotectedvirtual |
Post-solve actions.
Reimplemented from bitpit::SystemSolver.
Definition at line 480 of file system_solvers_split.cpp.
◆ preKrylovSetupActions()
|
overrideprotectedvirtual |
Perform actions before Krylov subspace method setup.
Reimplemented from bitpit::SystemSolver.
Definition at line 1578 of file system_solvers_split.cpp.
◆ resetKSPOptions()
|
protectedvirtual |
Reset the specified KSP options.
- Parameters
-
options are the options that will be reset
Reimplemented from bitpit::SystemSolver.
Definition at line 425 of file system_solvers_large.cpp.
◆ resetKSPStatus() [1/2]
|
overrideprotectedvirtual |
Reset the status of the KSP.
Reimplemented from bitpit::SystemSolver.
Definition at line 630 of file system_solvers_split.cpp.
◆ resetKSPStatus() [2/2]
|
protectedvirtual |
Reset the status of the specified KSP.
- Parameters
-
status on output will contain the status of the KSP
Reimplemented from bitpit::SystemSolver.
Definition at line 430 of file system_solvers_large.cpp.
◆ resetSplitKSPStatuses()
|
protectedvirtual |
Reset the status of the split KSPs.
Definition at line 639 of file system_solvers_split.cpp.
◆ restoreInfo()
|
overrideprotectedvirtual |
Restore system information.
- Parameters
-
systemStream is the stream from which system information is read
Reimplemented from bitpit::SystemSolver.
Definition at line 683 of file system_solvers_split.cpp.
◆ setupKrylov() [1/2]
|
overrideprotectedvirtual |
Set up the Krylov subspace method used to solve the system.
Reimplemented from bitpit::SystemSolver.
Definition at line 1529 of file system_solvers_split.cpp.
◆ setupKrylov() [2/2]
|
protectedvirtual |
Set up the Krylov subspace method using the given options.
This function is in charge of only setting the properties of the Krylov subspace method that will be used to solve the system. There is a dedicated function to set up the preconditioner.
- Parameters
-
ksp is the KSP whose Krylov subspace method will be setup options are the options that will be used to set up the KSP
Reimplemented from bitpit::SystemSolver.
Definition at line 416 of file system_solvers_large.cpp.
◆ setupPreconditioner() [1/2]
|
overrideprotectedvirtual |
Set up the preconditioner.
Reimplemented from bitpit::SystemSolver.
Definition at line 1477 of file system_solvers_split.cpp.
◆ setupPreconditioner() [2/2]
|
protectedvirtual |
Set up the specified preconditioner using the given options.
- Parameters
-
pc is the preconditioner to set up options are the options that will be used to set up the preconditioner
Reimplemented from bitpit::SystemSolver.
Definition at line 411 of file system_solvers_large.cpp.
◆ setupSplitKrylovs()
|
protectedvirtual |
Set up the Krylov subspace method used to solve the system.
Definition at line 1551 of file system_solvers_split.cpp.
◆ setupSplitPreconditioners()
|
protectedvirtual |
Set up the Krylov subspace method used to solve the system.
Definition at line 1508 of file system_solvers_split.cpp.
◆ SystemSolver() [1/6]
bitpit::SystemSolver::SystemSolver | ( | bool | debug = false | ) |
Constructor.
- Parameters
-
debug if set to true, debug information will be printed
Definition at line 269 of file system_solvers_large.cpp.
◆ SystemSolver() [2/6]
bitpit::SystemSolver::SystemSolver | ( | bool | flatten, |
bool | transpose, | ||
bool | debug ) |
Constuctor
- Parameters
-
flatten if set to true, block size will not be taken into account when allocating the internal storage of the system matrix. Even when the internal storage is flat, block size information are available and can be used by the solver to speed up the solution of the system. However, since the internal storage doesn't take blocks into account, some low level operations (e.g., matrix-matrix multiplications) cannot use block information. Some algorithms for the solution of the system (e.g., multigrid) requires a flat storage. transpose if set to true, transposed system will be solved debug if set to true, debug information will be printed
Definition at line 271 of file system_solvers_large.cpp.
◆ SystemSolver() [3/6]
bitpit::SystemSolver::SystemSolver | ( | bool | transpose, |
bool | debug ) |
Constuctor
- Parameters
-
transpose if set to true, transposed system will be solved debug if set to true, debug information will be printed
Definition at line 270 of file system_solvers_large.cpp.
◆ SystemSolver() [4/6]
bitpit::SystemSolver::SystemSolver | ( | const std::string & | prefix, |
bool | debug = false ) |
Constructor.
- Parameters
-
prefix is the prefix string to prepend to all option requests debug if set to true, debug information will be printed
Definition at line 272 of file system_solvers_large.cpp.
◆ SystemSolver() [5/6]
bitpit::SystemSolver::SystemSolver | ( | const std::string & | prefix, |
bool | flatten, | ||
bool | transpose, | ||
bool | debug ) |
Constuctor
- Parameters
-
prefix is the prefix string to prepend to all option requests flatten if set to true, block size will not be taken into account when allocating the internal storage of the system matrix. Even when the internal storage is flat, block size information are available and can be used by the solver to speed up the solution of the system. However, since the internal storage doesn't take blocks into account, some low level operations (e.g., matrix-matrix multiplications) cannot use block information. Some algorithms for the solution of the system (e.g., multigrid) requires a flat storage. transpose if set to true, transposed system will be solved debug if set to true, debug information will be printed
Definition at line 274 of file system_solvers_large.cpp.
◆ SystemSolver() [6/6]
bitpit::SystemSolver::SystemSolver | ( | const std::string & | prefix, |
bool | transpose, | ||
bool | debug ) |
Constructor.
- Parameters
-
prefix is the prefix string to prepend to all option requests debug if set to true, debug information will be printed
Definition at line 273 of file system_solvers_large.cpp.
◆ update() [1/7]
void bitpit::SplitSystemSolver::update | ( | const Assembler & | assembler | ) |
Update all the rows of the system.
Only the values of the system matrix can be updated, once the system is assembled its pattern cannot be modified.
- Parameters
-
assembler is the matrix assembler for the rows that will be updated
Definition at line 793 of file system_solvers_split.cpp.
◆ update() [2/7]
|
protected |
Update all the rows of the system.
Only the values of the system matrix can be updated, once the system is assembled its pattern cannot be modified.
- Parameters
-
assembler is the matrix assembler for the rows that will be updated
Definition at line 290 of file system_solvers_large.cpp.
◆ update() [3/7]
void bitpit::SplitSystemSolver::update | ( | const SparseMatrix & | elements | ) |
Update all the rows of the system.
Only the values of the system matrix can be updated, once the system is assembled its pattern cannot be modified.
- Parameters
-
elements are the elements that will be used to update the rows
Definition at line 758 of file system_solvers_split.cpp.
◆ update() [4/7]
void bitpit::SplitSystemSolver::update | ( | long | nRows, |
const long * | rows, | ||
const Assembler & | assembler ) |
Update the system.
Only the values of the system matrix can be updated, once the system is assembled its pattern cannot be modified.
- Parameters
-
nRows is the number of rows that will be updated rows are the indices of the rows that will be updated assembler is the matrix assembler for the rows that will be updated
Definition at line 808 of file system_solvers_split.cpp.
◆ update() [5/7]
|
protected |
Update the system.
Only the values of the system matrix can be updated, once the system is assembled its pattern cannot be modified.
- Parameters
-
nRows is the number of rows that will be updated rows are the indices of the rows that will be updated assembler is the matrix assembler for the rows that will be updated
Definition at line 291 of file system_solvers_large.cpp.
◆ update() [6/7]
void bitpit::SplitSystemSolver::update | ( | long | nRows, |
const long * | rows, | ||
const SparseMatrix & | elements ) |
Update the system.
Only the values of the system matrix can be updated, once the system is assembled its pattern cannot be modified.
- Parameters
-
nRows is the number of rows that will be updated rows are the indices of the rows that will be updated elements are the elements that will be used to update the rows
Definition at line 773 of file system_solvers_split.cpp.
◆ update() [7/7]
|
protected |
Update the system.
Only the values of the system matrix can be updated, once the system is assembled its pattern cannot be modified.
- Parameters
-
nRows is the number of rows that will be updated rows are the indices of the rows that will be updated assembler is the matrix assembler for the rows that will be updated
Definition at line 403 of file system_solvers_large.tpp.
◆ vectorsCreate()
|
overrideprotectedvirtual |
Create RHS and solution vectors.
Vectors will be created, but they will not be initialized.
Note that MatNest does not depend on VecNest, and in general, Jed Brown (one of the main developers of PETSc) recommend not using VecNest. For this reason, VecNest will not be used.
Reimplemented from bitpit::SystemSolver.
Definition at line 1353 of file system_solvers_split.cpp.
◆ vectorsCreateSplitPermutations()
|
protected |
Create the split permutation for RHS and solution vectors.
Given a vector whose elements are ordered according to the full matrix, the split permutations allow to obtain a vector whose elements are ordered according to the split matrices.
Definition at line 1398 of file system_solvers_split.cpp.
◆ vectorsDestroy()
|
overrideprotectedvirtual |
Destroy RHS and solution vectors.
Reimplemented from bitpit::SystemSolver.
Definition at line 1436 of file system_solvers_split.cpp.
◆ vectorsReorder()
|
overrideprotectedvirtual |
Reorder RHS and solution vectors to match the order of the system matrix.
- Parameters
-
invert is a flag for inverting the ordering
Reimplemented from bitpit::SystemSolver.
Definition at line 1367 of file system_solvers_split.cpp.
◆ vectorsRestore()
|
overrideprotectedvirtual |
Restore RHS and solution vectors.
- Parameters
-
systemStream is the stream from which system information is read directory is the directory from which the vector data file will be read prefix is the prefix added to the name of the file containing vectors data
Reimplemented from bitpit::SystemSolver.
Definition at line 1382 of file system_solvers_split.cpp.
Friends And Related Symbol Documentation
◆ SystemSolver
|
friend |
Definition at line 73 of file system_solvers_split.hpp.
Member Data Documentation
◆ m_splitAs
|
protected |
Definition at line 112 of file system_solvers_split.hpp.
◆ m_splitKSPOptions
|
protected |
Definition at line 114 of file system_solvers_split.hpp.
◆ m_splitKSPStatuses
|
protected |
Definition at line 115 of file system_solvers_split.hpp.
◆ m_splitStrategy
|
protected |
Definition at line 110 of file system_solvers_split.hpp.
The documentation for this class was generated from the following files:
- src/LA/system_solvers_split.hpp
- src/LA/system_solvers_split.cpp
