strumpack::HSS::HSSMatrix< scalar_t > Class Template Reference

Class to represent a sequential/threaded Hierarchically Semi-Separable matrix. More...

#include <HSSMatrix.hpp>

Inheritance diagram for strumpack::HSS::HSSMatrix< scalar_t >:
Collaboration diagram for strumpack::HSS::HSSMatrix< scalar_t >:

Public Member Functions

 HSSMatrix ()
 
 HSSMatrix (const DenseM_t &A, const opts_t &opts)
 
 HSSMatrix (std::size_t m, std::size_t n, const opts_t &opts)
 
 HSSMatrix (const structured::ClusterTree &t, const opts_t &opts)
 
 HSSMatrix (kernel::Kernel< real_t > &K, const opts_t &opts)
 
 HSSMatrix (const HSSMatrix< scalar_t > &other)
 
HSSMatrix< scalar_t > & operator= (const HSSMatrix< scalar_t > &other)
 
 HSSMatrix (HSSMatrix< scalar_t > &&other)=default
 
HSSMatrix< scalar_t > & operator= (HSSMatrix< scalar_t > &&other)=default
 
std::unique_ptr< HSSMatrixBase< scalar_t > > clone () const override
 
const HSSMatrix< scalar_t > * child (int c) const
 
HSSMatrix< scalar_t > * child (int c)
 
void compress (const DenseM_t &A, const opts_t &opts)
 
void compress (const std::function< void(DenseM_t &Rr, DenseM_t &Rc, DenseM_t &Sr, DenseM_t &Sc)> &Amult, const std::function< void(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J, DenseM_t &B)> &Aelem, const opts_t &opts)
 
void compress_with_coordinates (const DenseMatrix< real_t > &coords, const std::function< void(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J, DenseM_t &B)> &Aelem, const opts_t &opts)
 
void reset () override
 
void factor () override
 
void partial_factor ()
 
void solve (DenseM_t &b) const override
 
void forward_solve (WorkSolve< scalar_t > &w, const DenseM_t &b, bool partial) const override
 
void backward_solve (WorkSolve< scalar_t > &w, DenseM_t &x) const override
 
DenseM_t apply (const DenseM_t &b) const
 
void mult (Trans op, const DenseM_t &x, DenseM_t &y) const override
 
DenseM_t applyC (const DenseM_t &b) const
 
scalar_t get (std::size_t i, std::size_t j) const
 
DenseM_t extract (const std::vector< std::size_t > &I, const std::vector< std::size_t > &J) const
 
void extract_add (const std::vector< std::size_t > &I, const std::vector< std::size_t > &J, DenseM_t &B) const
 
std::size_t rank () const override
 
std::size_t memory () const override
 
std::size_t nonzeros () const override
 
std::size_t levels () const override
 
void print_info (std::ostream &out=std::cout, std::size_t roff=0, std::size_t coff=0) const override
 
DenseM_t dense () const
 
void shift (scalar_t sigma) override
 
void draw (std::ostream &of, std::size_t rlo=0, std::size_t clo=0) const override
 
void write (const std::string &fname) const
 
const HSSFactors< scalar_t > & ULV ()
 
- Public Member Functions inherited from strumpack::HSS::HSSMatrixBase< scalar_t >
 HSSMatrixBase (std::size_t m, std::size_t n, bool active)
 
virtual ~HSSMatrixBase ()=default
 
 HSSMatrixBase (const HSSMatrixBase< scalar_t > &other)
 
HSSMatrixBase< scalar_t > & operator= (const HSSMatrixBase< scalar_t > &other)
 
 HSSMatrixBase (HSSMatrixBase &&h)=default
 
HSSMatrixBaseoperator= (HSSMatrixBase &&h)=default
 
virtual std::unique_ptr< HSSMatrixBase< scalar_t > > clone () const =0
 
std::pair< std::size_t, std::size_t > dims () const
 
std::size_t rows () const override
 
std::size_t cols () const override
 
bool leaf () const
 
virtual std::size_t factor_nonzeros () const
 
const HSSMatrixBase< scalar_t > & child (int c) const
 
HSSMatrixBase< scalar_t > & child (int c)
 
bool is_compressed () const
 
bool is_untouched () const
 
bool active () const
 
virtual std::size_t levels () const =0
 
