Abstract base class for matrices of doubles – adds abstract interfaces for solving, LU decomposition and multiplication by vectors. More...
#include <matrices.h>
Public Member Functions | |
DoubleMatrixBase () | |
(Empty) constructor. More... | |
DoubleMatrixBase (const DoubleMatrixBase &matrix)=delete | |
Broken copy constructor. More... | |
void | operator= (const DoubleMatrixBase &)=delete |
Broken assignment operator. More... | |
virtual unsigned long | nrow () const =0 |
Return the number of rows of the matrix. More... | |
virtual unsigned long | ncol () const =0 |
Return the number of columns of the matrix. More... | |
virtual | ~DoubleMatrixBase () |
virtual (empty) destructor More... | |
virtual double | operator() (const unsigned long &i, const unsigned long &j) const =0 |
Round brackets to give access as a(i,j) for read only (we're not providing a general interface for component-wise write access since not all matrix formats allow efficient direct access!) More... | |
LinearSolver *& | linear_solver_pt () |
Return a pointer to the linear solver object. More... | |
LinearSolver *const & | linear_solver_pt () const |
Return a pointer to the linear solver object (const version) More... | |
void | solve (DoubleVector &rhs) |
Complete LU solve (replaces matrix by its LU decomposition and overwrites RHS with solution). The default should not need to be over-written. More... | |
void | solve (const DoubleVector &rhs, DoubleVector &soln) |
Complete LU solve (Nothing gets overwritten!). The default should not need to be overwritten. More... | |
void | solve (Vector< double > &rhs) |
Complete LU solve (replaces matrix by its LU decomposition and overwrites RHS with solution). The default should not need to be over-written. More... | |
void | solve (const Vector< double > &rhs, Vector< double > &soln) |
Complete LU solve (Nothing gets overwritten!). The default should not need to be overwritten. More... | |
virtual void | residual (const DoubleVector &x, const DoubleVector &b, DoubleVector &residual_) |
Find the residual, i.e. r=b-Ax the residual. More... | |
virtual double | max_residual (const DoubleVector &x, const DoubleVector &rhs) |
Find the maximum residual r=b-Ax – generic version, can be overloaded for specific derived classes where the max. can be determined "on the fly". More... | |
virtual void | multiply (const DoubleVector &x, DoubleVector &soln) const =0 |
Multiply the matrix by the vector x: soln=Ax. More... | |
virtual void | multiply_transpose (const DoubleVector &x, DoubleVector &soln) const =0 |
Multiply the transposed matrix by the vector x: soln=A^T x. More... | |
Protected Attributes | |
LinearSolver * | Linear_solver_pt |
LinearSolver * | Default_linear_solver_pt |
Abstract base class for matrices of doubles – adds abstract interfaces for solving, LU decomposition and multiplication by vectors.
Definition at line 260 of file matrices.h.
|
inline |
(Empty) constructor.
Definition at line 271 of file matrices.h.
|
delete |
Broken copy constructor.
|
inlinevirtual |
virtual (empty) destructor
Definition at line 286 of file matrices.h.
|
inline |
Return a pointer to the linear solver object.
Definition at line 296 of file matrices.h.
References Linear_solver_pt.
Referenced by oomph::HelmholtzMGPreconditioner< DIM >::disable_smoother_and_superlu_doc_time().
|
inline |
Return a pointer to the linear solver object (const version)
Definition at line 302 of file matrices.h.
References Linear_solver_pt.
|
inlinevirtual |
Find the maximum residual r=b-Ax – generic version, can be overloaded for specific derived classes where the max. can be determined "on the fly".
Definition at line 348 of file matrices.h.
References oomph::DoubleVector::max(), and residual().
|
pure virtual |
Multiply the matrix by the vector x: soln=Ax.
Implemented in oomph::SumOfMatrices, oomph::CCDoubleMatrix, oomph::DenseDoubleMatrix, and oomph::CRDoubleMatrix.
Referenced by oomph::SumOfMatrices::multiply(), residual(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), oomph::DampedJacobi< MATRIX >::solve_helper(), and oomph::GMRES< MATRIX >::solve_helper().
|
pure virtual |
Multiply the transposed matrix by the vector x: soln=A^T x.
Implemented in oomph::SumOfMatrices, oomph::CCDoubleMatrix, oomph::DenseDoubleMatrix, and oomph::CRDoubleMatrix.
|
pure virtual |
Return the number of columns of the matrix.
Implemented in oomph::SumOfMatrices, oomph::CCDoubleMatrix, oomph::DenseDoubleMatrix, and oomph::CRDoubleMatrix.
Referenced by oomph::SumOfMatrices::add_matrix(), oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::SuperLUSolver::factorise_serial(), oomph::SumOfMatrices::ncol(), oomph::HyprePreconditioner::setup(), oomph::TrilinosPreconditionerBase::setup(), oomph::DenseLU::solve(), oomph::SuperLUSolver::solve(), oomph::MumpsSolver::solve(), oomph::HypreSolver::solve(), oomph::TrilinosAztecOOSolver::solve(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), and oomph::SuperLUSolver::solve_transpose().
|
pure virtual |
Return the number of rows of the matrix.
Implemented in oomph::SumOfMatrices, oomph::CCDoubleMatrix, oomph::DenseDoubleMatrix, and oomph::CRDoubleMatrix.
Referenced by oomph::SumOfMatrices::add_matrix(), oomph::DenseLU::factorise(), oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::SuperLUSolver::factorise_serial(), oomph::SumOfMatrices::nrow(), oomph::MatrixBasedDiagPreconditioner::setup(), oomph::HyprePreconditioner::setup(), oomph::IdentityPreconditioner::setup(), oomph::TrilinosPreconditionerBase::setup(), oomph::DenseLU::solve(), oomph::SuperLUSolver::solve(), oomph::MumpsSolver::solve(), oomph::HypreSolver::solve(), oomph::TrilinosAztecOOSolver::solve(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), oomph::AugmentedProblemGMRES::solve_helper(), and oomph::SuperLUSolver::solve_transpose().
|
pure virtual |
Round brackets to give access as a(i,j) for read only (we're not providing a general interface for component-wise write access since not all matrix formats allow efficient direct access!)
Implemented in oomph::SumOfMatrices, oomph::CCDoubleMatrix, oomph::DenseDoubleMatrix, and oomph::CRDoubleMatrix.
|
delete |
Broken assignment operator.
|
inlinevirtual |
Find the residual, i.e. r=b-Ax the residual.
Definition at line 326 of file matrices.h.
References i, multiply(), oomph::DistributableLinearAlgebraObject::nrow_local(), and oomph::DoubleVector::values_pt().
Referenced by max_residual(), oomph::GMRESBlockPreconditioner::preconditioner_solve(), oomph::GS< MATRIX >::solve_helper(), oomph::GS< CRDoubleMatrix >::solve_helper(), and oomph::DampedJacobi< MATRIX >::solve_helper().
void oomph::DoubleMatrixBase::solve | ( | const DoubleVector & | rhs, |
DoubleVector & | soln | ||
) |
Complete LU solve (Nothing gets overwritten!). The default should not need to be overwritten.
Complete LU solve (Nothing gets overwritten!). This generic version should never need to be overwritten.
Definition at line 74 of file matrices.cc.
References Linear_solver_pt, and oomph::LinearSolver::solve().
Complete LU solve (Nothing gets overwritten!). The default should not need to be overwritten.
Complete LU solve (Nothing gets overwritten!). This generic version should never need to be overwritten.
Definition at line 116 of file matrices.cc.
References Linear_solver_pt, and oomph::LinearSolver::solve().
void oomph::DoubleMatrixBase::solve | ( | DoubleVector & | rhs | ) |
Complete LU solve (replaces matrix by its LU decomposition and overwrites RHS with solution). The default should not need to be over-written.
Complete LU solve (overwrites RHS with solution). This is the generic version which should not need to be over-written.
Definition at line 50 of file matrices.cc.
References Linear_solver_pt, and oomph::LinearSolver::solve().
Referenced by oomph::Newmark< NSTEPS >::assign_initial_data_values(), oomph::Newmark< NSTEPS >::assign_initial_data_values_stage2(), oomph::BlackBoxFDNewtonSolver::black_box_fd_newton_solve(), oomph::RefineableTetMeshBase::compute_volume_target(), oomph::HelmholtzMGPreconditioner< DIM >::direct_solve(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::find_distance_to_free_surface(), oomph::PVDEquationsBase< DIM >::get_principal_stress(), and oomph::FiniteElement::locate_zeta().
void oomph::DoubleMatrixBase::solve | ( | Vector< double > & | rhs | ) |
Complete LU solve (replaces matrix by its LU decomposition and overwrites RHS with solution). The default should not need to be over-written.
Complete LU solve (overwrites RHS with solution). This is the generic version which should not need to be over-written.
Definition at line 92 of file matrices.cc.
References Linear_solver_pt, and oomph::LinearSolver::solve().
|
protected |
Definition at line 267 of file matrices.h.
Referenced by oomph::CRDoubleMatrix::build(), oomph::CCDoubleMatrix::CCDoubleMatrix(), oomph::CRDoubleMatrix::CRDoubleMatrix(), oomph::DenseDoubleMatrix::DenseDoubleMatrix(), oomph::CRDoubleMatrix::lubksub(), oomph::DenseDoubleMatrix::lubksub(), oomph::CCDoubleMatrix::lubksub(), oomph::CRDoubleMatrix::ludecompose(), oomph::DenseDoubleMatrix::ludecompose(), oomph::CCDoubleMatrix::ludecompose(), oomph::CCDoubleMatrix::~CCDoubleMatrix(), oomph::CRDoubleMatrix::~CRDoubleMatrix(), and oomph::DenseDoubleMatrix::~DenseDoubleMatrix().
|
protected |
Definition at line 264 of file matrices.h.
Referenced by oomph::CCDoubleMatrix::CCDoubleMatrix(), oomph::CRDoubleMatrix::clear(), oomph::CRDoubleMatrix::CRDoubleMatrix(), oomph::DenseDoubleMatrix::DenseDoubleMatrix(), linear_solver_pt(), and solve().