25#ifndef __BITPIT_SYSTEM_SOLVERS_LARGE_HPP__
26#define __BITPIT_SYSTEM_SOLVERS_LARGE_HPP__
28#include <unordered_map>
29#include <unordered_set>
34#include "bitpit_IO.hpp"
36#include "system_matrix.hpp"
54 : overlap(PETSC_DEFAULT),
levels(PETSC_DEFAULT),
65 KSPConvergedReason convergence;
68 : error(0), its(-1), convergence(KSP_DIVERGED_BREAKDOWN)
119template<
typename RowRankStorage,
typename ColRankStorage>
130 const RowRankStorage *m_rowRankStorage;
131 const ColRankStorage *m_colRankStorage;
145#if BITPIT_ENABLE_MPI==1
146 virtual bool isPartitioned()
const = 0;
147 virtual const MPI_Comm & getCommunicator()
const = 0;
152 virtual int getBlockSize()
const = 0;
154 virtual long getRowCount()
const = 0;
155 virtual long getColCount()
const = 0;
157 virtual long getRowElementCount()
const = 0;
158 virtual long getColElementCount()
const = 0;
160#if BITPIT_ENABLE_MPI==1
161 virtual long getRowGlobalCount()
const = 0;
162 virtual long getColGlobalCount()
const = 0;
164 virtual long getRowGlobalElementCount()
const = 0;
165 virtual long getColGlobalElementCount()
const = 0;
167 virtual long getRowGlobalOffset()
const = 0;
168 virtual long getColGlobalOffset()
const = 0;
170 virtual long getRowGlobalElementOffset()
const = 0;
171 virtual long getColGlobalElementOffset()
const = 0;
174 virtual long getRowNZCount(
long rowIndex)
const = 0;
175 virtual long getMaxRowNZCount()
const = 0;
191#if BITPIT_ENABLE_MPI==1
206#if BITPIT_ENABLE_MPI==1
250 static PetscErrorCode displayLogView();
252 bool m_externalMPIInitialization;
253 bool m_externalPETScInitialization;
255 bool m_optionsEditable;
256 std::vector<std::string> m_options;
258 bool m_logViewEnabled;
260 void enableLogView();
277 SystemSolver(
const std::string &prefix,
bool debug =
false);
278 SystemSolver(
const std::string &prefix,
bool transpose,
bool debug);
279 SystemSolver(
const std::string &prefix,
bool flatten,
bool transpose,
bool debug);
283 virtual void clear();
294 void update(
long nRows,
const long *rows,
const SparseMatrix &elements);
296 void update(
long nRows,
const long *rows,
const Assembler &assembler);
310#if BITPIT_ENABLE_MPI==1
320 void solve(
const std::vector<double> &rhs, std::vector<double> *solution);
322 void dumpSystem(
const std::string &header,
const std::string &directory,
const std::string &prefix =
"")
const;
323#if BITPIT_ENABLE_MPI==1
324 void restoreSystem(MPI_Comm communicator,
bool redistribute,
const std::string &directory,
325 const std::string &prefix =
"");
327 void restoreSystem(
const std::string &directory,
const std::string &prefix =
"");
338 virtual void exportMatrix(
const std::string &filePath, FileFormat exportFormat = FILE_BINARY)
const;
346 virtual void exportRHS(
const std::string &filePath, FileFormat exportFormat = FILE_BINARY)
const;
347 virtual void importRHS(
const std::string &filePath);
354 virtual void exportSolution(
const std::string &filePath, FileFormat exportFormat = FILE_BINARY)
const;
376 bool m_convergenceMonitorEnabled;
387 virtual void matrixFill(
const std::string &filePath);
388 virtual void matrixDump(std::ostream &systemStream,
const std::string &directory,
const std::string &prefix)
const;
389#if BITPIT_ENABLE_MPI==1
390 virtual void matrixRestore(std::istream &systemStream,
const std::string &directory,
const std::string &prefix,
bool redistribute);
392 virtual void matrixRestore(std::istream &systemStream,
const std::string &directory,
const std::string &prefix);
397 virtual void vectorsFill(
const std::vector<double> &rhs,
const std::vector<double> &solution);
398 virtual void vectorsFill(
const std::string &rhsFilePath,
const std::string &solutionFilePath);
400 virtual void vectorsDump(std::ostream &systemStream,
const std::string &directory,
const std::string &prefix)
const;
401 virtual void vectorsRestore(std::istream &systemStream,
const std::string &directory,
const std::string &prefix);
407 template<
typename DerivedSystemSolver>
410 template<
typename DerivedSystemSolver>
411 void update(
long nRows,
const long *rows,
const typename DerivedSystemSolver::Assembler &assembler);
443 virtual void dumpInfo(std::ostream &stream)
const;
446 void createMatrix(
int rowBlockSize,
int colBlockSize, Mat *matrix)
const;
447 void createMatrix(
int rowBlockSize,
int colBlockSize,
int nNestRows,
int nNestCols, Mat *subMatrices, Mat *matrix)
const;
448 void fillMatrix(Mat matrix,
const std::string &filePath)
const;
449 void dumpMatrix(Mat matrix,
const std::string &directory,
const std::string &name)
const;
450#if BITPIT_ENABLE_MPI==1
451 void restoreMatrix(
const std::string &directory,
const std::string &name,
bool redistribute, Mat *matrix)
const;
453 void restoreMatrix(
const std::string &directory,
const std::string &name, Mat *matrix)
const;
455 void exportMatrix(Mat matrix,
const std::string &filePath, FileFormat fileFormat)
const;
459 void createVector(
int blockSize,
int nestSize, Vec *subVectors, Vec *vector)
const;
460 void fillVector(Vec vector,
const std::string &filePath)
const;
461 void fillVector(Vec vector,
const std::vector<double> &data)
const;
462 void reorderVector(Vec vector, IS permutations,
bool invert)
const;
463 void dumpVector(Vec vector,
const std::string &directory,
const std::string &name)
const;
464 void dumpVector(Vec vector, VectorSide side,
const std::string &directory,
const std::string &name)
const;
465 void restoreVector(
const std::string &directory,
const std::string &name, Mat matrix, VectorSide side, Vec *vector)
const;
466#if BITPIT_ENABLE_MPI==1
467 void restoreVector(
const std::string &directory,
const std::string &name,
bool redistribute, Vec *vector)
const;
469 void restoreVector(
const std::string &directory,
const std::string &name, Vec *vector)
const;
471 void restoreVector(
const std::string &directory,
const std::string &name, std::size_t localSize, std::size_t globalSize, Vec *vector)
const;
472 void exportVector(Vec vector,
const std::string &filePath, FileFormat fileFormat)
const;
473 void exportVector(Vec vector, std::vector<double> *data)
const;
476 std::string
getInfoFilePath(
const std::string &directory,
const std::string &name)
const;
477 std::string
getDataFilePath(
const std::string &directory,
const std::string &name)
const;
478 std::string
getFilePath(
const std::string &directory,
const std::string &name)
const;
480#if BITPIT_ENABLE_MPI==1
487 static int m_nInstances;
489 std::string m_prefix;
493#if BITPIT_ENABLE_MPI==1
494 MPI_Comm m_communicator;
499 bool m_forceConsistency;
501#if BITPIT_ENABLE_MPI==1
502 void setCommunicator(MPI_Comm communicator);
503 void freeCommunicator();
506 void resetPermutations();
508 void removeNullSpaceFromRHS();
510 void openBinaryArchive(
const std::string &directory,
const std::string &prefix,
512 void openBinaryArchive(
const std::string &header,
const std::string &directory,
513 const std::string &prefix,
int block,
OBinaryArchive *archive)
const;
517 int getBinaryArchiveBlock()
const;
523#include "system_solvers_large.tpp"
Base class for binary archives.
The NaturalSystemMatrixOrdering class defines allows to use a matrix natural ordering.
long getColPermutationRank(long col) const override
long getRowPermutationRank(long row) const override
The PetscManager class handles the interaction with PETSc library.
bool initialize(bool debug)
bool areOptionsEditable() const
void addInitOptions(int argc, char **args)
void addInitOption(const std::string &option)
The ProxySystemMatrixOrdering class defines allows to use a matrix ordering defined in some external ...
long getRowPermutationRank(long row) const override
long getColPermutationRank(long col) const override
ProxySystemMatrixOrdering(const RowRankStorage *rowRankStorage, const ColRankStorage *colRankStorage)
Metafunction for generating a list of elements that can be either stored in an external vectror or,...
The SystemMatrixAssembler class provides an interface for defining system matrix assemblers.
The SystemMatrixOrdering class provides an interface for defining classes that allows to reorder the ...
virtual long getColPermutationRank(long col) const =0
virtual long getRowPermutationRank(long row) const =0
The SystemSolver class provides methods for building and solving large linear systems.
void dumpMatrix(Mat matrix, const std::string &directory, const std::string &name) const
virtual void vectorsCreate()
void matrixUpdate(long nRows, const long *rows, const Assembler &assembler)
double * getSolutionRawPtr()
long getRowGlobalElementCount() const
bool isPartitioned() const
virtual int getDumpVersion() const
virtual void exportMatrix(const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
virtual void matrixFill(const std::string &filePath)
virtual void initializeKSPOptions()
virtual void prePreconditionerSetupActions()
void restoreSolutionRawPtr(double *raw_solution)
virtual void matrixRestore(std::istream &systemStream, const std::string &directory, const std::string &prefix, bool redistribute)
void fillMatrix(Mat matrix, const std::string &filePath) const
void restoreRHSRawReadPtr(const double *raw_rhs) const
void dumpVector(Vec vector, const std::string &directory, const std::string &name) const
virtual void postPreconditionerSetupActions()
virtual void importMatrix(const std::string &filePath)
virtual void initializeKSPStatus()
void reorderVector(Vec vector, IS permutations, bool invert) const
virtual void finalizeKSP()
std::string getInfoFilePath(const std::string &directory, const std::string &name) const
virtual void vectorsDump(std::ostream &systemStream, const std::string &directory, const std::string &prefix) const
virtual void resetKSPOptions(KSPOptions *options) const
void exportVector(Vec vector, const std::string &filePath, FileFormat fileFormat) const
void dumpSystem(const std::string &header, const std::string &directory, const std::string &prefix="") const
virtual void destroyKSPOptions()
virtual void clearWorkspace()
void setReordering(long nRows, long nCols, const SystemMatrixOrdering &reordering)
virtual void prepareKSP()
virtual void restoreInfo(std::istream &stream)
virtual void setNullSpace()
virtual void preKSPSolveActions()
virtual void setupKrylov()
virtual void postKrylovSetupActions()
long getColGlobalCount() const
void restoreSystem(MPI_Comm communicator, bool redistribute, const std::string &directory, const std::string &prefix="")
const double * getSolutionRawReadPtr() const
bool getTranspose() const
virtual void exportRHS(const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
virtual int getBlockSize() const
void restoreRHSRawPtr(double *raw_rhs)
void setTranspose(bool transpose)
const KSPStatus & getKSPStatus() const
void restoreMatrix(const std::string &directory, const std::string &name, bool redistribute, Mat *matrix) const
void matrixAssembly(const Assembler &assembler)
void createVector(int blockSize, Vec *vector) const
void destroyMatrix(Mat *matrix) const
std::string getFilePath(const std::string &directory, const std::string &name) const
long getColElementCount() const
virtual void matrixDump(std::ostream &systemStream, const std::string &directory, const std::string &prefix) const
virtual void dumpInfo(std::ostream &stream) const
void assembly(const SparseMatrix &matrix)
virtual void destroyKSP()
virtual void postKSPSolveActions()
const MPI_Comm & getCommunicator() const
virtual void fillKSPStatus()
std::string getDataFilePath(const std::string &directory, const std::string &name) const
const double * getRHSRawReadPtr() const
virtual void preKrylovSetupActions()
virtual void importSolution(const std::string &filePath)
virtual void resetKSPStatus()
SystemSolver(bool debug=false)
virtual void importRHS(const std::string &filePath)
void destroyVector(Vec *vector) const
void createMatrix(int rowBlockSize, int colBlockSize, Mat *matrix) const
void fillVector(Vec vector, const std::string &filePath) const
virtual void matrixDestroy()
virtual void destroyKSPStatus()
virtual void vectorsFill(const std::vector< double > &rhs, const std::vector< double > &solution)
virtual void vectorsDestroy()
virtual void setupPreconditioner()
void enableForceConsistency(bool enable)
long getRowGlobalCount() const
void restoreSolutionRawReadPtr(const double *raw_solution) const
bool isForceConsistencyEnabled() const
KSPOptions & getKSPOptions()
long getRowElementCount() const
virtual void vectorsReorder(bool invert)
virtual void vectorsRestore(std::istream &systemStream, const std::string &directory, const std::string &prefix)
long getColGlobalElementCount() const
void restoreVector(const std::string &directory, const std::string &name, Mat matrix, VectorSide side, Vec *vector) const
virtual void exportSolution(const std::string &filePath, FileFormat exportFormat=FILE_BINARY) const
The SystemSparseMatrixAssembler class defines an assembler for building the system matrix form a spar...
long getColGlobalCount() const override
long getRowGlobalCount() const override
void getRowValues(long rowIndex, ConstProxyVector< double > *values) const override
const MPI_Comm & getCommunicator() const override
long getRowElementCount() const override
long getColGlobalElementCount() const override
long getMaxRowNZCount() const override
long getRowGlobalElementCount() const override
void getRowData(long rowIndex, ConstProxyVector< long > *pattern, ConstProxyVector< double > *values) const override
long getRowNZCount(long rowIndex) const override
long getRowGlobalOffset() const override
long getRowCount() const override
bool isPartitioned() const override
AssemblyOptions getOptions() const override
long getColGlobalOffset() const override
long getColGlobalElementOffset() const override
void getRowPattern(long rowIndex, ConstProxyVector< long > *pattern) const override
long getColCount() const override
long getRowGlobalElementOffset() const override
long getColElementCount() const override
int getBlockSize() const override
SystemSparseMatrixAssembler(const SparseMatrix *matrix)
#define BITPIT_DEPRECATED(func)
PetscScalar subrtol
Deprecated, ASM ILU levels should be set using the "levels" member.
PetscInt levels
Overlap between a pair of subdomains for the partitioned preconditioner.
PetscInt restart
Tells the iterative solver that the initial guess is nonzero.
PetscInt maxits
Number of iterations at which the GMRES method restarts.
PetscScalar atol
Relative convergence tolerance, relative decrease in the preconditioned residual norm.
PetscInt sublevels
Absolute convergence tolerance, absolute size of the preconditioned residual norm.
PetscScalar rtol
Maximum number of iterations.
PetscBool initial_non_zero
Number of levels of fill used by the preconditioner.
KSPOptions()
Deprecated, it has never had any effect on the solution of the system.
bool full
Controls if the assembler is providing all the non-zero values of a row.