virtual void print_info (std::ostream &out=std::cout, std::size_t roff=0, std::size_t coff=0) const =0
 
void set_openmp_task_depth (int depth)
 
virtual void shift (scalar_t sigma) override=0
 
virtual void draw (std::ostream &of, std::size_t rlo, std::size_t clo) const
 
virtual void forward_solve (WorkSolveMPI< scalar_t > &w, const DistM_t &b, bool partial) const
 
virtual void backward_solve (WorkSolveMPI< scalar_t > &w, DistM_t &x) const
 
virtual const BLACSGridgrid () const
 
virtual const BLACSGridgrid (const BLACSGrid *local_grid) const
 
virtual const BLACSGridgrid_local () const
 
virtual int Ptotal () const
 
virtual int Pactive () const
 
virtual void to_block_row (const DistM_t &A, DenseM_t &sub_A, DistM_t &leaf_A) const
 
virtual void allocate_block_row (int d, DenseM_t &sub_A, DistM_t &leaf_A) const
 
virtual void from_block_row (DistM_t &A, const DenseM_t &sub_A, const DistM_t &leaf_A, const BLACSGrid *lg) const
 
- Public Member Functions inherited from strumpack::structured::StructuredMatrix< scalar_t >
virtual ~StructuredMatrix ()=default
 
virtual std::size_t rows () const =0
 
virtual std::size_t cols () const =0
 
virtual std::size_t memory () const =0
 
virtual std::size_t nonzeros () const =0
 
virtual std::size_t rank () const =0
 
virtual std::size_t local_rows () const
 
virtual std::size_t begin_row () const
 
virtual std::size_t end_row () const
 
virtual const std::vector< int > & dist () const
 
virtual const std::vector< int > & rdist () const
 
virtual const std::vector< int > & cdist () const
 
virtual void mult (Trans op, const DenseMatrix< scalar_t > &x, DenseMatrix< scalar_t > &y) const
 
void mult (Trans op, int m, const scalar_t *x, int ldx, scalar_t *y, int ldy) const
 
virtual void mult (Trans op, const DistributedMatrix< scalar_t > &x, DistributedMatrix< scalar_t > &y) const
 
virtual void factor ()
 
virtual void solve (DenseMatrix< scalar_t > &b) const
 
virtual void solve (int nrhs, scalar_t *b, int ldb) const
 
virtual void solve (DistributedMatrix< scalar_t > &b) const
 
virtual void shift (scalar_t s)
 

Static Public Member Functions

static HSSMatrix< scalar_t > read (const std::string &fname)
 

Friends

class HSSMatrixMPI< scalar_t >
 
template<typename T >
void apply_HSS (Trans ta, const HSSMatrix< T > &a, const DenseMatrix< T > &b, T beta, DenseMatrix< T > &c)
 
template<typename T >
void draw (const HSSMatrix< T > &H, const std::string &name)
 

Detailed Description

template<typename scalar_t>
class strumpack::HSS::HSSMatrix< scalar_t >

Class to represent a sequential/threaded Hierarchically Semi-Separable matrix.

This is for non-symmetric matrices, but can be used with symmetric matrices as well. This class inherits from HSSMatrixBase. An HSS matrix is represented recursively, using 2 children which are also HSSMatrix objects, except at the lowest level (the leafs), where the HSSMatrix has no children. Hence, the tree representing the HSS matrix is always a binary tree, but not necessarily a complete binary tree.

Template Parameters
scalar_tCan be float, double, std:complex<float> or std::complex<double>.
See also
HSSMatrixMPI, HSSMatrixBase

Constructor & Destructor Documentation

◆ HSSMatrix() [1/7]

template<typename scalar_t >
strumpack::HSS::HSSMatrix< scalar_t >::HSSMatrix ( )

Default constructor, constructs an empty 0 x 0 matrix.

◆ HSSMatrix() [2/7]

template<typename scalar_t >
strumpack::HSS::HSSMatrix< scalar_t >::HSSMatrix ( const DenseM_t A,
const opts_t opts 
)

Construct an HSS representation for a dense matrix A. The HSS tree will be constructed by splitting the row/column set evenly and recursively until the leafs in the HSS tree are smaller than opts.leaf_size(). Alternative constructors can be used to specify a specific HSS partitioning tree. Internally, this will call the appropriate compress routine to construct the HSS representation, using the options (such as compression tolerance, adaptive compression scheme etc) as specified in the HSSOptions opts object.

