A class for compressed row matrices. This is a distributable object. More...
#include <matrices.h>
Classes | |
struct | CRDoubleMatrixComparisonHelper |
Create a struct to provide a comparison function for std::sort. More... | |
Public Member Functions | |
CRDoubleMatrix () | |
Default constructor. More... | |
CRDoubleMatrix (const LinearAlgebraDistribution *distribution_pt, const unsigned &ncol, const Vector< double > &value, const Vector< int > &column_index, const Vector< int > &row_start) | |
Constructor: vector of values, vector of column indices, vector of row starts and number of rows and columns. More... | |
CRDoubleMatrix (const LinearAlgebraDistribution *distribution_pt) | |
Constructor: just stores the distribution but does not build the matrix. More... | |
CRDoubleMatrix (const CRDoubleMatrix &matrix) | |
Copy constructor. More... | |
void | operator= (const CRDoubleMatrix &) |
Broken assignment operator. More... | |
virtual | ~CRDoubleMatrix () |
Destructor. More... | |
const Vector< int > | get_index_of_diagonal_entries () const |
Access function: returns the vector Index_of_diagonal_entries. The i-th entry of the vector contains the index of the last entry below or on the diagonal. If there are no entries below or on the diagonal then the corresponding entry is -1. If, however, there are no entries in the row then the entry is irrelevant and is kept as the initialised value; 0. More... | |
bool | entries_are_sorted (const bool &doc_unordered_entries=false) const |
Runs through the column index vector and checks if the entries follow the regular lexicographical ordering of matrix entries, i.e. it will check (at the i-th row of the matrix) if the entries in the column index vector associated with this row are in increasing order. More... | |
void | sort_entries () |
Sorts the entries associated with each row of the matrix in the column index vector and the value vector into ascending order and sets up the Index_of_diagonal_entries vector. More... | |
void | build (const LinearAlgebraDistribution *distribution_pt, const unsigned &ncol, const Vector< double > &value, const Vector< int > &column_index, const Vector< int > &row_start) |
build method: vector of values, vector of column indices, vector of row starts and number of rows and columns. More... | |
void | build (const LinearAlgebraDistribution *distribution_pt) |
rebuild the matrix - assembles an empty matrix will a defined distribution More... | |
void | build (const unsigned &ncol, const Vector< double > &value, const Vector< int > &column_index, const Vector< int > &row_start) |
keeps the existing distribution and just matrix that is stored More... | |
void | build_without_copy (const unsigned &ncol, const unsigned &nnz, double *value, int *column_index, int *row_start) |
keeps the existing distribution and just matrix that is stored without copying the matrix data More... | |
void | redistribute (const LinearAlgebraDistribution *const &dist_pt) |
void | clear () |
clear More... | |
unsigned long | nrow () const |
Return the number of rows of the matrix. More... | |
unsigned long | ncol () const |
Return the number of columns of the matrix. More... | |
void | output_bottom_right_zero_helper (std::ostream &outfile) const |
Output the "bottom right" entry regardless of it being zero or not (this allows automatic detection of matrix size in e.g. matlab, python). More... | |
void | sparse_indexed_output_helper (std::ostream &outfile) const |
Indexed output function to print a matrix to the stream outfile as i,j,a(i,j) for a(i,j)!=0 only. More... | |
void | sparse_indexed_output_with_offset (std::string filename) |
Indexed output function to print a matrix to a file as i,j,a(i,j) for a(i,j)!=0 only. Specify filename. This uses acual global row numbers. More... | |
double | operator() (const unsigned long &i, const unsigned long &j) const |
int * | row_start () |
Access to C-style row_start array. More... | |
const int * | row_start () const |
Access to C-style row_start array (const version) More... | |
int * | column_index () |
Access to C-style column index array. More... | |
const int * | column_index () const |
Access to C-style column index array (const version) More... | |
double * | value () |
Access to C-style value array. More... | |
const double * | value () const |
Access to C-style value array (const version) More... | |
unsigned long | nnz () const |
Return the number of nonzero entries (the local nnz) More... | |
virtual void | ludecompose () |
LU decomposition using SuperLU if matrix is not distributed or distributed onto a single processor. More... | |
virtual void | lubksub (DoubleVector &rhs) |
LU back solve for given RHS. More... | |
void | multiply (const DoubleVector &x, DoubleVector &soln) const |
Multiply the matrix by the vector x: soln=Ax. More... | |
void | multiply_transpose (const DoubleVector &x, DoubleVector &soln) const |
Multiply the transposed matrix by the vector x: soln=A^T x. More... | |
void | multiply (const CRDoubleMatrix &matrix_in, CRDoubleMatrix &result) const |
Function to multiply this matrix by the CRDoubleMatrix matrix_in. In a serial matrix, there are 4 methods available: Method 1: First runs through this matrix and matrix_in to find the storage requirements for result - arrays of the correct size are then allocated before performing the calculation. Minimises memory requirements but more costly. Method 2: Grows storage for values and column indices of result 'on the fly' using an array of maps. Faster but more memory intensive. Method 3: Grows storage for values and column indices of result 'on the fly' using a vector of vectors. Not particularly impressive on the platforms we tried... Method 4: Trilinos Epetra Matrix Matrix multiply. Method 5: Trilinox Epetra Matrix Matrix Mulitply (ml based) If Trilinos is installed then Method 4 is employed by default, otherwise Method 2 is employed by default. In a distributed matrix, only Trilinos Epetra Matrix Matrix multiply is available. More... | |
void | matrix_reduction (const double &alpha, CRDoubleMatrix &reduced_matrix) |
For every row, find the maximum absolute value of the entries in this row. Set all values that are less than alpha times this maximum to zero and return the resulting matrix in reduced_matrix. Note: Diagonal entries are retained regardless of their size. More... | |
unsigned & | serial_matrix_matrix_multiply_method () |
Access function to Serial_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for serial matrices. Method 1: First runs through this matrix and matrix_in to find the storage requirements for result - arrays of the correct size are then allocated before performing the calculation. Minimises memory requirements but more costly. Method 2: Grows storage for values and column indices of result 'on the fly' using an array of maps. Faster but more memory intensive. Method 3: Grows storage for values and column indices of result 'on the fly' using a vector of vectors. Not particularly impressive on the platforms we tried... Method 4: Trilinos Epetra Matrix Matrix multiply. Method 5: Trilinos Epetra Matrix Matrix multiply (ML based). More... | |
const unsigned & | serial_matrix_matrix_multiply_method () const |
Read only access function (const version) to Serial_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for serial matrices. Method 1: First runs through this matrix and matrix_in to find the storage requirements for result - arrays of the correct size are then allocated before performing the calculation. Minimises memory requirements but more costly. Method 2: Grows storage for values and column indices of result 'on the fly' using an array of maps. Faster but more memory intensive. Method 3: Grows storage for values and column indices of result 'on the fly' using a vector of vectors. Not particularly impressive on the platforms we tried... Method 4: Trilinos Epetra Matrix Matrix multiply. Method 5: Trilinos Epetra Matrix Matrix multiply (ML based). More... | |
unsigned & | distributed_matrix_matrix_multiply_method () |
Access function to Distributed_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for distributed matrices. Method 1: Trilinos Epetra Matrix Matrix multiply. Method 2: Trilinos Epetra Matrix Matrix multiply (ML based). More... | |
const unsigned & | distributed_matrix_matrix_multiply_method () const |
Read only access function (const version) to Distributed_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for distributed matrices. Method 1: Trilinos Epetra Matrix Matrix multiply. Method 2: Trilinos Epetra Matrix Matrix multiply (ML based). More... | |
bool | built () const |
access function to the Built flag - indicates whether the matrix has been build - i.e. the distribution has been defined and the matrix assembled. More... | |
CRDoubleMatrix * | global_matrix () const |
if this matrix is distributed then a the equivalent global matrix is built using new and returned. The calling method is responsible for the destruction of the new matrix. More... | |
void | get_matrix_transpose (CRDoubleMatrix *result) const |
Returns the transpose of this matrix. More... | |
double | inf_norm () const |
returns the inf-norm of this matrix More... | |
Vector< double > | diagonal_entries () const |
returns a Vector of diagonal entries of this matrix. This only works with square matrices. This condition may be relaxed in the future if need be. More... | |
void | add (const CRDoubleMatrix &matrix_in, CRDoubleMatrix &result_matrix) const |
element-wise addition of this matrix with matrix_in. More... | |
![]() | |
Matrix () | |
(Empty) constructor More... | |
Matrix (const Matrix &matrix) | |
Broken copy constructor. More... | |
void | operator= (const Matrix &) |
Broken assignment operator. More... | |
virtual | ~Matrix () |
Virtual (empty) destructor. More... | |
double | operator() (const unsigned long &i, const unsigned long &j) const |
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!) The function uses the MATRIX_TYPE template parameter to call the get_entry() function which must be defined in all derived classes that are to be fully instantiated. More... | |
double & | operator() (const unsigned long &i, const unsigned long &j) |
Round brackets to give access as a(i,j) for read-write access. The function uses the MATRIX_TYPE template parameter to call the entry() function which must be defined in all derived classes that are to be fully instantiated. If the particular Matrix does not allow write access, the function should break with an error message. More... | |
virtual void | output (std::ostream &outfile) const |
Output function to print a matrix row-by-row, in the form a(0,0) a(0,1) ... a(1,0) a(1,1) ... ... to the stream outfile. Broken virtual since it might not be sensible to implement this for some sparse matrices. More... | |
void | sparse_indexed_output (std::ostream &outfile, const unsigned &precision=0, const bool &output_bottom_right_zero=false) const |
Indexed output function to print a matrix to the stream outfile as i,j,a(i,j) for a(i,j)!=0 only with specified precision (if precision=0 then nothing is changed). If optional boolean flag is set to true we also output the "bottom right" entry regardless of it being zero or not (this allows automatic detection of matrix size in e.g. matlab, python). More... | |
void | sparse_indexed_output (std::string filename, const unsigned &precision=0, const bool &output_bottom_right_zero=false) const |
Indexed output function to print a matrix to the file named filename as i,j,a(i,j) for a(i,j)!=0 only with specified precision. If optional boolean flag is set to true we also output the "bottom right" entry regardless of it being zero or not (this allows automatic detection of matrix size in e.g. matlab, python). More... | |
![]() | |
DoubleMatrixBase () | |
(Empty) constructor. More... | |
DoubleMatrixBase (const DoubleMatrixBase &matrix) | |
Broken copy constructor. More... | |
void | operator= (const DoubleMatrixBase &) |
Broken assignment operator. More... | |
virtual | ~DoubleMatrixBase () |
virtual (empty) destructor 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... | |
![]() | |
DistributableLinearAlgebraObject () | |
Default constructor - create a distribution. More... | |
DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix) | |
Broken copy constructor. More... | |
void | operator= (const DistributableLinearAlgebraObject &) |
Broken assignment operator. More... | |
virtual | ~DistributableLinearAlgebraObject () |
Destructor. More... | |
LinearAlgebraDistribution * | distribution_pt () const |
access to the LinearAlgebraDistribution More... | |
unsigned | nrow () const |
access function to the number of global rows. More... | |
unsigned | nrow_local () const |
access function for the num of local rows on this processor. More... | |
unsigned | nrow_local (const unsigned &p) const |
access function for the num of local rows on this processor. More... | |
unsigned | first_row () const |
access function for the first row on this processor More... | |
unsigned | first_row (const unsigned &p) const |
access function for the first row on this processor More... | |
bool | distributed () const |
distribution is serial or distributed More... | |
bool | distribution_built () const |
void | build_distribution (const LinearAlgebraDistribution *const dist_pt) |
setup the distribution of this distributable linear algebra object More... | |
void | build_distribution (const LinearAlgebraDistribution &dist) |
setup the distribution of this distributable linear algebra object More... | |
Public Attributes | |
struct oomph::CRDoubleMatrix::CRDoubleMatrixComparisonHelper | Comparison_struct |
Private Attributes | |
Vector< int > | Index_of_diagonal_entries |
Vector whose i'th entry contains the index of the last entry below or on the diagonal of the i'th row of the matrix. More... | |
unsigned | Serial_matrix_matrix_multiply_method |
Flag to determine which matrix-matrix multiplication method is used (for serial (or global) matrices) More... | |
unsigned | Distributed_matrix_matrix_multiply_method |
Flag to determine which matrix-matrix multiplication method is used (for distributed matrices) More... | |
CRMatrix< double > | CR_matrix |
Storage for the Matrix in CR Format. More... | |
bool | Built |
Flag to indicate whether the matrix has been built - i.e. the distribution has been setup AND the matrix has been assembled. More... | |
Additional Inherited Members | |
![]() | |
void | range_check (const unsigned long &i, const unsigned long &j) const |
Range check to catch when an index is out of bounds, if so, it issues a warning message and dies by throwing an OomphLibError . More... | |
![]() | |
void | clear_distribution () |
clear the distribution of this distributable linear algebra object More... | |
![]() | |
LinearSolver * | Linear_solver_pt |
LinearSolver * | Default_linear_solver_pt |
A class for compressed row matrices. This is a distributable object.
Definition at line 872 of file matrices.h.
oomph::CRDoubleMatrix::CRDoubleMatrix | ( | ) |
Default constructor.
Definition at line 1238 of file matrices.cc.
References oomph::DoubleMatrixBase::Default_linear_solver_pt, and oomph::DoubleMatrixBase::Linear_solver_pt.
oomph::CRDoubleMatrix::CRDoubleMatrix | ( | const LinearAlgebraDistribution * | distribution_pt, |
const unsigned & | ncol, | ||
const Vector< double > & | value, | ||
const Vector< int > & | column_index, | ||
const Vector< int > & | row_start | ||
) |
Constructor: vector of values, vector of column indices, vector of row starts and number of rows and columns.
Constructor: Takes the distribution and the number of columns, as well as the vector of values, vector of column indices,vector of row starts.
Definition at line 1336 of file matrices.cc.
References oomph::DoubleMatrixBase::Default_linear_solver_pt, oomph::DoubleMatrixBase::Linear_solver_pt, oomph::CCDoubleMatrix::ncol(), and oomph::LinearAlgebraDistribution::nrow_local().
oomph::CRDoubleMatrix::CRDoubleMatrix | ( | const LinearAlgebraDistribution * | distribution_pt | ) |
Constructor: just stores the distribution but does not build the matrix.
Constructor: just stores the distribution but does not build the matrix
Definition at line 1310 of file matrices.cc.
References oomph::DoubleMatrixBase::Default_linear_solver_pt, and oomph::DoubleMatrixBase::Linear_solver_pt.
oomph::CRDoubleMatrix::CRDoubleMatrix | ( | const CRDoubleMatrix & | matrix | ) |
Copy constructor.
Definition at line 1258 of file matrices.cc.
References oomph::CCMatrix< double >::build_without_copy(), column_index(), oomph::DoubleMatrixBase::Default_linear_solver_pt, oomph::DistributableLinearAlgebraObject::distribution_pt(), oomph::DoubleMatrixBase::Linear_solver_pt, ncol(), oomph::SparseMatrix< double, CCMatrix< double > >::nnz(), nnz(), oomph::DistributableLinearAlgebraObject::nrow_local(), row_start(), and value().
|
virtual |
Destructor.
Definition at line 1367 of file matrices.cc.
References oomph::DoubleMatrixBase::Default_linear_solver_pt.
void oomph::CRDoubleMatrix::add | ( | const CRDoubleMatrix & | matrix_in, |
CRDoubleMatrix & | result_matrix | ||
) | const |
element-wise addition of this matrix with matrix_in.
Element-wise addition of this matrix with matrix_in.
Definition at line 3495 of file matrices.cc.
References build(), oomph::LinearAlgebraDistribution::built(), built(), column_index(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, ncol(), oomph::CCDoubleMatrix::ncol(), nrow(), oomph::CCDoubleMatrix::nrow(), row_start(), oomph::SparseMatrix< double, CCMatrix< double > >::value(), and value().
void oomph::CRDoubleMatrix::build | ( | const LinearAlgebraDistribution * | distribution_pt, |
const unsigned & | ncol, | ||
const Vector< double > & | value, | ||
const Vector< int > & | column_index, | ||
const Vector< int > & | row_start | ||
) |
build method: vector of values, vector of column indices, vector of row starts and number of rows and columns.
build method: Takes the distribution and the number of columns, as well as the vector of values, vector of column indices,vector of row starts.
Definition at line 1692 of file matrices.cc.
References oomph::CCMatrix< double >::build(), and oomph::DoubleMatrixBase::Default_linear_solver_pt.
Referenced by add(), oomph::ComplexDampedJacobi< MATRIX >::complex_solve_helper(), oomph::CRDoubleMatrixHelpers::concatenate(), oomph::CRDoubleMatrixHelpers::concatenate_without_communication(), oomph::CRDoubleMatrixHelpers::create_uniformly_distributed_matrix(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::Problem::get_eigenproblem_matrices(), oomph::Problem::get_jacobian(), get_matrix_transpose(), oomph::BlockPreconditioner< CRDoubleMatrix >::internal_get_block(), oomph::TrilinosEpetraHelpers::multiply(), multiply(), oomph::LagrangeEnforcedFlowPreconditioner::setup(), oomph::HelmholtzMGPreconditioner< DIM >::setup_mg_structures(), oomph::GS< CRDoubleMatrix >::solve(), and oomph::GS< CRDoubleMatrix >::solve_helper().
void oomph::CRDoubleMatrix::build | ( | const LinearAlgebraDistribution * | distribution_pt | ) |
rebuild the matrix - assembles an empty matrix will a defined distribution
Rebuild the matrix - assembles an empty matrix with a defined distribution.
Definition at line 1377 of file matrices.cc.
void oomph::CRDoubleMatrix::build | ( | const unsigned & | ncol, |
const Vector< double > & | value, | ||
const Vector< int > & | column_index, | ||
const Vector< int > & | row_start | ||
) |
keeps the existing distribution and just matrix that is stored
method to rebuild the matrix, but not the distribution
Definition at line 1714 of file matrices.cc.
void oomph::CRDoubleMatrix::build_without_copy | ( | const unsigned & | ncol, |
const unsigned & | nnz, | ||
double * | value, | ||
int * | column_index, | ||
int * | row_start | ||
) |
keeps the existing distribution and just matrix that is stored without copying the matrix data
method to rebuild the matrix, but not the distribution
Definition at line 1731 of file matrices.cc.
Referenced by oomph::NavierStokesSchurComplementPreconditioner::assemble_inv_press_and_veloc_mass_matrix_diagonal(), oomph::PressureBasedSolidLSCPreconditioner::assemble_mass_matrix_diagonal(), oomph::CRDoubleMatrixHelpers::concatenate_without_communication(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::Problem::get_eigenproblem_matrices(), oomph::Problem::get_jacobian(), global_matrix(), oomph::BlockPreconditioner< CRDoubleMatrix >::internal_get_block(), oomph::TrilinosEpetraHelpers::multiply(), and multiply().
|
inline |
access function to the Built flag - indicates whether the matrix has been build - i.e. the distribution has been defined and the matrix assembled.
Definition at line 1177 of file matrices.h.
References oomph::CRDoubleMatrixHelpers::inf_norm().
Referenced by add(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::CRDoubleMatrixHelpers::create_uniformly_distributed_matrix(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::BlockPreconditioner< CRDoubleMatrix >::get_concatenated_block(), oomph::TrilinosEpetraHelpers::multiply(), multiply(), oomph::DoubleVector::norm(), and oomph::BlockPreconditioner< CRDoubleMatrix >::set_replacement_dof_block().
void oomph::CRDoubleMatrix::clear | ( | ) |
clear
Clean method.
Definition at line 1677 of file matrices.cc.
References oomph::LinearSolver::clean_up_memory(), and oomph::DoubleMatrixBase::Linear_solver_pt.
Referenced by oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::HypreInterface::hypre_matrix_setup(), and oomph::FSIPreconditioner::setup().
|
inline |
Access to C-style column index array.
Definition at line 1056 of file matrices.h.
Referenced by add(), oomph::BlockPreconditioner< CRDoubleMatrix >::block_setup(), oomph::CRDoubleMatrixHelpers::concatenate(), CRDoubleMatrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::SuperLUSolver::factorise_serial(), oomph::BlockPreconditioner< CRDoubleMatrix >::internal_get_block(), multiply(), oomph::ILUZeroPreconditioner< CRDoubleMatrix >::setup(), oomph::PseudoElasticPreconditionerSubsidiaryPreconditionerOld::setup(), oomph::HelmholtzMGPreconditioner< DIM >::setup_coarsest_level_structures(), and oomph::GS< CRDoubleMatrix >::solve_helper().
|
inline |
Access to C-style column index array (const version)
Definition at line 1059 of file matrices.h.
Vector< double > oomph::CRDoubleMatrix::diagonal_entries | ( | ) | const |
returns a Vector of diagonal entries of this matrix. This only works with square matrices. This condition may be relaxed in the future if need be.
Return the diagonal entries of the matrix. This only works with square matrices. This condition may be relaxed in the future if need be.
Definition at line 3442 of file matrices.cc.
References i, oomph::CCDoubleMatrix::ncol(), and oomph::CCDoubleMatrix::nrow().
Referenced by oomph::LagrangeEnforcedFlowPreconditioner::setup().
|
inline |
Access function to Distributed_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for distributed matrices. Method 1: Trilinos Epetra Matrix Matrix multiply. Method 2: Trilinos Epetra Matrix Matrix multiply (ML based).
Definition at line 1158 of file matrices.h.
Referenced by oomph::CRDoubleMatrixHelpers::deep_copy().
|
inline |
Read only access function (const version) to Distributed_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for distributed matrices. Method 1: Trilinos Epetra Matrix Matrix multiply. Method 2: Trilinos Epetra Matrix Matrix multiply (ML based).
Definition at line 1169 of file matrices.h.
bool oomph::CRDoubleMatrix::entries_are_sorted | ( | const bool & | doc_unordered_entries = false | ) | const |
Runs through the column index vector and checks if the entries follow the regular lexicographical ordering of matrix entries, i.e. it will check (at the i-th row of the matrix) if the entries in the column index vector associated with this row are in increasing order.
Runs through the column index vector and checks if the entries are arranged arbitrarily or if they follow the regular lexicographical of matrices. If a boolean argument is provided with the assignment TRUE then information on the first entry which is not in the correct position will also be given.
Definition at line 1390 of file matrices.cc.
References i, oomph::CCDoubleMatrix::nrow(), and oomph::oomph_info.
|
inline |
Access function: returns the vector Index_of_diagonal_entries. The i-th entry of the vector contains the index of the last entry below or on the diagonal. If there are no entries below or on the diagonal then the corresponding entry is -1. If, however, there are no entries in the row then the entry is irrelevant and is kept as the initialised value; 0.
Definition at line 912 of file matrices.h.
References oomph::Global_string_for_annotation::string().
void oomph::CRDoubleMatrix::get_matrix_transpose | ( | CRDoubleMatrix * | result | ) | const |
Returns the transpose of this matrix.
Compute transpose of matrix.
Definition at line 3250 of file matrices.cc.
References oomph::LinearAlgebraDistribution::build(), build(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::CCDoubleMatrix::ncol(), oomph::SparseMatrix< double, CCMatrix< double > >::nnz(), oomph::CCDoubleMatrix::nrow(), and oomph::SparseMatrix< double, CCMatrix< double > >::value().
CRDoubleMatrix * oomph::CRDoubleMatrix::global_matrix | ( | ) | const |
if this matrix is distributed then a the equivalent global matrix is built using new and returned. The calling method is responsible for the destruction of the new matrix.
if this matrix is distributed then the equivalent global matrix is built using new and returned. The calling method is responsible for the destruction of the new matrix.
Definition at line 2458 of file matrices.cc.
References build_without_copy(), i, oomph::CCDoubleMatrix::ncol(), oomph::SparseMatrix< double, CCMatrix< double > >::nnz(), oomph::CCDoubleMatrix::nrow(), and oomph::SparseMatrix< double, CCMatrix< double > >::value().
Referenced by oomph::ILUZeroPreconditioner< CRDoubleMatrix >::setup().
double oomph::CRDoubleMatrix::inf_norm | ( | ) | const |
returns the inf-norm of this matrix
Compute infinity (maximum) norm of matrix.
Definition at line 3392 of file matrices.cc.
References i, and oomph::SparseMatrix< double, CCMatrix< double > >::value().
Referenced by oomph::PseudoElasticPreconditionerScalingHelperOld::s_inf_norm().
|
virtual |
LU back solve for given RHS.
Do back-substitution.
Definition at line 1772 of file matrices.cc.
References oomph::DoubleVector::built(), oomph::DoubleMatrixBase::Default_linear_solver_pt, and oomph::DistributableLinearAlgebraObject::distribution_pt().
|
virtual |
LU decomposition using SuperLU if matrix is not distributed or distributed onto a single processor.
Do LU decomposition.
Definition at line 1749 of file matrices.cc.
References oomph::DoubleMatrixBase::Default_linear_solver_pt.
void oomph::CRDoubleMatrix::matrix_reduction | ( | const double & | alpha, |
CRDoubleMatrix & | reduced_matrix | ||
) |
For every row, find the maximum absolute value of the entries in this row. Set all values that are less than alpha times this maximum to zero and return the resulting matrix in reduced_matrix. Note: Diagonal entries are retained regardless of their size.
For every row, find the maximum absolute value of the entries in this row. Set all values that are less than alpha times this maximum to zero and return the resulting matrix in reduced_matrix. Note: Diagonal entries are retained regardless of their size.
Definition at line 2392 of file matrices.cc.
References oomph::CCMatrix< double >::build(), i, oomph::CCDoubleMatrix::ncol(), and oomph::SparseMatrix< double, CCMatrix< double > >::value().
|
virtual |
Multiply the matrix by the vector x: soln=Ax.
Multiply the matrix by the vector x.
Implements oomph::DoubleMatrixBase.
Definition at line 1805 of file matrices.cc.
References oomph::DoubleVector::build(), oomph::DoubleVector::built(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DoubleVector::initialise(), oomph::TrilinosEpetraHelpers::multiply(), oomph::CCDoubleMatrix::ncol(), oomph::LinearAlgebraDistribution::nrow(), oomph::CCDoubleMatrix::nrow(), oomph::SparseMatrix< double, CCMatrix< double > >::value(), and oomph::DoubleVector::values_pt().
Referenced by oomph::ProblemBasedShiftInvertOperator::apply(), oomph::MatrixVectorProduct::multiply(), oomph::DoubleVector::norm(), oomph::LagrangeEnforcedFlowPreconditioner::setup(), oomph::PressureBasedSolidLSCPreconditioner::setup(), and oomph::NavierStokesSchurComplementPreconditioner::setup().
void oomph::CRDoubleMatrix::multiply | ( | const CRDoubleMatrix & | matrix_in, |
CRDoubleMatrix & | result | ||
) | const |
Function to multiply this matrix by the CRDoubleMatrix matrix_in. In a serial matrix, there are 4 methods available: Method 1: First runs through this matrix and matrix_in to find the storage requirements for result - arrays of the correct size are then allocated before performing the calculation. Minimises memory requirements but more costly. Method 2: Grows storage for values and column indices of result 'on the fly' using an array of maps. Faster but more memory intensive. Method 3: Grows storage for values and column indices of result 'on the fly' using a vector of vectors. Not particularly impressive on the platforms we tried... Method 4: Trilinos Epetra Matrix Matrix multiply. Method 5: Trilinox Epetra Matrix Matrix Mulitply (ml based) If Trilinos is installed then Method 4 is employed by default, otherwise Method 2 is employed by default. In a distributed matrix, only Trilinos Epetra Matrix Matrix multiply is available.
Function to multiply this matrix by the CRDoubleMatrix matrix_in. In a serial matrix, there are 4 methods available: Method 1: First runs through this matrix and matrix_in to find the storage requirements for result - arrays of the correct size are then allocated before performing the calculation. Minimises memory requirements but more costly. Method 2: Grows storage for values and column indices of result 'on the fly' using an array of maps. Faster but more memory intensive. Method 3: Grows storage for values and column indices of result 'on the fly' using a vector of vectors. Not particularly impressive on the platforms we tried... Method 4: Trilinos Epetra Matrix Matrix multiply. Method 5: Trilinox Epetra Matrix Matrix Mulitply (ml based) If Trilinos is installed then Method 4 is employed by default, otherwise Method 2 is employed by default. In a distributed matrix, only Trilinos Epetra Matrix Matrix multiply is available.
Definition at line 2017 of file matrices.cc.
References build(), build_without_copy(), built(), column_index(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DistributableLinearAlgebraObject::distribution_built(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::SparseMatrix< double, CCMatrix< double > >::M, oomph::TrilinosEpetraHelpers::multiply(), oomph::SparseMatrix< double, CCMatrix< double > >::N, ncol(), oomph::SparseMatrix< double, CCMatrix< double > >::Nnz, oomph::CCDoubleMatrix::nrow(), row_start(), oomph::SparseMatrix< double, CCMatrix< double > >::Value, oomph::SparseMatrix< double, CCMatrix< double > >::value(), and value().
|
virtual |
Multiply the transposed matrix by the vector x: soln=A^T x.
Implements oomph::DoubleMatrixBase.
Definition at line 1907 of file matrices.cc.
References oomph::DoubleVector::build(), oomph::DoubleVector::built(), oomph::LinearAlgebraDistribution::communicator_pt(), oomph::DistributableLinearAlgebraObject::distribution_pt(), i, oomph::DoubleVector::initialise(), oomph::TrilinosEpetraHelpers::multiply(), oomph::CCDoubleMatrix::ncol(), oomph::LinearAlgebraDistribution::nrow(), oomph::CCDoubleMatrix::nrow(), oomph::SparseMatrix< double, CCMatrix< double > >::value(), and oomph::DoubleVector::values_pt().
Referenced by oomph::MatrixVectorProduct::multiply_transpose(), and oomph::SuperLUSolver::solve().
|
inlinevirtual |
Return the number of columns of the matrix.
Implements oomph::DoubleMatrixBase.
Definition at line 998 of file matrices.h.
Referenced by add(), CRDoubleMatrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::BlockPreconditioner< CRDoubleMatrix >::get_concatenated_block(), oomph::TrilinosEpetraHelpers::multiply(), multiply(), oomph::BlockPreconditioner< CRDoubleMatrix >::set_replacement_dof_block(), oomph::MatrixVectorProduct::setup(), and oomph::HelmholtzMGPreconditioner< DIM >::setup_coarsest_level_structures().
|
inline |
Return the number of nonzero entries (the local nnz)
Definition at line 1068 of file matrices.h.
References oomph::TrilinosEpetraHelpers::multiply().
Referenced by oomph::BlockPreconditioner< CRDoubleMatrix >::block_setup(), CRDoubleMatrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::SuperLUSolver::factorise_serial(), oomph::LagrangeEnforcedFlowPreconditioner::setup(), oomph::HyprePreconditioner::setup(), and oomph::HelmholtzMGPreconditioner< DIM >::setup_coarsest_level_structures().
|
inlinevirtual |
Return the number of rows of the matrix.
Implements oomph::DoubleMatrixBase.
Definition at line 992 of file matrices.h.
References oomph::DistributableLinearAlgebraObject::nrow().
Referenced by add(), oomph::CRDoubleMatrixHelpers::concatenate(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::SuperLUSolver::factorise_distributed(), oomph::Problem::get_eigenproblem_matrices(), oomph::HypreInterface::hypre_matrix_setup(), oomph::MGSolver< DIM >::modify_restriction_matrices(), oomph::TrilinosEpetraHelpers::multiply(), oomph::ILUZeroPreconditioner< CRDoubleMatrix >::setup(), oomph::HyprePreconditioner::setup(), and oomph::HelmholtzMGPreconditioner< DIM >::setup_coarsest_level_structures().
|
inlinevirtual |
Overload the round-bracket access operator for read-only access. In a distributed matrix i refers to the local row index.
Implements oomph::DoubleMatrixBase.
Definition at line 1045 of file matrices.h.
|
inline |
Broken assignment operator.
Definition at line 898 of file matrices.h.
References oomph::BrokenCopy::broken_assign().
|
inlinevirtual |
Output the "bottom right" entry regardless of it being zero or not (this allows automatic detection of matrix size in e.g. matlab, python).
Implements oomph::Matrix< double, CRDoubleMatrix >.
Definition at line 1006 of file matrices.h.
void oomph::CRDoubleMatrix::redistribute | ( | const LinearAlgebraDistribution *const & | dist_pt | ) |
The contents of the matrix are redistributed to match the new distribution. In a non-MPI build this method does nothing. NOTE 1: The current distribution and the new distribution must have the same number of global rows. NOTE 2: The current distribution and the new distribution must have the same Communicator.
Definition at line 2580 of file matrices.cc.
References oomph::CCMatrix< double >::build(), oomph::CCMatrix< double >::build_without_copy(), oomph::LinearAlgebraDistribution::communicator_pt(), oomph::LinearAlgebraDistribution::distributed(), oomph::LinearAlgebraDistribution::first_row(), i, oomph::CCDoubleMatrix::ncol(), oomph::SparseMatrix< double, CCMatrix< double > >::nnz(), oomph::LinearAlgebraDistribution::nrow(), oomph::CCDoubleMatrix::nrow(), oomph::LinearAlgebraDistribution::nrow_local(), and oomph::SparseMatrix< double, CCMatrix< double > >::value().
Referenced by oomph::CRDoubleMatrixHelpers::create_uniformly_distributed_matrix(), oomph::Problem::get_eigenproblem_matrices(), and oomph::Problem::get_jacobian().
|
inline |
Access to C-style row_start array.
Definition at line 1050 of file matrices.h.
Referenced by add(), oomph::CRDoubleMatrixHelpers::concatenate(), CRDoubleMatrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::SuperLUSolver::factorise_serial(), oomph::CRDoubleMatrixHelpers::gershgorin_eigenvalue_estimate(), oomph::CRDoubleMatrixHelpers::inf_norm(), oomph::BlockPreconditioner< CRDoubleMatrix >::internal_get_block(), oomph::MGSolver< DIM >::modify_restriction_matrices(), multiply(), oomph::MatrixBasedLumpedPreconditioner< oomph::CRDoubleMatrix >::setup(), oomph::ILUZeroPreconditioner< CRDoubleMatrix >::setup(), oomph::PseudoElasticPreconditionerSubsidiaryPreconditionerOld::setup(), oomph::HelmholtzMGPreconditioner< DIM >::setup_coarsest_level_structures(), and oomph::GS< CRDoubleMatrix >::solve_helper().
|
inline |
Access to C-style row_start array (const version)
Definition at line 1053 of file matrices.h.
|
inline |
Access function to Serial_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for serial matrices. Method 1: First runs through this matrix and matrix_in to find the storage requirements for result - arrays of the correct size are then allocated before performing the calculation. Minimises memory requirements but more costly. Method 2: Grows storage for values and column indices of result 'on the fly' using an array of maps. Faster but more memory intensive. Method 3: Grows storage for values and column indices of result 'on the fly' using a vector of vectors. Not particularly impressive on the platforms we tried... Method 4: Trilinos Epetra Matrix Matrix multiply. Method 5: Trilinos Epetra Matrix Matrix multiply (ML based).
Definition at line 1127 of file matrices.h.
Referenced by oomph::CRDoubleMatrixHelpers::deep_copy().
|
inline |
Read only access function (const version) to Serial_matrix_matrix_multiply_method, the flag which determines the matrix matrix multiplication method used for serial matrices. Method 1: First runs through this matrix and matrix_in to find the storage requirements for result - arrays of the correct size are then allocated before performing the calculation. Minimises memory requirements but more costly. Method 2: Grows storage for values and column indices of result 'on the fly' using an array of maps. Faster but more memory intensive. Method 3: Grows storage for values and column indices of result 'on the fly' using a vector of vectors. Not particularly impressive on the platforms we tried... Method 4: Trilinos Epetra Matrix Matrix multiply. Method 5: Trilinos Epetra Matrix Matrix multiply (ML based).
Definition at line 1148 of file matrices.h.
void oomph::CRDoubleMatrix::sort_entries | ( | ) |
Sorts the entries associated with each row of the matrix in the column index vector and the value vector into ascending order and sets up the Index_of_diagonal_entries vector.
This helper function sorts the entries in the column index vector and the value vector. During the construction of the matrix the entries were most likely assigned in an arbitrary order. As a result, it cannot be assumed that the entries in the column index vector corresponding to each row of the matrix have been arranged in increasing order. During the setup an additional vector will be set up; Index_of_diagonal_entries. The i-th entry of this vector contains the index of the last entry below or on the diagonal. If there are no entries below or on the diagonal then the corresponding entry is -1. If, however, there are no entries in the row then the entry is irrelevant and is kept as the initialised value; 0.
Definition at line 1473 of file matrices.cc.
References i, oomph::CCDoubleMatrix::nrow(), and oomph::SparseMatrix< double, CCMatrix< double > >::value().
Referenced by oomph::GS< CRDoubleMatrix >::setup_helper().
|
inlinevirtual |
Indexed output function to print a matrix to the stream outfile as i,j,a(i,j) for a(i,j)!=0 only.
Implements oomph::Matrix< double, CRDoubleMatrix >.
Definition at line 1013 of file matrices.h.
|
inline |
Indexed output function to print a matrix to a file as i,j,a(i,j) for a(i,j)!=0 only. Specify filename. This uses acual global row numbers.
Definition at line 1021 of file matrices.h.
References oomph::LinearAlgebraDistribution::first_row(), and i.
Referenced by oomph::NavierStokesSchurComplementPreconditioner::setup().
|
inline |
Access to C-style value array.
Definition at line 1062 of file matrices.h.
Referenced by add(), oomph::CRDoubleMatrixHelpers::concatenate(), CRDoubleMatrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::CRDoubleMatrixHelpers::deep_copy(), oomph::MumpsSolver::factorise(), oomph::SuperLUSolver::factorise_distributed(), oomph::SuperLUSolver::factorise_serial(), oomph::CRDoubleMatrixHelpers::gershgorin_eigenvalue_estimate(), oomph::CRDoubleMatrixHelpers::inf_norm(), oomph::BlockPreconditioner< CRDoubleMatrix >::internal_get_block(), oomph::MGSolver< DIM >::modify_restriction_matrices(), multiply(), oomph::MatrixBasedLumpedPreconditioner< oomph::CRDoubleMatrix >::setup(), oomph::ILUZeroPreconditioner< CRDoubleMatrix >::setup(), oomph::PseudoElasticPreconditionerSubsidiaryPreconditionerOld::setup(), oomph::HelmholtzMGPreconditioner< DIM >::setup_coarsest_level_structures(), and oomph::GS< CRDoubleMatrix >::solve_helper().
|
inline |
Access to C-style value array (const version)
Definition at line 1065 of file matrices.h.
|
private |
Flag to indicate whether the matrix has been built - i.e. the distribution has been setup AND the matrix has been assembled.
Definition at line 1217 of file matrices.h.
struct oomph::CRDoubleMatrix::CRDoubleMatrixComparisonHelper oomph::CRDoubleMatrix::Comparison_struct |
|
private |
Storage for the Matrix in CR Format.
Definition at line 1213 of file matrices.h.
|
private |
Flag to determine which matrix-matrix multiplication method is used (for distributed matrices)
Definition at line 1210 of file matrices.h.
|
private |
Vector whose i'th entry contains the index of the last entry below or on the diagonal of the i'th row of the matrix.
Definition at line 1202 of file matrices.h.
|
private |
Flag to determine which matrix-matrix multiplication method is used (for serial (or global) matrices)
Definition at line 1206 of file matrices.h.