Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
bitpit::SystemSolver Class Reference

The SystemSolver class provides methods for building and solving large linear systems. More...

Inheritance diagram for bitpit::SystemSolver:
Inheritance graph
[legend]
Collaboration diagram for bitpit::SystemSolver:
Collaboration graph
[legend]

Public Types

typedef SystemMatrixAssembler Assembler
 
enum  FileFormat { FILE_BINARY , FILE_ASCII }
 

Public Member Functions

 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 exportMatrix (const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
 
virtual void exportRHS (const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
 
virtual void exportSolution (const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
 
virtual int getBlockSize () const
 
long getColCount () const
 
long getColElementCount () const
 
long getColGlobalCount () const
 
long getColGlobalElementCount () const
 
KSPOptionsgetKSPOptions ()
 
const KSPOptionsgetKSPOptions () const
 
const KSPStatusgetKSPStatus () 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
 
bool isSetUp () 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 setUp ()
 
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 Types

enum  VectorSide { VECTOR_SIDE_RIGHT , VECTOR_SIDE_LEFT }
 

Protected Member Functions

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 ()
 
virtual void destroyKSPOptions ()
 
virtual void destroyKSPStatus ()
 
void destroyMatrix (Mat *matrix) const
 
void destroyVector (Vec *vector) const
 
virtual void dumpInfo (std::ostream &stream) 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
 
virtual void fillKSPStatus ()
 
virtual void fillKSPStatus (KSP ksp, KSPStatus *status) 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
 
virtual void initializeKSPOptions ()
 
virtual void initializeKSPStatus ()
 
void matrixAssembly (const Assembler &assembler)
 
virtual void matrixDestroy ()
 
virtual void matrixDump (std::ostream &systemStream, const std::string &directory, const std::string &prefix) const
 
virtual void matrixFill (const std::string &filePath)
 
virtual void matrixRestore (std::istream &systemStream, const std::string &directory, const std::string &prefix, bool redistribute)
 
void matrixUpdate (long nRows, const long *rows, const Assembler &assembler)
 
virtual void postKrylovSetupActions ()
 
virtual void postKSPSolveActions ()
 
virtual void postPreconditionerSetupActions ()
 
virtual void preKrylovSetupActions ()
 
virtual void preKSPSolveActions ()
 
virtual void prepareKSP ()
 
virtual void prePreconditionerSetupActions ()
 
void reorderVector (Vec vector, IS permutations, bool invert) const
 
virtual void resetKSPOptions (KSPOptions *options) const
 
virtual void resetKSPStatus ()
 
virtual void resetKSPStatus (KSPStatus *status) const
 
virtual void restoreInfo (std::istream &stream)
 
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 setupKrylov ()
 
virtual void setupKrylov (KSP ksp, const KSPOptions &options) const
 
virtual void setupPreconditioner ()
 
virtual void setupPreconditioner (PC pc, const KSPOptions &options) const
 
virtual void solveKSP ()
 
template<typename DerivedSystemSolver >
void update (long nRows, const long *rows, const typename DerivedSystemSolver::Assembler &assembler)
 
virtual void vectorsCreate ()
 
virtual void vectorsDestroy ()
 
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)
 
virtual void vectorsReorder (bool invert)
 
virtual void vectorsRestore (std::istream &systemStream, const std::string &directory, const std::string &prefix)
 

Protected Attributes

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
 

Detailed Description

The SystemSolver class provides methods for building and solving large linear systems.

Rather than working with individual elements in the system matrix, it is possible to employ blocks of elements. The size of the blocks can be defined during assembly. When a size different that one is provided, the matrix will store elements by fixed-sized dense nb × nb blocks, where nb is the size of the blocks. Blocking may be advantageous when solving PDE-based simulations that leads to matrices with a naturally blocked structure (with a block size equal to the number of degrees of freedom per cell).

When blocking is used, row and column indexes will count the number of blocks in the row/column direction, not the number of rows/columns of the matrix.

Examples
LA_example_00001.cpp, and RBF_example_00001.cpp.

Definition at line 264 of file system_solvers_large.hpp.

Member Typedef Documentation

◆ Assembler

Definition at line 267 of file system_solvers_large.hpp.

Member Enumeration Documentation

◆ FileFormat

enum bitpit::SystemSolver::FileFormat

Definition at line 269 of file system_solvers_large.hpp.

◆ VectorSide

enum bitpit::SystemSolver::VectorSide
protected

Definition at line 361 of file system_solvers_large.hpp.

Constructor & Destructor Documentation

◆ SystemSolver() [1/6]

bitpit::SystemSolver::SystemSolver ( bool debug = false)

Constructor.

Parameters
debugif set to true, debug information will be printed
Examples
LA_example_00001.cpp.

Definition at line 770 of file system_solvers_large.cpp.

◆ SystemSolver() [2/6]

bitpit::SystemSolver::SystemSolver ( bool transpose,
bool debug )

Constuctor

Parameters
transposeif set to true, transposed system will be solved
debugif set to true, debug information will be printed

Definition at line 781 of file system_solvers_large.cpp.

◆ SystemSolver() [3/6]

bitpit::SystemSolver::SystemSolver ( bool flatten,
bool transpose,
bool debug )

Constuctor

Parameters
flattenif 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.
transposeif set to true, transposed system will be solved
debugif set to true, debug information will be printed

Definition at line 799 of file system_solvers_large.cpp.

◆ SystemSolver() [4/6]

bitpit::SystemSolver::SystemSolver ( const std::string & prefix,
bool debug = false )

Constructor.

Parameters
prefixis the prefix string to prepend to all option requests
debugif set to true, debug information will be printed

Definition at line 810 of file system_solvers_large.cpp.

◆ SystemSolver() [5/6]

bitpit::SystemSolver::SystemSolver ( const std::string & prefix,
bool transpose,
bool debug )

Constructor.

Parameters
prefixis the prefix string to prepend to all option requests
debugif set to true, debug information will be printed

Definition at line 821 of file system_solvers_large.cpp.

◆ SystemSolver() [6/6]

bitpit::SystemSolver::SystemSolver ( const std::string & prefix,
bool flatten,
bool transpose,
bool debug )

Constuctor

Parameters
prefixis the prefix string to prepend to all option requests
flattenif 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.
transposeif set to true, transposed system will be solved
debugif set to true, debug information will be printed

Definition at line 839 of file system_solvers_large.cpp.

◆ ~SystemSolver()

bitpit::SystemSolver::~SystemSolver ( )
virtual

Destructor.

Definition at line 863 of file system_solvers_large.cpp.

Member Function Documentation

◆ assembly() [1/5]

void bitpit::SystemSolver::assembly ( const Assembler & assembler)

Assembly the system.

After assembying th system solver, its options will be reset.

Parameters
assembleris the matrix assembler

Definition at line 949 of file system_solvers_large.cpp.

◆ assembly() [2/5]

void bitpit::SystemSolver::assembly ( const Assembler & assembler,
const SystemMatrixOrdering & reordering )

Assembly the system.

After assembying th system solver, its options will be reset.

Parameters
assembleris the matrix assembler
reorderingis the reordering that will be applied when assembling the system

Definition at line 962 of file system_solvers_large.cpp.

◆ assembly() [3/5]

void bitpit::SystemSolver::assembly ( const SparseMatrix & matrix)

Assembly the system.

After assembying th system solver, its options will be reset.

Parameters
matrixis the matrix
Examples
RBF_example_00001.cpp.

Definition at line 917 of file system_solvers_large.cpp.

◆ assembly() [4/5]

void bitpit::SystemSolver::assembly ( const SparseMatrix & matrix,
const SystemMatrixOrdering & reordering )

Assembly the system.

After assembying th system solver, its options will be reset.

Parameters
matrixis the matrix
reorderingis the reordering that will be applied when assembling the system

Definition at line 930 of file system_solvers_large.cpp.

◆ assembly() [5/5]

template<typename DerivedSystemSolver >
void bitpit::SystemSolver::assembly ( const typename DerivedSystemSolver::Assembler & assembler,
const SystemMatrixOrdering & reordering )
protected

Assembly the system.

After assembying th system solver, its options will be reset.

Parameters
assembleris the matrix assembler
reorderingis the reordering that will be applied when assemblying the system

Definition at line 95 of file system_solvers_large.tpp.

◆ clear()

void bitpit::SystemSolver::clear ( )
virtual

Clear the system

Reimplemented in bitpit::DiscretizationStencilSolver< stencil_t, solver_kernel_t >.

Definition at line 875 of file system_solvers_large.cpp.

◆ clearReordering()

void bitpit::SystemSolver::clearReordering ( )
protected

Clear the reordering that will be applied when assembling the matrix.

The function will clear any reordering previously set. With no reordering defined, the matrix will be assembled using its natural ordering.

Definition at line 3146 of file system_solvers_large.cpp.

◆ clearWorkspace()

void bitpit::SystemSolver::clearWorkspace ( )
virtual

Clear and release the memory of all data structures needed for the solution of the system.

These data structures will be re-created the next time the system will be solved.

Definition at line 899 of file system_solvers_large.cpp.

◆ createKSP()

void bitpit::SystemSolver::createKSP ( )
protectedvirtual

Create the KSP.

Definition at line 3239 of file system_solvers_large.cpp.

◆ createMatrix() [1/2]

void bitpit::SystemSolver::createMatrix ( int rowBlockSize,
int colBlockSize,
int nNestRows,
int nNestCols,
Mat * subMatrices,
Mat * matrix ) const
protected

Create a nest matrix.

Parameters
rowBlockSizeis the row block size of the matrix
colBlockSizeis the column block size of the matrix
nNestRowsis the number of rows in the nest, i.e. the number of submatrices along the rows of the nest
nNestColsis the number of columns in the nest, i.e. the number of submatrices along the columns of the nest
subMatricesare the submatrices stored in row-major order, empty submatrices can be passed using nullptr
vectoron output will contain the newly created matrix

Definition at line 2321 of file system_solvers_large.cpp.

◆ createMatrix() [2/2]

void bitpit::SystemSolver::createMatrix ( int rowBlockSize,
int colBlockSize,
Mat * matrix ) const
protected

Create a matrix.

Parameters
rowBlockSizeis the row block size of the matrix
colBlockSizeis the column block size of the matrix
matrixon output will contain the newly created matrix

Definition at line 2266 of file system_solvers_large.cpp.

◆ createVector() [1/2]

void bitpit::SystemSolver::createVector ( int blockSize,
int nestSize,
Vec * subVectors,
Vec * vector ) const
protected

Create a nest vector.

Parameters
blockSizeis the block size of the vector
nestSizeis the number of subvectors that will be contained in the nest
subVectorsare the subvectors, empty subvectors can be passed using nullptr
vectoron output will contain the newly created vector

Definition at line 2626 of file system_solvers_large.cpp.

◆ createVector() [2/2]

void bitpit::SystemSolver::createVector ( int blockSize,
Vec * vector ) const
protected

Create a vector.

Parameters
blockSizeis the block size of the vector
vectoron output will contain the newly created vector

Definition at line 2593 of file system_solvers_large.cpp.

◆ destroyKSP()

void bitpit::SystemSolver::destroyKSP ( )
protectedvirtual

Destroy the KSP.

Definition at line 3257 of file system_solvers_large.cpp.

◆ destroyKSPOptions()

void bitpit::SystemSolver::destroyKSPOptions ( )
protectedvirtual

Destroy the options associated with the KSP.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3465 of file system_solvers_large.cpp.

◆ destroyKSPStatus()

void bitpit::SystemSolver::destroyKSPStatus ( )
protectedvirtual

Destroy the status of the KSP.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3538 of file system_solvers_large.cpp.

◆ destroyMatrix()

void bitpit::SystemSolver::destroyMatrix ( Mat * matrix) const
protected

Destroy the specified matrix.

Parameters
matrixis the matrix that will be destroyed

Definition at line 2579 of file system_solvers_large.cpp.

◆ destroyVector()

void bitpit::SystemSolver::destroyVector ( Vec * vector) const
protected

Destroy the specified vector.

Parameters
vectoris the vector that will be destroyed

Definition at line 3010 of file system_solvers_large.cpp.

◆ dumpInfo()

void bitpit::SystemSolver::dumpInfo ( std::ostream & systemStream) const
protectedvirtual

Dump system information.

Parameters
systemStreamis the stream in which system information is written

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 2230 of file system_solvers_large.cpp.

◆ dumpMatrix()

void bitpit::SystemSolver::dumpMatrix ( Mat matrix,
const std::string & directory,
const std::string & name ) const
protected

Dump the specified matrix.

In addition to the file that contains the contant of the matrix, some other files will be created:

  • a file with the information needed for re-creating the matrix (called "<NAME>.info");

a file with the information about the partitioning (called "<NAME>.partitioning"), this file is created only if the system is partitioned and bitpit MPI support is enabled.

Parameters
matrixis the matrix that will be dumped
directoryis the directory in which the matrix data file will be written
nameis the name of the matrix that will be dumped

Definition at line 2390 of file system_solvers_large.cpp.

◆ dumpSystem()

void bitpit::SystemSolver::dumpSystem ( const std::string & header,
const std::string & directory,
const std::string & prefix = "" ) const

Dump the system to files.

Only the contents of the system will be dumped, this include the matrix, the RHS vector, the solution vector, and the information needed to properly restore the aforementioned data structures (e.g., the block size or the transpose flag).

Parameters
headeris the header that will be written in the system information archive
directoryis the directory where the files will be saved
prefixis the prefix that will be added to the files

Definition at line 2132 of file system_solvers_large.cpp.

◆ dumpVector()

void bitpit::SystemSolver::dumpVector ( Vec vector,
const std::string & directory,
const std::string & name ) const
protected

Dump the specified vector.

In addition to the file that contains the contant of the vector, some other files will be created:

  • a file with the information needed for re-creating the vector (called "<NAME>.info");

a file with the information about the partitioning (called "<NAME>.partitioning"), this file is created only if the system is partitioned and bitpit MPI support is enabled.

Parameters
vectoris the vector that will be dumped
directoryis the directory in which the vector data file will be written
nameis the name of the vector that will be dumped

Definition at line 2711 of file system_solvers_large.cpp.

◆ enableForceConsistency()

void bitpit::SystemSolver::enableForceConsistency ( bool enable)

Enable or disable forcing right hand side consistency.

Parameters
enableif set to true, right hand side consistency will be forced before solving the system

Definition at line 3617 of file system_solvers_large.cpp.

◆ exportMatrix() [1/2]

void bitpit::SystemSolver::exportMatrix ( const std::string & filePath,
FileFormat fileFormat = FILE_BINARY ) const
virtual

Export the matrix to the specified file.

Parameters
filePathis the path of the file
fileFormatis 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 in bitpit::SplitSystemSolver, and bitpit::SplitSystemSolver.

Definition at line 1980 of file system_solvers_large.cpp.

◆ exportMatrix() [2/2]

void bitpit::SystemSolver::exportMatrix ( Mat matrix,
const std::string & filePath,
FileFormat fileFormat ) const
protected

Export the given matrix to the specified file.

Parameters
matrixis the matrix that will be filled
filePathis the path of the file that will contain matrix data
fileFormatis 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 2531 of file system_solvers_large.cpp.

◆ exportRHS()

void bitpit::SystemSolver::exportRHS ( const std::string & filePath,
FileFormat fileFormat = FILE_BINARY ) const
virtual

Export the RHS vector to the specified file.

Parameters
filePathis the path of the file
fileFormatis the file format that will be used for exporting the RHS vector, note that the ASCII format may not be able to handle large vectors

Definition at line 2022 of file system_solvers_large.cpp.

◆ exportSolution()

void bitpit::SystemSolver::exportSolution ( const std::string & filePath,
FileFormat fileFormat = FILE_BINARY ) const
virtual

Export the solution vector to the specified file.

Parameters
filePathis the path of the file
fileFormatis the file format that will be used for exporting the solution vector, note that the ASCII format may not be able to handle large vectors

Definition at line 2075 of file system_solvers_large.cpp.

◆ exportVector() [1/2]

void bitpit::SystemSolver::exportVector ( Vec vector,
const std::string & filePath,
FileFormat fileFormat ) const
protected

Export the specified vector to the file given file.

Parameters
vectoris the vector that will be exported
filePathis the path of the file
fileFormatis the file format that will be used for exporting the vector, note that the ASCII format may not be able to handle large vectors

Definition at line 2943 of file system_solvers_large.cpp.

◆ exportVector() [2/2]

void bitpit::SystemSolver::exportVector ( Vec vector,
std::vector< double > * data ) const
protected

Export the specified vector into the given container.

Parameters
vectoris the vector that will be exported
dataon output it will contain the data of the vector of the linear system

Definition at line 2992 of file system_solvers_large.cpp.

◆ fillKSPStatus() [1/2]

void bitpit::SystemSolver::fillKSPStatus ( )
protectedvirtual

Fill the status of the KSP.

Reimplemented in bitpit::SplitSystemSolver, and bitpit::SplitSystemSolver.

Definition at line 3497 of file system_solvers_large.cpp.

◆ fillKSPStatus() [2/2]

void bitpit::SystemSolver::fillKSPStatus ( KSP ksp,
KSPStatus * status ) const
protectedvirtual

Fill the status of the specified KSP.

Parameters
kspis the KSP the status will be fetched from
statuson output will contain the status of the KSP

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3508 of file system_solvers_large.cpp.

◆ fillMatrix()

void bitpit::SystemSolver::fillMatrix ( Mat matrix,
const std::string & filePath ) const
protected

Fill the given matrix reading its contents from the specified file.

The input file should contain a compatible matrix stored in PETSc binary format. If the data file cannot be read an exception is thrown.

Parameters
matrixis the matrix that will be filled
filePathis the path of the file that contains matrix data

Definition at line 2348 of file system_solvers_large.cpp.

◆ fillVector() [1/2]

void bitpit::SystemSolver::fillVector ( Vec vector,
const std::string & filePath ) const
protected

Fill the given vector reading its contents from the specified file.

The input file should contain a compatible vector stored in PETSc binary format. If the data file cannot be read an exception is thrown.

Parameters
vectoris the vector that will be filled
filePathis the path of the file that contains vector data

Definition at line 2650 of file system_solvers_large.cpp.

◆ fillVector() [2/2]

void bitpit::SystemSolver::fillVector ( Vec vector,
const std::vector< double > & data ) const
protected

Fill the specified vector with the given data.

Parameters
vectoris the vector that will be filled
datais the data that will be copied into the vector

Definition at line 2679 of file system_solvers_large.cpp.

◆ finalizeKSP()

void bitpit::SystemSolver::finalizeKSP ( )
protectedvirtual

Finalize the KSP after the solution of the system.

Definition at line 3231 of file system_solvers_large.cpp.

◆ getBlockSize()

int bitpit::SystemSolver::getBlockSize ( ) const
virtual

Get the block size of the system.

Returns
The block size of the system.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1022 of file system_solvers_large.cpp.

◆ getColCount()

long bitpit::SystemSolver::getColCount ( ) const

Get the number of columns of the system.

If the matrix is a block matrix (i.e., the block size is greater than one), this function will return the number of block columns, where a block column is defined as a group of block-size matrix columns.

Returns
The number of columns of the system.

Definition at line 1065 of file system_solvers_large.cpp.

◆ getColElementCount()

long bitpit::SystemSolver::getColElementCount ( ) const

Get the number of elements in the columns of the system.

This function will return the effective number of columns of the system matrix. This value is the same as the local size used in creating the x vector for the matrix-vector product y = Ax.

Returns
The number of elements in the columns of the system.

Definition at line 1108 of file system_solvers_large.cpp.

◆ getColGlobalCount()

long bitpit::SystemSolver::getColGlobalCount ( ) const

Get number of global (block) columns.

If the matrix is a block matrix (i.e., the block size is greater than one), this function will return the global number of block columns, where a block column is defined as a group of block-size matrix columns.

Returns
The number of global (block) columns.

Definition at line 1152 of file system_solvers_large.cpp.

◆ getColGlobalElementCount()

long bitpit::SystemSolver::getColGlobalElementCount ( ) const

Get the number of global elements in the columns of the system.

This function will return the effective global number of columns of the system matrix.

Returns
The number of global elements in the columns of the system.

Definition at line 1193 of file system_solvers_large.cpp.

◆ getCommunicator()

const MPI_Comm & bitpit::SystemSolver::getCommunicator ( ) const
protected

Gets the MPI communicator associated to the system.

Returns
The MPI communicator associated to the system.

Definition at line 3549 of file system_solvers_large.cpp.

◆ getDataFilePath()

std::string bitpit::SystemSolver::getDataFilePath ( const std::string & directory,
const std::string & name ) const
protected

Get the path of the specified data file.

Parameters
directoryis the directory that contains the file
nameis the name, without extension, of the file
Returns
The path of the data file.

Definition at line 3037 of file system_solvers_large.cpp.

◆ getDumpVersion()

int bitpit::SystemSolver::getDumpVersion ( ) const
protectedvirtual

Get the version associated with the binary dumps.

Returns
The version associated with the binary dumps.

Definition at line 1334 of file system_solvers_large.cpp.

◆ getFilePath()

std::string bitpit::SystemSolver::getFilePath ( const std::string & directory,
const std::string & name ) const
protected

Get the path of the specified file.

Parameters
directoryis the directory that contains the file
nameis the name of the file
Returns
The path of the file.

Definition at line 3049 of file system_solvers_large.cpp.

◆ getInfoFilePath()

std::string bitpit::SystemSolver::getInfoFilePath ( const std::string & directory,
const std::string & name ) const
protected

Get the path of the specified info file.

Parameters
directoryis the directory that contains the file
nameis the name, without extension, of the file
Returns
The path of the info file.

Definition at line 3025 of file system_solvers_large.cpp.

◆ getKSPOptions() [1/2]

KSPOptions & bitpit::SystemSolver::getKSPOptions ( )

Get a reference to the options associated with the Krylov solver.

The options associated with the Krylov solver can only be accessed after assembling the system.

Returns
A reference to the options associated with the Krylov solver.

Definition at line 3419 of file system_solvers_large.cpp.

◆ getKSPOptions() [2/2]

const KSPOptions & bitpit::SystemSolver::getKSPOptions ( ) const

Get a constant reference to the options associated with the Krylov solver.

The options associated with the Krylov solver can only be accessed after assembling the system.

Returns
A constant reference to the options associated with the Krylov solver.

Definition at line 3435 of file system_solvers_large.cpp.

◆ getKSPStatus()

const KSPStatus & bitpit::SystemSolver::getKSPStatus ( ) const

Get a constant reference to the status of the Krylov solver.

The status of the Krylov solver can only be accessed after assembling the system.

Returns
A constant reference to the status of the Krylov solver.

Definition at line 3477 of file system_solvers_large.cpp.

◆ getRHSRawPtr() [1/2]

double * bitpit::SystemSolver::getRHSRawPtr ( )

Get a raw pointer to the solution vector.

Returns
A raw pointer to the solution vector.
Examples
RBF_example_00001.cpp.

Definition at line 1862 of file system_solvers_large.cpp.

◆ getRHSRawPtr() [2/2]

const double * bitpit::SystemSolver::getRHSRawPtr ( ) const

Get a constant raw pointer to the solution vector.

Returns
A constant raw pointer to the solution vector.

Definition at line 1875 of file system_solvers_large.cpp.

◆ getRHSRawReadPtr()

const double * bitpit::SystemSolver::getRHSRawReadPtr ( ) const

Get a constant raw pointer to the solution vector.

Returns
A constant raw pointer to the solution vector.

Definition at line 1885 of file system_solvers_large.cpp.

◆ getRowCount()

long bitpit::SystemSolver::getRowCount ( ) const

Get the number of rows of the system.

If the matrix is a block matrix (i.e., the block size is greater than one), this function will return the number of block rows, where a block row is defined as a group of block-size matrix rows.

Returns
The number of rows of the system.

Definition at line 1043 of file system_solvers_large.cpp.

◆ getRowElementCount()

long bitpit::SystemSolver::getRowElementCount ( ) const

Get the number of elements in the rows of the system.

This function will return the effective number of rows of the system matrix. This value is the same as the local size used in creating the y vector for the matrix-vector product y = Ax.

Returns
The number of elements in the rows of the system.

Definition at line 1087 of file system_solvers_large.cpp.

◆ getRowGlobalCount()

long bitpit::SystemSolver::getRowGlobalCount ( ) const

Get number of global (block) rows.

If the matrix is a block matrix (i.e., the block size is greater than one), this function will return the global number of block rows, where a block row is defined as a group of block-size matrix rows.

Returns
The number of global rows

Definition at line 1130 of file system_solvers_large.cpp.

◆ getRowGlobalElementCount()

long bitpit::SystemSolver::getRowGlobalElementCount ( ) const

Get the number of global elements in the rows of the system.

This function will return the effective global number of rows of the system matrix.

Returns
The number of global elements in the rows of the system.

Definition at line 1173 of file system_solvers_large.cpp.

◆ getSolutionRawPtr() [1/2]

double * bitpit::SystemSolver::getSolutionRawPtr ( )

Get a raw pointer to the solution vector.

Returns
A raw pointer to the solution vector.
Examples
RBF_example_00001.cpp.

Definition at line 1920 of file system_solvers_large.cpp.

◆ getSolutionRawPtr() [2/2]

const double * bitpit::SystemSolver::getSolutionRawPtr ( ) const

Get a constant raw pointer to the solution vector.

Returns
A constant raw pointer to the solution vector.

Definition at line 1933 of file system_solvers_large.cpp.

◆ getSolutionRawReadPtr()

const double * bitpit::SystemSolver::getSolutionRawReadPtr ( ) const

Get a constant raw pointer to the solution vector.

Returns
A constant raw pointer to the solution vector.
Examples
RBF_example_00001.cpp.

Definition at line 1943 of file system_solvers_large.cpp.

◆ getTranspose()

bool bitpit::SystemSolver::getTranspose ( ) const

Get the transpose flag.

Returns
Returns true if the transposed system will be solved, false otherwise.

Definition at line 165 of file system_solvers_large.cpp.

◆ importMatrix()

void bitpit::SystemSolver::importMatrix ( const std::string & filePath)
virtual

Import the matrix reading its contents form the specified file.

The input file should contain a compatible matrix stored in PETSc binary format. It's up to the caller of this routine to make sure the loaded matrix is compatible with the system. If the matrix file cannot be read an exception is thrown.

Parameters
filePathis the path of the file

Definition at line 1994 of file system_solvers_large.cpp.

◆ importRHS()

void bitpit::SystemSolver::importRHS ( const std::string & filePath)
virtual

Import the RHS vector reading its contents from the specified file.

The input file should contain a compatible vector stored in PETSc binary format. If the size of the vector stored in the file is not compatible with the matrix, an exception is thrown. An exception is also raised if the file cannot be read.

It is possible to fill the RHS vector only after the system has been assembled.

Parameters
filePathis the path of the file that contains RHS vector data

Definition at line 2038 of file system_solvers_large.cpp.

◆ importSolution()

void bitpit::SystemSolver::importSolution ( const std::string & filePath)
virtual

Import the solution vector reading its contents from the specified file.

The input file should contain a compatible vector stored in PETSc binary format. If the size of the vector stored in the file is not compatible with the matrix, an exception is thrown. An exception is also raised if the file cannot be read.

It is possible to fill the solution vector only after the system has been assembled.

Parameters
filePathis the path of the file that contains solution vector data

Definition at line 2091 of file system_solvers_large.cpp.

◆ initializeKSPOptions()

void bitpit::SystemSolver::initializeKSPOptions ( )
protectedvirtual

Initialize the options associated with the KSP.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3447 of file system_solvers_large.cpp.

◆ initializeKSPStatus()

void bitpit::SystemSolver::initializeKSPStatus ( )
protectedvirtual

Initialize the status of the KSP.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3489 of file system_solvers_large.cpp.

◆ isAssembled()

bool bitpit::SystemSolver::isAssembled ( ) const

Check if the system is assembled.

Returns
Returns true if the system is assembled, false otherwise.

Definition at line 1221 of file system_solvers_large.cpp.

◆ isForceConsistencyEnabled()

bool bitpit::SystemSolver::isForceConsistencyEnabled ( ) const

Check if right hand side consistency is forced before solving the system.

In order to force right hand side consistency, null space components are removed from the right hand side before solving the system.

Definition at line 3606 of file system_solvers_large.cpp.

◆ isPartitioned()

bool bitpit::SystemSolver::isPartitioned ( ) const

Checks if the matrix is partitioned.

Returns
Returns true if the patch is partitioned, false otherwise.
Examples
LA_example_00001.cpp.

Definition at line 1210 of file system_solvers_large.cpp.

◆ isSetUp()

bool bitpit::SystemSolver::isSetUp ( ) const

Check if the system is set up.

This function is deprecated, it will always return true.

Returns
Returns true if the system is set up, false otherwise.

Definition at line 1233 of file system_solvers_large.cpp.

◆ matrixAssembly()

void bitpit::SystemSolver::matrixAssembly ( const Assembler & assembler)
protected

Assemble the matrix.

Parameters
assembleris the matrix assembler

Definition at line 1346 of file system_solvers_large.cpp.

◆ matrixDestroy()

void bitpit::SystemSolver::matrixDestroy ( )
protectedvirtual

Destroy the matrix.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1743 of file system_solvers_large.cpp.

◆ matrixDump()

void bitpit::SystemSolver::matrixDump ( std::ostream & systemStream,
const std::string & directory,
const std::string & prefix ) const
protectedvirtual

Dump the matrix.

Parameters
systemStreamis the stream in which system information is written
directoryis the directory in which the matrix data file will be written
prefixis the prefix added to the name of the file containing matrix data

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1703 of file system_solvers_large.cpp.

◆ matrixFill()

void bitpit::SystemSolver::matrixFill ( const std::string & filePath)
protectedvirtual

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
filePathis the path of the file

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1491 of file system_solvers_large.cpp.

◆ matrixRestore()

void bitpit::SystemSolver::matrixRestore ( std::istream & systemStream,
const std::string & directory,
const std::string & prefix,
bool redistribute )
protectedvirtual

Restore the matrix.

Parameters
systemStreamis the stream from which system information is read
directoryis the directory from which the matrix data file will be read
prefixis the prefix that will be was added to the files during
redistributeif 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 in bitpit::SplitSystemSolver.

Definition at line 1724 of file system_solvers_large.cpp.

◆ matrixUpdate()

void bitpit::SystemSolver::matrixUpdate ( long nRows,
const long * rows,
const Assembler & assembler )
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
nRowsis the number of rows that will be updated
rowsare 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).
assembleris the matrix assembler for the rows that will be updated

Definition at line 1518 of file system_solvers_large.cpp.

◆ postKrylovSetupActions()

void bitpit::SystemSolver::postKrylovSetupActions ( )
protectedvirtual

Perform actions after Krylov subspace method setup.

Definition at line 3408 of file system_solvers_large.cpp.

◆ postKSPSolveActions()

void bitpit::SystemSolver::postKSPSolveActions ( )
protectedvirtual

Post-solve actions.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1317 of file system_solvers_large.cpp.

◆ postPreconditionerSetupActions()

void bitpit::SystemSolver::postPreconditionerSetupActions ( )
protectedvirtual

Perform actions after preconditioner setup.

Definition at line 3351 of file system_solvers_large.cpp.

◆ preKrylovSetupActions()

void bitpit::SystemSolver::preKrylovSetupActions ( )
protectedvirtual

Perform actions before Krylov subspace method setup.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3389 of file system_solvers_large.cpp.

◆ preKSPSolveActions()

void bitpit::SystemSolver::preKSPSolveActions ( )
protectedvirtual

Pre-solve actions.

Definition at line 1300 of file system_solvers_large.cpp.

◆ prepareKSP()

void bitpit::SystemSolver::prepareKSP ( )
protectedvirtual

Prepare the KSP before the solution of the system.

Definition at line 3170 of file system_solvers_large.cpp.

◆ prePreconditionerSetupActions()

void bitpit::SystemSolver::prePreconditionerSetupActions ( )
protectedvirtual

Perform actions before preconditioner setup.

Definition at line 3343 of file system_solvers_large.cpp.

◆ reorderVector()

void bitpit::SystemSolver::reorderVector ( Vec vector,
IS permutations,
bool invert ) const
protected

Reorder the specified vector.

Parameters
invertis a flag for inverting the ordering

Definition at line 2919 of file system_solvers_large.cpp.

◆ resetKSPOptions()

void bitpit::SystemSolver::resetKSPOptions ( KSPOptions * options) const
protectedvirtual

Reset the specified KSP options.

Parameters
optionsare the options that will be reset

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3457 of file system_solvers_large.cpp.

◆ resetKSPStatus() [1/2]

void bitpit::SystemSolver::resetKSPStatus ( )
protectedvirtual

Reset the status of the KSP.

Reimplemented in bitpit::SplitSystemSolver, and bitpit::SplitSystemSolver.

Definition at line 3518 of file system_solvers_large.cpp.

◆ resetKSPStatus() [2/2]

void bitpit::SystemSolver::resetKSPStatus ( KSPStatus * status) const
protectedvirtual

Reset the status of the specified KSP.

Parameters
statuson output will contain the status of the KSP

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3528 of file system_solvers_large.cpp.

◆ restoreInfo()

void bitpit::SystemSolver::restoreInfo ( std::istream & systemStream)
protectedvirtual

Restore system information.

Parameters
systemStreamis the stream from which system information is read

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 2245 of file system_solvers_large.cpp.

◆ restoreMatrix()

void bitpit::SystemSolver::restoreMatrix ( const std::string & directory,
const std::string & name,
bool redistribute,
Mat * matrix ) const
protected

Restore the specified matrix.

Parameters
directoryis the directory from which the matrix data file will be read
nameis the name of the matrix that will be dumped
redistributeif 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
[out]matrixon output will contain the restored matrix

Definition at line 2465 of file system_solvers_large.cpp.

◆ restoreRHSRawPtr()

void bitpit::SystemSolver::restoreRHSRawPtr ( double * raw_rhs)

Restores the solution vector after getRHSRawPtr() has been called.

Parameters
raw_rhsis the location of pointer to array obtained from getRHSRawPtr()
Examples
RBF_example_00001.cpp.

Definition at line 1899 of file system_solvers_large.cpp.

◆ restoreRHSRawReadPtr()

void bitpit::SystemSolver::restoreRHSRawReadPtr ( const double * raw_rhs) const

Restores the solution vector after getRHSRawReadPtr() has been called.

Parameters
raw_rhsis the location of pointer to array obtained from getRHSRawReadPtr()

Definition at line 1910 of file system_solvers_large.cpp.

◆ restoreSolutionRawPtr()

void bitpit::SystemSolver::restoreSolutionRawPtr ( double * raw_solution)

Restores the solution vector after getSolutionRawPtr() has been called.

Parameters
raw_solutionis the location of pointer to array obtained from getSolutionRawPtr()
Examples
RBF_example_00001.cpp.

Definition at line 1957 of file system_solvers_large.cpp.

◆ restoreSolutionRawReadPtr()

void bitpit::SystemSolver::restoreSolutionRawReadPtr ( const double * raw_solution) const

Restores the solution vector after getSolutionRawReadPtr() has been called.

Parameters
raw_solutionis the location of pointer to array obtained from getSolutionRawReadPtr()

Definition at line 1968 of file system_solvers_large.cpp.

◆ restoreSystem()

void bitpit::SystemSolver::restoreSystem ( MPI_Comm communicator,
bool redistribute,
const std::string & directory,
const std::string & prefix = "" )

Restore the system from files.

Only the contents of the system will be restored, this include the matrix, the RHS vector, the solution vector, and the information needed to properly create the aforementioned data structures (e.g., the block size or the transpose flag).

Matrix file is mandatory, whereas solution and RHS files are not. If no solution or RHS files are found they will be re-created from scratch (but they will not be initialized).

Parameters
communicatoris the MPI communicator
redistributeif set to true, the system will be redistributed among the available processes, allowing to restore the system with a different number of processes than those used to dump it
directoryis the directory where the files will be read from
prefixis the prefix that will be was added to the files during

Definition at line 2180 of file system_solvers_large.cpp.

◆ restoreVector() [1/3]

void bitpit::SystemSolver::restoreVector ( const std::string & directory,
const std::string & name,
bool redistribute,
Vec * vector ) const
protected

Restore the specified vector.

Parameters
directoryis the directory from which the vector data file will be read
nameis the name of the vector that will be dumped
redistributeif set to true, the vector will be redistributed among the available processes, allowing to restore the vector with a different number of processes than those used to dump it
[out]vectoron output will contain the restored vector

Definition at line 2828 of file system_solvers_large.cpp.

◆ restoreVector() [2/3]

void bitpit::SystemSolver::restoreVector ( const std::string & directory,
const std::string & name,
Mat matrix,
VectorSide side,
Vec * vector ) const
protected

Restore the specified vector.

Parameters
directoryis the directory from which the vector data file will be read
nameis the name of the vector that will be dumped
matrixis the matrix the vectors should be compatible with
sidespecifies whether the vector can be multiplied by the matrix, or whether the matrix-vector product can be stored in it
[out]vectoron output will contain the restored vector

Definition at line 2775 of file system_solvers_large.cpp.

◆ restoreVector() [3/3]

void bitpit::SystemSolver::restoreVector ( const std::string & directory,
const std::string & name,
std::size_t localSize,
std::size_t globalSize,
Vec * vector ) const
protected

Restore the specified vector.

Parameters
directoryis the directory from which the vector data file will be read
nameis the name of the vector that will be dumped
localSizeis the local size of the restored vector
globalSizeis the global size of the restored vector
[out]vectoron output will contain the restored vector

Definition at line 2868 of file system_solvers_large.cpp.

◆ setNullSpace()

void bitpit::SystemSolver::setNullSpace ( )
virtual

Attaches a null space to the system matrix.

Definition at line 3059 of file system_solvers_large.cpp.

◆ setReordering()

void bitpit::SystemSolver::setReordering ( long nRows,
long nCols,
const SystemMatrixOrdering & reordering )
protected

Set the reordering that will be applied when assembling the matrix.

Reordering will be applied when the system is assembled and its sole purpose is to speed up the resolution of the system (e.g., reorder can be used to reduce the fill-in of the LU factorization).

Reordering is only applied internally, all public functions expects row and column indices to be in natural matrix order (i.e., not reordered).

If the system is partitioned, each process can reorder only it's local part of the matrix.

Parameters
nRowsis the number of rows of the system matrix
nColsis the number of columns of the system matrix
reorderingis the reordering that will be applied

Definition at line 3096 of file system_solvers_large.cpp.

◆ setTranspose()

void bitpit::SystemSolver::setTranspose ( bool transpose)

Set the transpose flag.

If the system is already assembled and the transpose flag needs to be changed, both the solution vector and the RHS one will be destroyed and re-created.

If the transpose flag needs to be changed, the workspaces associated wit the system will be cleared.

Parameters
transposeif set to true, transposed system will be solved

Definition at line 181 of file system_solvers_large.cpp.

◆ setUp()

void bitpit::SystemSolver::setUp ( )

Setup the system.

This function is deprecated, there is no need to call it.

Definition at line 3162 of file system_solvers_large.cpp.

◆ setupKrylov() [1/2]

void bitpit::SystemSolver::setupKrylov ( )
protectedvirtual

Set up the Krylov subspace method used to solve the system.

Reimplemented in bitpit::SplitSystemSolver, and bitpit::SplitSystemSolver.

Definition at line 3359 of file system_solvers_large.cpp.

◆ setupKrylov() [2/2]

void bitpit::SystemSolver::setupKrylov ( KSP ksp,
const KSPOptions & options ) const
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
kspis the KSP whose Krylov subspace method will be setup
optionsare the options that will be used to set up the KSP

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 3374 of file system_solvers_large.cpp.

◆ setupPreconditioner() [1/2]

void bitpit::SystemSolver::setupPreconditioner ( )
protectedvirtual

Set up the preconditioner.

Reimplemented in bitpit::SplitSystemSolver, and bitpit::SplitSystemSolver.

Examples
LA_example_00001.cpp.

Definition at line 3267 of file system_solvers_large.cpp.

◆ setupPreconditioner() [2/2]

void bitpit::SystemSolver::setupPreconditioner ( PC pc,
const KSPOptions & options ) const
protectedvirtual

Set up the specified preconditioner using the given options.

Parameters
pcis the preconditioner to set up
optionsare the options that will be used to set up the preconditioner

Reimplemented in bitpit::SplitSystemSolver, AMGSystemSolver, and AMGSplitSystemSolver.

Definition at line 3280 of file system_solvers_large.cpp.

◆ solve() [1/2]

void bitpit::SystemSolver::solve ( )

Solve the system

Examples
RBF_example_00001.cpp.

Definition at line 1241 of file system_solvers_large.cpp.

◆ solve() [2/2]

void bitpit::SystemSolver::solve ( const std::vector< double > & rhs,
std::vector< double > * solution )

Solve the system

Parameters
rhsis the right-hand-side of the system
solutionin input should contain the initial solution, on output it contains the solution of the linear system

Definition at line 1265 of file system_solvers_large.cpp.

◆ solveKSP()

void bitpit::SystemSolver::solveKSP ( )
protectedvirtual

Solve KSP.

Definition at line 1280 of file system_solvers_large.cpp.

◆ unsetNullSpace()

void bitpit::SystemSolver::unsetNullSpace ( )

Removes the null space from the system matrix.

Definition at line 3074 of file system_solvers_large.cpp.

◆ update() [1/4]

void bitpit::SystemSolver::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
assembleris the matrix assembler for the rows that will be updated

Definition at line 997 of file system_solvers_large.cpp.

◆ update() [2/4]

void bitpit::SystemSolver::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
nRowsis the number of rows that will be updated
rowsare the indices of the rows that will be updated
assembleris the matrix assembler for the rows that will be updated

Definition at line 1012 of file system_solvers_large.cpp.

◆ update() [3/4]

void bitpit::SystemSolver::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
nRowsis the number of rows that will be updated
rowsare the indices of the rows that will be updated
elementsare the elements that will be used to update the rows

Definition at line 977 of file system_solvers_large.cpp.

◆ update() [4/4]

template<typename DerivedSystemSolver >
void bitpit::SystemSolver::update ( long nRows,
const long * rows,
const typename DerivedSystemSolver::Assembler & assembler )
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
nRowsis the number of rows that will be updated
rowsare the indices of the rows that will be updated
assembleris the matrix assembler for the rows that will be updated

Definition at line 139 of file system_solvers_large.tpp.

◆ vectorsCreate()

void bitpit::SystemSolver::vectorsCreate ( )
protectedvirtual

Create RHS and solution vectors.

Vectors will be created, but they will not be initialized.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1753 of file system_solvers_large.cpp.

◆ vectorsDestroy()

void bitpit::SystemSolver::vectorsDestroy ( )
protectedvirtual

Destroy RHS and solution vectors.

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1851 of file system_solvers_large.cpp.

◆ vectorsDump()

void bitpit::SystemSolver::vectorsDump ( std::ostream & systemStream,
const std::string & directory,
const std::string & prefix ) const
protectedvirtual

Dump RHS and solution vectors.

Parameters
systemStreamis the stream from which system information is read
directoryis the directory in which the vectors data file will be written
prefixis the prefix added to the names of the files containing vector data

Definition at line 1823 of file system_solvers_large.cpp.

◆ vectorsFill() [1/2]

void bitpit::SystemSolver::vectorsFill ( const std::string & rhsFilePath,
const std::string & solutionFilePath )
protectedvirtual

Fill RHS and solution vectors reading their contents from the specified file.

Parameters
rhsFilePathis the file path containing the content of the RHS vector, if the path is empty the function will leave the RHS vector unaltered
solutionFilePathis the file path containing the content of the solution vector, if the path is empty the function will leave the solution vector unaltered

Definition at line 1789 of file system_solvers_large.cpp.

◆ vectorsFill() [2/2]

void bitpit::SystemSolver::vectorsFill ( const std::vector< double > & rhs,
const std::vector< double > & solution )
protectedvirtual

Fill RHS and solution vectors using the given data.

Parameters
rhscontains the that that will be copied into the RHS vector, if the vector is empty no data will be copied and function will leave the RHS vector unaltered
solutioncontains the that that will be copied into the solution vector, if the vector is empty no data will be copied and function will leave the RHS vector unaltered

Definition at line 1770 of file system_solvers_large.cpp.

◆ vectorsReorder()

void bitpit::SystemSolver::vectorsReorder ( bool invert)
protectedvirtual

Reorder RHS and solution vectors to match the order of the system matrix.

Parameters
invertis a flag for inverting the ordering

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1805 of file system_solvers_large.cpp.

◆ vectorsRestore()

void bitpit::SystemSolver::vectorsRestore ( std::istream & systemStream,
const std::string & directory,
const std::string & prefix )
protectedvirtual

Restore RHS and solution vectors.

Parameters
systemStreamis the stream that contains system information
directoryis the directory from which the vector data file will be read
prefixis the prefix added to the name of the file containing vectors data

Reimplemented in bitpit::SplitSystemSolver.

Definition at line 1839 of file system_solvers_large.cpp.

Member Data Documentation

◆ m_A

Mat bitpit::SystemSolver::m_A
protected

Definition at line 369 of file system_solvers_large.hpp.

◆ m_colReordering

IS bitpit::SystemSolver::m_colReordering
protected

Definition at line 374 of file system_solvers_large.hpp.

◆ m_convergenceMonitorEnabled

bool bitpit::SystemSolver::m_convergenceMonitorEnabled
protected

Definition at line 376 of file system_solvers_large.hpp.

◆ m_flatten

bool bitpit::SystemSolver::m_flatten
protected

Definition at line 366 of file system_solvers_large.hpp.

◆ m_KSP

KSP bitpit::SystemSolver::m_KSP
protected

Definition at line 378 of file system_solvers_large.hpp.

◆ m_KSPDirty

bool bitpit::SystemSolver::m_KSPDirty
protected

Definition at line 379 of file system_solvers_large.hpp.

◆ m_KSPOptions

KSPOptions bitpit::SystemSolver::m_KSPOptions
protected

Definition at line 380 of file system_solvers_large.hpp.

◆ m_KSPStatus

KSPStatus bitpit::SystemSolver::m_KSPStatus
protected

Definition at line 381 of file system_solvers_large.hpp.

◆ m_rhs

Vec bitpit::SystemSolver::m_rhs
protected

Definition at line 370 of file system_solvers_large.hpp.

◆ m_rowReordering

IS bitpit::SystemSolver::m_rowReordering
protected

Definition at line 373 of file system_solvers_large.hpp.

◆ m_solution

Vec bitpit::SystemSolver::m_solution
protected

Definition at line 371 of file system_solvers_large.hpp.

◆ m_transpose

bool bitpit::SystemSolver::m_transpose
protected

Definition at line 367 of file system_solvers_large.hpp.


The documentation for this class was generated from the following files:
--- layout: doxygen_footer ---