Parameters
Adense matrix (unmodified) to compress as HSS
optsobject containing a number of options for HSS compression
See also
DenseMatrix
HSSOptions

◆ HSSMatrix() [3/7]

template<typename scalar_t >
strumpack::HSS::HSSMatrix< scalar_t >::HSSMatrix ( std::size_t  m,
std::size_t  n,
const opts_t opts 
)

Construct an HSS representation for an m x n matrix. The HSS tree will be constructed by splitting the row/column set evenly and recursively until the leafs in the HSS tree are smaller than opts.leaf_size(). After construction, the HSS matrix will be empty, and can be filled by calling one of the compress member routines. Alternative constructors can be used to specify a specific HSS partitioning tree.

Parameters
mnumber of rows in the constructed HSS matrix
nnumber of rows in the constructed HSS matrix
optsobject containing a number of options for HSS compression
See also
compress, HSSOptions

◆ HSSMatrix() [4/7]

template<typename scalar_t >
strumpack::HSS::HSSMatrix< scalar_t >::HSSMatrix ( const structured::ClusterTree t,
const opts_t opts 
)

Construct an HSS matrix using a specified HSS tree. After construction, the HSS matrix will be empty, and can be filled by calling one of the compress member routines.

Parameters
ttree specifying the HSS matrix partitioning
optsobject containing a number of options for HSS compression
See also
compress, HSSOptions

◆ HSSMatrix() [5/7]

template<typename scalar_t >
strumpack::HSS::HSSMatrix< scalar_t >::HSSMatrix ( kernel::Kernel< real_t > &  K,
const opts_t opts 
)

Construct an HSS approximation for the kernel matrix K.

Parameters
KKernel matrix object. The data associated with this kernel will be permuted according to the clustering algorithm selected by the HSSOptions objects. The permutation will be stored in the kernel object.
optsobject containing a number of HSS options

◆ HSSMatrix() [6/7]

template<typename scalar_t >
strumpack::HSS::HSSMatrix< scalar_t >::HSSMatrix ( const HSSMatrix< scalar_t > &  other)

Copy constructor. Copying an HSSMatrix can be an expensive operation.

Parameters
otherHSS matrix to be copied

◆ HSSMatrix() [7/7]

template<typename scalar_t >
strumpack::HSS::HSSMatrix< scalar_t >::HSSMatrix ( HSSMatrix< scalar_t > &&  other)
default

Move constructor.

Parameters
otherHSS matrix to be moved from, will be emptied

Member Function Documentation

◆ apply()

template<typename scalar_t >
DenseM_t strumpack::HSS::HSSMatrix< scalar_t >::apply ( const DenseM_t b) const

Multiply this HSS matrix with a dense matrix (vector), ie, compute x = this * b.

Parameters
bMatrix to multiply with, from the left.
Returns
The result of this * b.
See also
applyC, HSS::apply_HSS

◆ applyC()

template<typename scalar_t >
DenseM_t strumpack::HSS::HSSMatrix< scalar_t >::applyC ( const DenseM_t b) const

Multiply the transpose or complex conjugate of this HSS matrix with a dense matrix (vector), ie, compute x = this^C * b.

Parameters
bMatrix to multiply with, from the left.
Returns
The result of this^C * b.
See also
apply, HSS::apply_HSS

◆ backward_solve()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::backward_solve ( WorkSolve< scalar_t > &  w,
DenseM_t x 
) const
overridevirtual

Perform only the backward phase of the ULV linear solve. This is for advanced use only, typically to be used in combination with partial_factor. You should really just use factor/solve when possible.

Parameters
wtemporary working storage, to pass information from forward_solve to backward_solve
bon input, the vector obtained from forward_solve, on output the solution of A x = b (with A this HSS matrix). The vector b should be b.rows() == cols().
See also
factor, partial_factor, backward_solve

Reimplemented from strumpack::HSS::HSSMatrixBase< scalar_t >.

◆ child() [1/2]

template<typename scalar_t >
HSSMatrix< scalar_t > * strumpack::HSS::HSSMatrix< scalar_t >::child ( int  c)
inline

