30 #ifndef OOMPH_COMPLEX_MATRICES_HEADER
31 #define OOMPH_COMPLEX_MATRICES_HEADER
35 #include <oomph-lib-config.h>
68 virtual unsigned long nrow()
const = 0;
71 virtual unsigned long ncol()
const = 0;
79 virtual std::complex<double>
operator()(
const unsigned long&
i,
80 const unsigned long& j)
const = 0;
87 "ludecompose() has not been written for this matrix class\n",
88 OOMPH_CURRENT_FUNCTION,
89 OOMPH_EXCEPTION_LOCATION);
100 "lubksub() has not been written for this matrix class\n",
101 OOMPH_CURRENT_FUNCTION,
102 OOMPH_EXCEPTION_LOCATION);
109 virtual void solve(
Vector<std::complex<double>>& rhs);
113 virtual void solve(
const Vector<std::complex<double>>& rhs,
114 Vector<std::complex<double>>& soln);
118 const Vector<std::complex<double>>& b,
125 const Vector<std::complex<double>>& rhs)
127 unsigned long n = rhs.size();
131 for (
unsigned long i = 0;
i < n;
i++)
133 ans = std::max(ans, std::abs(res[
i]));
140 Vector<std::complex<double>>& soln) = 0;
144 Vector<std::complex<double>>& soln) = 0;
202 const unsigned long& m,
203 const std::complex<double>& initial_val)
204 :
DenseMatrix<std::complex<double>>(n, m, initial_val),
219 inline unsigned long nrow()
const
225 inline unsigned long ncol()
const
233 const unsigned long& j)
const
241 const unsigned long& j)
261 const Vector<std::complex<double>>& rhs,
266 Vector<std::complex<double>>& soln);
270 Vector<std::complex<double>>& soln);
307 const unsigned long& n,
308 const unsigned long& m)
349 inline unsigned long nrow()
const
355 inline unsigned long ncol()
const
362 const unsigned long& j)
const
378 const Vector<std::complex<double>>& b,
383 Vector<std::complex<double>>& soln);
388 Vector<std::complex<double>>& soln);
445 public CCMatrix<std::complex<double>>
475 const unsigned long& n,
476 const unsigned long& m)
512 inline unsigned long nrow()
const
518 inline unsigned long ncol()
const
526 const unsigned long& j)
const
542 const Vector<std::complex<double>>& b,
548 Vector<std::complex<double>>& soln);
553 Vector<std::complex<double>>& soln);
A class for compressed column matrices that store doubles.
virtual ~CCComplexMatrix()
Destructor: Kill the LU factors if they have been setup.
std::complex< double > operator()(const unsigned long &i, const unsigned long &j) const
Overload the round-bracket access operator to provide read-only (const) access to the data.
void multiply(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the matrix by the vector x: soln=Ax.
int Info
Info flag for the SuperLU solver.
void lubksub(Vector< std::complex< double >> &rhs)
LU back solve for given RHS.
unsigned long nrow() const
Return the number of rows of the matrix.
CCComplexMatrix()
Default constructor.
unsigned long ncol() const
Return the number of columns of the matrix.
int ludecompose()
LU decomposition using SuperLU.
void clean_up_memory()
LU clean up (perhaps this should happen in the destructor)
bool Doc_stats_during_solve
Flag to indicate if stats are to be displayed during solution of linear system with SuperLU.
void operator=(const CCComplexMatrix &)=delete
Broken assignment operator.
void residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &b, Vector< std::complex< double >> &residual)
Find the residulal to x of Ax=b, ie r=b-Ax.
CCComplexMatrix(const CCComplexMatrix &matrix)=delete
Broken copy constructor.
void * F_factors
Storage for the LU factors as required by SuperLU.
void multiply_transpose(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the transposed matrix by the vector x: soln=A^T x.
void disable_doc_stats()
the solve
CCComplexMatrix(const Vector< std::complex< double >> &value, const Vector< int > &row_index, const Vector< int > &column_start, const unsigned long &n, const unsigned long &m)
Constructor: Pass vector of values, vector of row indices, vector of column starts and number of rows...
void enable_doc_stats()
Set flag to indicate that stats are to be displayed during solution of linear system with SuperLU.
/////////////////////////////////////////////////////////////// /////////////////////////////////////...
int * column_start()
Access to C-style column_start array.
std::complex< double > get_entry(const unsigned long &i, const unsigned long &j) const
Access function that will be called by the read-only round-bracket operator (const)
int * row_index()
Access to C-style row index array.
A class for compressed row matrices.
unsigned long nrow() const
Return the number of rows of the matrix.
std::complex< double > operator()(const unsigned long &i, const unsigned long &j) const
Overload the round-bracket access operator for read-only access.
void * F_factors
Storage for the LU factors as required by SuperLU.
SerialMatrixMultiplyMethod
Serial matrix multiply method enumeration.
void enable_doc_stats()
Set flag to indicate that stats are to be displayed during solution of linear system with SuperLU.
void multiply(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the matrix by the vector x: soln=Ax.
virtual ~CRComplexMatrix()
Destructor: Kill the LU decomposition if it has been computed.
SerialMatrixMultiplyMethod & serial_matrix_matrix_multiply_method()
Access function to Serial_matrix_matrix_multiply_method, the flag which determines the matrix matrix ...
unsigned long ncol() const
Return the number of columns of the matrix.
bool Doc_stats_during_solve
Flag to indicate if stats are to be displayed during solution of linear system with SuperLU.
SerialMatrixMultiplyMethod Serial_matrix_matrix_multiply_method
A switch variable for selecting the matrix multiply method for serial (non-parallel) runs....
void add(const CRDoubleMatrix &matrix_in, CRComplexMatrix &result_matrix) const
Element-wise addition of this matrix with matrix_in.
SerialMatrixMultiplyMethod serial_matrix_matrix_multiply_method() const
Read only access function (const version) to Serial_matrix_matrix_multiply_method,...
int ludecompose()
LU decomposition using SuperLU.
CRComplexMatrix(const CRComplexMatrix &matrix)=delete
Broken copy constructor.
void disable_doc_stats()
the solve
CRComplexMatrix()
Default constructor.
void lubksub(Vector< std::complex< double >> &rhs)
LU back solve for given RHS.
void operator=(const CRComplexMatrix &)=delete
Broken assignment operator.
int Info
Info flag for the SuperLU solver.
void residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &b, Vector< std::complex< double >> &residual)
Find the residual to x of Ax=b, i.e. r=b-Ax.
void multiply_transpose(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the transposed matrix by the vector x: soln=A^T x.
CRComplexMatrix(const Vector< std::complex< double >> &value, const Vector< int > &column_index, const Vector< int > &row_start, const unsigned long &n, const unsigned long &m)
Constructor: Pass vector of values, vector of column indices, vector of row starts and number of colu...
void clean_up_memory()
LU clean up (perhaps this should happen in the destructor)
A class for compressed row matrices. This is a distributable object.
A class for compressed row matrices, a sparse storage format Once again the recursive template trick ...
int * row_start()
Access to C-style row_start array.
int * column_index()
Access to C-style column index array.
std::complex< double > get_entry(const unsigned long &i, const unsigned long &j) const
Access function that will be called by the read-only round-bracket operator (const)
Abstract base class for matrices of complex doubles – adds abstract interfaces for solving,...
virtual int ludecompose()
LU decomposition of the matrix using the appropriate linear solver. Return the sign of the determinan...
virtual ~ComplexMatrixBase()
virtual (empty) destructor
virtual double max_residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &rhs)
Find the maximum residual r=b-Ax – generic version, can be overloaded for specific derived classes wh...
virtual std::complex< 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 co...
virtual unsigned long ncol() const =0
Return the number of columns of the matrix.
virtual void multiply(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)=0
Multiply the matrix by the vector x: soln=Ax.
virtual void lubksub(Vector< std::complex< double >> &rhs)
LU backsubstitue a previously LU-decomposed matrix; The passed rhs will be over-written with the solu...
void operator=(const ComplexMatrixBase &)=delete
Broken assignment operator.
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.
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.
virtual unsigned long nrow() const =0
Return the number of rows of the matrix.
ComplexMatrixBase()
(Empty) constructor.
virtual void solve(Vector< std::complex< double >> &rhs)
Complete LU solve (replaces matrix by its LU decomposition and overwrites RHS with solution)....
ComplexMatrixBase(const ComplexMatrixBase &matrix)=delete
Broken copy constructor.
Class of matrices containing double complex, and stored as a DenseMatrix<complex<double> >,...
std::complex< double > & operator()(const unsigned long &i, const unsigned long &j)
Overload the non-const version of the round-bracket access operator for read-write access.
bool Overwrite_matrix_storage
Boolean flag used to decided whether the LU decomposition is stored separately, or not.
void lubksub(Vector< std::complex< double >> &rhs)
Overload the LU back substitution. Note that the rhs will be overwritten with the solution vector.
virtual ~DenseComplexMatrix()
Destructor, delete the storage for Index vector and LU storage (if any)
unsigned long ncol() const
Return the number of columns of the matrix.
Vector< long > * Index
Pointer to storage for the index of permutations in the LU solve.
void residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &rhs, Vector< std::complex< double >> &residual)
Find the residual of Ax=b, ie r=b-Ax for the "solution" x.
void multiply(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the matrix by the vector x: soln=Ax.
int ludecompose()
Overload the LU decomposition. For this storage scheme the matrix storage will be over-writeen by its...
DenseComplexMatrix(const unsigned long &n, const unsigned long &m, const std::complex< double > &initial_val)
Constructor to build a matrix with n rows and m columns, with initial value initial_val.
DenseComplexMatrix(const DenseComplexMatrix &matrix)=delete
Broken copy constructor.
void delete_lu_factors()
Function that deletes the storage for the LU_factors, if it is not the same as the matrix storage.
DenseComplexMatrix(const unsigned long &n, const unsigned long &m)
Constructor to build a matrix with n rows and m columns.
void operator=(const DenseComplexMatrix &)=delete
Broken assignment operator.
DenseComplexMatrix(const unsigned long &n)
Constructor to build a square n by n matrix.
DenseComplexMatrix()
Empty constructor, simply assign the lengths N and M to 0.
std::complex< double > operator()(const unsigned long &i, const unsigned long &j) const
Overload the const version of the round-bracket access operator for read-only access.
std::complex< double > * LU_factors
Pointer to storage for the LU decomposition.
unsigned long nrow() const
Return the number of rows of the matrix.
void multiply_transpose(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the transposed matrix by the vector x: soln=A^T x.
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
std::complex< double > & entry(const unsigned long &i, const unsigned long &j)
The access function that will be called by the read-write round-bracket operator.
std::complex< double > get_entry(const unsigned long &i, const unsigned long &j) const
The access function the will be called by the read-only (const version) round-bracket operator.
An OomphLibError object which should be thrown when an run-time error is encountered....
T * value()
Access to C-style value array.
A slight extension to the standard template vector class so that we can include "graceful" array rang...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...