oomph::ComplexMatrixBase Class Referenceabstract

#include <complex_matrices.h>

+ Inheritance diagram for oomph::ComplexMatrixBase:

Public Member Functions

 ComplexMatrixBase ()
 (Empty) constructor. More...
 
 ComplexMatrixBase (const ComplexMatrixBase &matrix)=delete
 Broken copy constructor. More...
 
void operator= (const ComplexMatrixBase &)=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 ~ComplexMatrixBase ()
 virtual (empty) destructor More...
 
virtual std::complex< doubleoperator() (const unsigned long &i, const unsigned long &j) const =0
 
virtual int ludecompose ()
 
virtual void lubksub (Vector< std::complex< double >> &rhs)
 
virtual void solve (Vector< std::complex< double >> &rhs)
 
virtual void solve (const Vector< std::complex< double >> &rhs, Vector< std::complex< double >> &soln)
 
virtual void residual (const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &b, Vector< std::complex< double >> &residual)=0
 Find the residual, i.e. r=b-Ax the residual. More...
 
virtual double max_residual (const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &rhs)
 
virtual void multiply (const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)=0
 Multiply the matrix by the vector x: soln=Ax. More...
 
virtual void multiply_transpose (const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)=0
 Multiply the transposed matrix by the vector x: soln=A^T x. More...
 

Detailed Description

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

Constructor & Destructor Documentation

◆ ComplexMatrixBase() [1/2]

oomph::ComplexMatrixBase::ComplexMatrixBase ( )
inline

(Empty) constructor.

59 {}

◆ ComplexMatrixBase() [2/2]

oomph::ComplexMatrixBase::ComplexMatrixBase ( const ComplexMatrixBase matrix)
delete

Broken copy constructor.

◆ ~ComplexMatrixBase()

virtual oomph::ComplexMatrixBase::~ComplexMatrixBase ( )
inlinevirtual

virtual (empty) destructor

74 {}

Member Function Documentation

◆ lubksub()

virtual void oomph::ComplexMatrixBase::lubksub ( Vector< std::complex< double >> &  rhs)
inlinevirtual

LU backsubstitue a previously LU-decomposed matrix; The passed rhs will be over-written with the solution vector

Reimplemented in oomph::CCComplexMatrix, oomph::CRComplexMatrix, and oomph::DenseComplexMatrix.

98  {
99  throw OomphLibError(
100  "lubksub() has not been written for this matrix class\n",
103  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by solve().

◆ ludecompose()

virtual int oomph::ComplexMatrixBase::ludecompose ( )
inlinevirtual

LU decomposition of the matrix using the appropriate linear solver. Return the sign of the determinant

Dummy return

Reimplemented in oomph::CCComplexMatrix, oomph::CRComplexMatrix, and oomph::DenseComplexMatrix.

85  {
86  throw OomphLibError(
87  "ludecompose() has not been written for this matrix class\n",
90 
92  return 1;
93  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by solve().

◆ max_residual()

virtual double oomph::ComplexMatrixBase::max_residual ( const Vector< std::complex< double >> &  x,
const Vector< std::complex< double >> &  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"

126  {
127  unsigned long n = rhs.size();
128  Vector<std::complex<double>> res(n);
129  residual(x, rhs, res);
130  double ans = 0.0;
131  for (unsigned long i = 0; i < n; i++)
132  {
133  ans = std::max(ans, std::abs(res[i]));
134  }
135  return ans;
136  }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
virtual void residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &b, Vector< std::complex< double >> &residual)=0
Find the residual, i.e. r=b-Ax the residual.
#define max(a, b)
Definition: datatypes.h:23
list x
Definition: plotDoE.py:28

References abs(), i, max, n, res, residual(), and plotDoE::x.

◆ multiply()

virtual void oomph::ComplexMatrixBase::multiply ( const Vector< std::complex< double >> &  x,
Vector< std::complex< double >> &  soln 
)
pure virtual

Multiply the matrix by the vector x: soln=Ax.

Implemented in oomph::CCComplexMatrix, oomph::CRComplexMatrix, and oomph::DenseComplexMatrix.

◆ multiply_transpose()

virtual void oomph::ComplexMatrixBase::multiply_transpose ( const Vector< std::complex< double >> &  x,
Vector< std::complex< double >> &  soln 
)
pure virtual

Multiply the transposed matrix by the vector x: soln=A^T x.

Implemented in oomph::CCComplexMatrix, oomph::CRComplexMatrix, and oomph::DenseComplexMatrix.

◆ ncol()

virtual unsigned long oomph::ComplexMatrixBase::ncol ( ) const
pure virtual

Return the number of columns of the matrix.

Implemented in oomph::CCComplexMatrix, oomph::CRComplexMatrix, and oomph::DenseComplexMatrix.

Referenced by solve().

◆ nrow()

virtual unsigned long oomph::ComplexMatrixBase::nrow ( ) const
pure virtual

Return the number of rows of the matrix.

Implemented in oomph::CCComplexMatrix, oomph::CRComplexMatrix, and oomph::DenseComplexMatrix.

Referenced by solve().

◆ operator()()

virtual std::complex<double> oomph::ComplexMatrixBase::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::CCComplexMatrix, oomph::CRComplexMatrix, and oomph::DenseComplexMatrix.

◆ operator=()

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

Broken assignment operator.

◆ residual()

virtual void oomph::ComplexMatrixBase::residual ( const Vector< std::complex< double >> &  x,
const Vector< std::complex< double >> &  b,
Vector< std::complex< double >> &  residual 
)
pure virtual

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

Implemented in oomph::DenseComplexMatrix, oomph::CCComplexMatrix, and oomph::CRComplexMatrix.

Referenced by max_residual().

◆ solve() [1/2]

void oomph::ComplexMatrixBase::solve ( const Vector< std::complex< double >> &  rhs,
Vector< std::complex< double >> &  soln 
)
virtual

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

73  {
74  // Set the solution vector equal to the rhs
75  // N.B. This won't work if we change to another vector format
76  soln = rhs;
77 
78  // Overwrite the solution vector (rhs is unchanged)
79  solve(soln);
80  }
virtual void solve(Vector< std::complex< double >> &rhs)
Definition: complex_matrices.cc:34

References solve().

◆ solve() [2/2]

void oomph::ComplexMatrixBase::solve ( Vector< std::complex< double >> &  rhs)
virtual

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.

35  {
36 #ifdef PARANOID
37  // Check Matrix is square
38  if (nrow() != ncol())
39  {
40  throw OomphLibError(
41  "This matrix is not square, the matrix MUST be square!",
44  }
45  // Check that size of rhs = nrow()
46  if (rhs.size() != nrow())
47  {
48  std::ostringstream error_message_stream;
49  error_message_stream << "The rhs vector is not the right size. It is "
50  << rhs.size() << ", it should be " << nrow()
51  << std::endl;
52 
53  throw OomphLibError(error_message_stream.str(),
56  }
57 #endif
58 
59  // Call the LU decomposition
60  ludecompose();
61 
62  // Call the back substitution
63  lubksub(rhs);
64  }
virtual int ludecompose()
Definition: complex_matrices.h:84
virtual unsigned long ncol() const =0
Return the number of columns of the matrix.
virtual void lubksub(Vector< std::complex< double >> &rhs)
Definition: complex_matrices.h:97
virtual unsigned long nrow() const =0
Return the number of rows of the matrix.

References lubksub(), ludecompose(), ncol(), nrow(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by solve().


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