Return a raw (non-owning) pointer to child c of this HSS matrix. A child of an HSS matrix is itself an HSS matrix. The value of c should be 0 or 1, and this HSS matrix should not be a leaf!

◆ child() [2/2]

template<typename scalar_t >
const HSSMatrix< scalar_t > * strumpack::HSS::HSSMatrix< scalar_t >::child ( int  c) const
inline

Return a const raw (non-owning) pointer to child c of this HSS matrix. A child of an HSS matrix is itself an HSS matrix. The value of c should be 0 or 1, and this HSS matrix should not be a leaf!

◆ clone()

template<typename scalar_t >
std::unique_ptr< HSSMatrixBase< scalar_t > > strumpack::HSS::HSSMatrix< scalar_t >::clone ( ) const
overridevirtual

Create a clone of this matrix. TODO remove this!!???

Implements strumpack::HSS::HSSMatrixBase< scalar_t >.

◆ compress() [1/2]

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::compress ( const DenseM_t A,
const opts_t opts 
)

Initialize this HSS matrix as the compressed HSS representation of a given dense matrix. The HSS matrix should have been constructed with the proper sizes, i.e., rows() == A.rows() and cols() == A.cols(). Internaly, this will call the appropriate compress routine to construct the HSS representation, using the options (such as compression tolerance, adaptive compression scheme etc) as specified in the HSSOptions opts object.

Parameters
Adense matrix (unmodified) to compress as HSS
optsobject containing a number of options for HSS compression
See also
DenseMatrix
HSSOptions

◆ compress() [2/2]

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::compress ( const std::function< void(DenseM_t &Rr, DenseM_t &Rc, DenseM_t &Sr, DenseM_t &Sc)> &  Amult,
const std::function< void(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J, DenseM_t &B)> &  Aelem,
const opts_t opts 
)

Initialize this HSS matrix as the compressed HSS representation. The compression uses the matrix-(multiple)vector multiplication routine Amult, and the element (sub-matrix) extraction routine Aelem, both provided by the user. The HSS matrix should have been constructed with the proper sizes, i.e., rows() == A.rows() and cols() == A.cols(). Internaly, this will call the appropriate compress routine to construct the HSS representation, using the options (such as compression tolerance, adaptive compression scheme etc) as specified in the HSSOptions opts object.

Parameters
Amultmatrix-(multiple)vector product routine. This can be a functor, or a lambda function for instance.
RrParameter to the matvec routine. Random matrix. This will be set by the compression routine.
RcParameter to the matvec routine. Random matrix. This will be set by the compression routine.
SrParameter to the matvec routine. Random sample matrix, to be computed by the matrix-(multiple)vector multiplication routine as A*Rr. This will aready be allocated by the compression routine and should be Sr.rows() == this->rows() and Sr.cols() == Rr.cols().
Scrandom sample matrix, to be computed by the matrix-(multiple)vector multiplication routine as A^T*Rc, or A^C*Rc. This will aready be allocated by the compression routine and should be Sc.rows() == this->cols() and Sc.cols() == Rc.cols().
Aelemelement extraction routine. This can be a functor, or a lambda function for instance.
IParameter in the element extraction routine. Set of row indices of elements to extract.
JParameter in the element extraction routine. Set of column indices of elements to extract.
BParameter in the element extraction routine. Matrix where to place extracted elements. This matrix will already be allocated. It will have B.rows() == I.size() and B.cols() == J.size().
optsobject containing a number of options for HSS compression
See also
DenseMatrix
HSSOptions

◆ compress_with_coordinates()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::compress_with_coordinates ( const DenseMatrix< real_t > &  coords,
const std::function< void(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J, DenseM_t &B)> &  Aelem,
const opts_t opts 
)

Initialize this HSS matrix as the compressed HSS representation. The compression uses nearest neighbor information (coordinates provided by the user). The HSS matrix should have been constructed with the proper sizes, i.e., rows() == A.rows() and cols() == A.cols().

Parameters
coordsmatrix with coordinates for the underlying geometry that defined the HSS matrix. This should be a d x n matrix (d rows, n columns), where d is the dimension of the coordinates and n is rows() and cols() of this matrix.
Aelemelement extraction routine. This can be a functor, or a lambda function for instance.
IParameter in the element extraction routine. Set of row indices of elements to extract.
JParameter in the element extraction routine. Set of column indices of elements to extract.
BParameter in the element extraction routine. Matrix where to place extracted elements. This matrix will already be allocated. It will have B.rows() == I.size() and B.cols() == J.size().
optsobject containing a number of options for HSS compression
See also
DenseMatrix
HSSOptions

