oomph::DoubleMatrixBase Class Referenceabstract

#include <matrices.h>

+ Inheritance diagram for oomph::DoubleMatrixBase:

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
 
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)
 
void solve (const DoubleVector &rhs, DoubleVector &soln)
 
void solve (Vector< double > &rhs)
 
void solve (const Vector< double > &rhs, Vector< double > &soln)
 
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)
 
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

LinearSolverLinear_solver_pt
 
LinearSolverDefault_linear_solver_pt
 

Detailed Description

Abstract base class for matrices of doubles – adds abstract interfaces for solving, LU decomposition and multiplication by vectors.

Constructor & Destructor Documentation

◆ DoubleMatrixBase() [1/2]

oomph::DoubleMatrixBase::DoubleMatrixBase ( )
inline

(Empty) constructor.

LinearSolver * Default_linear_solver_pt
Definition: matrices.h:267
LinearSolver * Linear_solver_pt
Definition: matrices.h:264

◆ DoubleMatrixBase() [2/2]

oomph::DoubleMatrixBase::DoubleMatrixBase ( const DoubleMatrixBase matrix)
delete

Broken copy constructor.

◆ ~DoubleMatrixBase()

virtual oomph::DoubleMatrixBase::~DoubleMatrixBase ( )
inlinevirtual

virtual (empty) destructor

286 {}

Member Function Documentation

◆ linear_solver_pt() [1/2]

LinearSolver*& oomph::DoubleMatrixBase::linear_solver_pt ( )
inline

◆ linear_solver_pt() [2/2]

LinearSolver* const& oomph::DoubleMatrixBase::linear_solver_pt ( ) const
inline

Return a pointer to the linear solver object (const version)

303  {
304  return Linear_solver_pt;
305  }

References Linear_solver_pt.

◆ max_residual()

virtual double oomph::DoubleMatrixBase::max_residual ( const DoubleVector x,
const DoubleVector rhs 
)
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"

349  {
351  residual(x, rhs, res);
352  return res.max();
353  }
std::vector< double > DoubleVector
loads clump configuration
Definition: ClumpInput.h:26
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
virtual void residual(const DoubleVector &x, const DoubleVector &b, DoubleVector &residual_)
Find the residual, i.e. r=b-Ax the residual.
Definition: matrices.h:326
list x
Definition: plotDoE.py:28

References res, residual(), and plotDoE::x.

◆ multiply()

◆ multiply_transpose()

virtual void oomph::DoubleMatrixBase::multiply_transpose ( const DoubleVector x,
DoubleVector soln 
) const
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.

◆ ncol()

◆ nrow()

◆ operator()()

virtual double oomph::DoubleMatrixBase::operator() ( const unsigned long &  i,
const unsigned long &  j 
) const
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.

◆ operator=()

void oomph::DoubleMatrixBase::operator= ( const DoubleMatrixBase )
delete

Broken assignment operator.

◆ residual()

virtual void oomph::DoubleMatrixBase::residual ( const DoubleVector x,
const DoubleVector b,
DoubleVector residual_ 
)
inlinevirtual

Find the residual, i.e. r=b-Ax the residual.

329  {
330  // compute residual = Ax
331  this->multiply(x, residual_);
332 
333  // set residual to -residual (-Ax)
334  unsigned nrow_local = residual_.nrow_local();
335  double* residual_pt = residual_.values_pt();
336  for (unsigned i = 0; i < nrow_local; i++)
337  {
338  residual_pt[i] = -residual_pt[i];
339  }
340 
341  // set residual = b + residuals
342  residual_ += b;
343  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Scalar * b
Definition: benchVecAdd.cpp:17
virtual void multiply(const DoubleVector &x, DoubleVector &soln) const =0
Multiply the matrix by the vector x: soln=Ax.

References b, 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().

◆ solve() [1/4]

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

75  {
76 #ifdef PARANOID
77  if (Linear_solver_pt == 0)
78  {
79  throw OomphLibError("Linear_solver_pt not set in matrix",
82  }
83 #endif
84  // Use the linear algebra interface to the linear solver
85  Linear_solver_pt->solve(this, rhs, soln);
86  }
virtual void solve(Problem *const &problem_pt, DoubleVector &result)=0
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Linear_solver_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::LinearSolver::solve().

◆ solve() [2/4]

void oomph::DoubleMatrixBase::solve ( const Vector< double > &  rhs,
Vector< double > &  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

117  {
118 #ifdef PARANOID
119  if (Linear_solver_pt == 0)
120  {
121  throw OomphLibError("Linear_solver_pt not set in matrix",
124  }
125 #endif
126  // Use the linear algebra interface to the linear solver
127  Linear_solver_pt->solve(this, rhs, soln);
128  }

References Linear_solver_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::LinearSolver::solve().

◆ solve() [3/4]

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.

51  {
52 #ifdef PARANOID
53  if (Linear_solver_pt == 0)
54  {
55  throw OomphLibError("Linear_solver_pt not set in matrix",
58  }
59 #endif
60 
61  // Copy rhs vector into local storage so it doesn't get overwritten
62  // if the linear solver decides to initialise the solution vector, say,
63  // which it's quite entitled to do!
64  DoubleVector actual_rhs(rhs);
65 
66  // Use the linear algebra interface to the linear solver
67  Linear_solver_pt->solve(this, actual_rhs, rhs);
68  }

References Linear_solver_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::LinearSolver::solve().

Referenced by oomph::BlackBoxFDNewtonSolver::black_box_fd_newton_solve(), oomph::HelmholtzMGPreconditioner< DIM >::direct_solve(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::find_distance_to_free_surface(), oomph::PVDEquationsBase< DIM >::get_principal_stress(), oomph::FiniteElement::locate_zeta(), and oomph::InterpolateFromIntegralPointsBase::pre_compute_ipt_to_node_mapping().

◆ solve() [4/4]

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.

93  {
94 #ifdef PARANOID
95  if (Linear_solver_pt == 0)
96  {
97  throw OomphLibError("Linear_solver_pt not set in matrix",
100  }
101 #endif
102 
103  // Copy rhs vector into local storage so it doesn't get overwritten
104  // if the linear solver decides to initialise the solution vector, say,
105  // which it's quite entitled to do!
106  Vector<double> actual_rhs(rhs);
107 
108  // Use the linear algebra interface to the linear solver
109  Linear_solver_pt->solve(this, actual_rhs, rhs);
110  }

References Linear_solver_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::LinearSolver::solve().

Member Data Documentation

◆ Default_linear_solver_pt

◆ Linear_solver_pt


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