◆ dense()

template<typename scalar_t >
DenseM_t strumpack::HSS::HSSMatrix< scalar_t >::dense ( ) const

Return a full/dense representation of this HSS matrix.

Returns
Dense matrix obtained from multiplying out the low-rank representation.

◆ draw()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::draw ( std::ostream &  of,
std::size_t  rlo = 0,
std::size_t  clo = 0 
) const
overridevirtual

Internal routine to draw this HSS matrix. Do not use this directly. Use HSS::draw.

See also
HSS::draw

Reimplemented from strumpack::HSS::HSSMatrixBase< scalar_t >.

◆ extract()

template<typename scalar_t >
DenseM_t strumpack::HSS::HSSMatrix< scalar_t >::extract ( const std::vector< std::size_t > &  I,
const std::vector< std::size_t > &  J 
) const

Compute a submatrix this(I,J) from this HSS matrix.

Parameters
Iset of row indices of elements to compute from this HSS matrix.
Iset of column indices of elements to compute from this HSS matrix.
Returns
Submatrix from this HSS matrix this(I,J)
See also
get, extract_add

◆ extract_add()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::extract_add ( const std::vector< std::size_t > &  I,
const std::vector< std::size_t > &  J,
DenseM_t B 
) const

Compute a submatrix this(I,J) from this HSS matrix and add it to a given matrix B.

Parameters
Iset of row indices of elements to compute from this HSS matrix.
Iset of column indices of elements to compute from this HSS matrix.
Theextracted submatrix this(I,J) will be added to this matrix. Should satisfy B.rows() == I.size() and B.cols() == J.size().
See also
get, extract

◆ factor()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::factor ( )
overridevirtual

Compute a ULV factorization of this matrix.

Reimplemented from strumpack::structured::StructuredMatrix< scalar_t >.

◆ forward_solve()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::forward_solve ( WorkSolve< scalar_t > &  w,
const DenseM_t b,
bool  partial 
) const
overridevirtual

Perform only the forward phase of the ULV linear solve. This is for advanced use only, typically to be used in combination with partial_factor. You should really just use factor/solve when possible.

Parameters
wtemporary working storage, to pass information from forward_solve to backward_solve
bon input, the right hand side vector, on output the intermediate solution. The vector b should be b.rows() == cols().
partialdenotes wether the matrix was fully or partially factored
See also
factor, partial_factor, backward_solve

Reimplemented from strumpack::HSS::HSSMatrixBase< scalar_t >.

◆ get()

template<typename scalar_t >
scalar_t strumpack::HSS::HSSMatrix< scalar_t >::get ( std::size_t  i,
std::size_t  j 
) const

Extract a single element this(i,j) from this HSS matrix. This is expensive and should not be used to compute multiple elements.

Parameters
iGlobal row index of element to compute.
jGlobal column index of element to compute.
Returns
this(i,j)
See also
extract

◆ levels()

template<typename scalar_t >
std::size_t strumpack::HSS::HSSMatrix< scalar_t >::levels ( ) const
overridevirtual

Return the number of levels in the HSS matrix.

Returns
Number of HSS levels (>= 1).

Implements strumpack::HSS::HSSMatrixBase< scalar_t >.

◆ memory()

template<typename scalar_t >
std::size_t strumpack::HSS::HSSMatrix< scalar_t >::memory ( ) const
overridevirtual

Return the total amount of memory used by this matrix, in bytes.

Returns
Memory usage in bytes.
See also
nonzeros

Implements strumpack::structured::StructuredMatrix< scalar_t >.

◆ mult()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::mult ( Trans  op,
const DenseM_t x,
DenseM_t y 
) const
overridevirtual

Multiply this HSS matrix with a dense matrix (vector), ie, compute y = op(this) * x. Overrides from the StructuredMatrix class method.

Parameters
opTranspose or complex conjugate
xright hand side matrix to multiply with, from the left, rows(x) == cols(op(this))
yresult of op(this) * b, cols(y) == cols(x), rows(r) = rows(op(this))
See also
applyC, HSS::apply_HSS

Reimplemented from strumpack::structured::StructuredMatrix< scalar_t >.

◆ nonzeros()

template<typename scalar_t >
std::size_t strumpack::HSS::HSSMatrix< scalar_t >::nonzeros ( ) const
overridevirtual

Return the total number of nonzeros stored by this matrix.

Returns
Nonzeros in the matrix representation.
See also
memory

Implements strumpack::structured::StructuredMatrix< scalar_t >.

◆ operator=() [1/2]

template<typename scalar_t >
HSSMatrix< scalar_t > & strumpack::HSS::HSSMatrix< scalar_t >::operator= ( const HSSMatrix< scalar_t > &  other)

Copy assignment operator. Copying an HSSMatrix can be an expensive operation.

Parameters
otherHSS matrix to be copied

◆ operator=() [2/2]

template<typename scalar_t >
HSSMatrix< scalar_t > & strumpack::HSS::HSSMatrix< scalar_t >::operator= ( HSSMatrix< scalar_t > &&  other)
default

Move assignment operator.

Parameters
otherHSS matrix to be moved from, will be emptied

◆ partial_factor()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::partial_factor ( )

Compute a partial ULV factorization of this matrix. Only the left child is factored. This is not similar to calling child(0)->factor(), except that the HSSFactors resulting from calling partial_factor can be used to compute the Schur complement.

See also
Schur_update, Schur_product_direct and Schur_product_indirect

◆ print_info()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::print_info ( std::ostream &  out = std::cout,
std::size_t  roff = 0,
std::size_t  coff = 0 
) const
overridevirtual

Print info about this HSS matrix, such as tree info, ranks, etc.

Parameters
outStream to print to, defaults to std::cout
roffRow offset of top left corner, defaults to 0. This is used to recursively print the tree, you can leave this at the default.
coffColumn offset of top left corner, defaults to 0. This is used to recursively print the tree, you can leave this at the default.

Implements strumpack::HSS::HSSMatrixBase< scalar_t >.

◆ rank()

template<typename scalar_t >
std::size_t strumpack::HSS::HSSMatrix< scalar_t >::rank ( ) const
overridevirtual

Return the maximum rank of this matrix over all low-rank compressed blocks.

Returns
Maximum rank.

Implements strumpack::structured::StructuredMatrix< scalar_t >.

◆ read()

template<typename scalar_t >
static HSSMatrix< scalar_t > strumpack::HSS::HSSMatrix< scalar_t >::read ( const std::string &  fname)
static

Read an HSSMatrix<scalar_t> from a binary file, called fname.

See also
write

◆ reset()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::reset ( )
override

Reset the matrix to an empty, 0 x 0 matrix, freeing up all it's memory.

◆ shift()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::shift ( scalar_t  sigma)
overridevirtual

Apply a shift to the diagonal of this matrix. Ie, this += sigma * I, with I the identity matrix. Call this after compression.

Parameters
sigmaShift to be applied to the diagonal.

Implements strumpack::HSS::HSSMatrixBase< scalar_t >.

◆ solve()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::solve ( DenseM_t b) const
overridevirtual

Solve a linear system with the ULV factorization of this HSSMatrix. The right hand side vector (or matrix) b is overwritten with the solution x of Ax=b.

Parameters
ULVULV factorization of this matrix
bon input, the right hand side vector, on output the solution of A x = b (with A this HSS matrix). The vector b should be b.rows() == cols().
See also
factor

Reimplemented from strumpack::structured::StructuredMatrix< scalar_t >.

◆ write()

template<typename scalar_t >
void strumpack::HSS::HSSMatrix< scalar_t >::write ( const std::string &  fname) const

Write this HSSMatrix<scalar_t> to a binary file, called fname.

See also
read

Friends And Related Function Documentation

◆ apply_HSS

template<typename scalar_t >
template<typename T >
void apply_HSS ( Trans  ta,
const HSSMatrix< T > &  a,
const DenseMatrix< T > &  b,
beta,
DenseMatrix< T > &  c 
)
friend
See also
HSS::apply_HSS

◆ draw

template<typename scalar_t >
template<typename T >
void draw ( const HSSMatrix< T > &  H,
const std::string &  name 
)
friend
See also
HSS::draw

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