H2Lib  3.0
Data Structures | Typedefs | Functions
harith

Algebraic operations with hierarchical matrices. More...

Data Structures

struct  _haccum
 Accumulator for H-matrix products. More...
 

Typedefs

typedef struct _haccum haccum
 Accumulator for H-matrix products.
 
typedef haccumphaccum
 Pointer to haccum object.
 
typedef const haccumpchaccum
 Pointer to constant haccum object.
 
typedef struct _hprodentry hprodentry
 Unhandled H-matrix products.
 
typedef hprodentryphprodentry
 Points to hprodentry object.
 

Functions

void trunc_rkmatrix (pctruncmode tm, real eps, prkmatrix r)
 Truncate an rkmatrix, $A \gets \operatorname{trunc}(A,\epsilon)$. More...
 
void add_amatrix_rkmatrix (field alpha, bool atrans, pcamatrix a, pctruncmode tm, real eps, prkmatrix b)
 Truncated addition of a matrix to a low-rank matrices in rkmatrix representation, $B \gets \operatorname{trunc}(B+\alpha A,\epsilon)$. More...
 
void add_rkmatrix_amatrix (field alpha, bool atrans, pcrkmatrix a, pamatrix b)
 Truncated addition of a full matrix to a low-rank matrix. More...
 
void add_hmatrix_amatrix (field alpha, bool atrans, pchmatrix a, pamatrix b)
 Addition of a hierarchical matrix to a matrix, $B \gets B+\alpha A$. More...
 
void add_rkmatrix (field alpha, pcrkmatrix src, pctruncmode tm, real eps, prkmatrix trg)
 Truncated addition of two low-rank matrices in rkmatrix representation, $A \gets \operatorname{trunc}(A+\alpha B,\epsilon)$. More...
 
void add_amatrix_hmatrix (field alpha, bool atrans, pcamatrix a, pctruncmode tm, real eps, phmatrix b)
 Truncated addition of two low-rank matrices in rkmatrix. More...
 
void add_lower_amatrix_hmatrix (field alpha, bool atrans, pcamatrix a, pctruncmode tm, real eps, phmatrix b)
 Locally truncated addition of a matrix in amatrix representation to the lower triangular part of a hierarchical matrix in hmatrix representation, $\operatorname{lower}(B) \gets \operatorname{blocktrunc}(\operatorname{lower}(B + \alpha A),\epsilon)$. More...
 
void add_rkmatrix_hmatrix (field alpha, pcrkmatrix a, pctruncmode tm, real eps, phmatrix b)
 Locally truncated addition of a low-rank matrix in rkmatrix representation to a hierarchical matrix in hmatrix representation, $B \gets \operatorname{blocktrunc}(B + \alpha A,\epsilon)$. More...
 
void add_lower_rkmatrix_hmatrix (field alpha, pcrkmatrix a, pctruncmode tm, real eps, phmatrix b)
 Locally truncated addition of a low-rank matrix in rkmatrix representation to the lower triangular part of a hierarchical matrix in hmatrix representation, $\operatorname{lower}(B) \gets \operatorname{blocktrunc}(\operatorname{lower}(B + \alpha A),\epsilon)$. More...
 
void add_hmatrix (field alpha, pchmatrix a, pctruncmode tm, real eps, phmatrix b)
 Locally truncated addition of a hierarchical matrix to a hierarchical matrix representation, $B \gets \operatorname{blocktrunc}(B + \alpha A,\epsilon)$. More...
 
phmatrix split_sub_amatrix (pcamatrix f, pccluster rc, pccluster cc, bool rsplit, bool csplit)
 Split an amatrix into submatrices. More...
 
phmatrix split_rkmatrix (pcrkmatrix r, pccluster rc, pccluster cc, bool rsplit, bool csplit, bool copy)
 Split an rkmatrix into submatrices. More...
 
phmatrix split_sub_rkmatrix (pcrkmatrix r, pccluster rc, pccluster cc, bool rsplit, bool csplit)
 Split an rkmatrix into submatrices referencing the original matrix. More...
 
void merge_rkmatrix (bool colmerge, pcrkmatrix src, pctruncmode tm, real eps, prkmatrix trg)
 Split a leaf hmatrix into submatrices. More...
 
prkmatrix merge_hmatrix_rkmatrix (pchmatrix s, pctruncmode tm, real eps)
 Merge a depth-one hmatrix containing only rkmatrix submatrices into a new rkmatrix, $A \gets \operatorname{succtrunc}\left(\begin{pmatrix} B_{11} & \ldots & B_{1m}\\ \vdots & \ddots & \vdots\\ B_{n1} & \ldots & B_{nm} \end{pmatrix},\epsilon\right)$. More...
 
void addmul_rkmatrix_amatrix_amatrix (field alpha, bool atrans, pcrkmatrix a, bool btrans, pcamatrix b, bool ctrans, pamatrix c)
 Multiply a matrix by an rkmatrix, $C \gets C + \alpha A B$. More...
 
void addmul_hmatrix_amatrix_amatrix (field alpha, bool atrans, pchmatrix a, bool btrans, pcamatrix bp, bool ctrans, pamatrix cp)
 Multiply a matrix by a hmatrix, $C \gets C + \alpha A B$. More...
 
void addmul_hmatrix (field alpha, bool xtrans, pchmatrix x, bool ytrans, pchmatrix y, pctruncmode tm, real eps, phmatrix z)
 Multiply two H-matrices, $Z \gets \operatorname{succtrunc}(Z + \alpha X Y,\epsilon)$. More...
 
void addmul_lower_hmatrix (field alpha, bool xtrans, pchmatrix x, bool ytrans, pchmatrix y, pctruncmode tm, real eps, phmatrix z)
 Multiply two H-matrices, computing only the lower triangular part of the result, $\operatorname{lower}(Z) \gets \operatorname{succtrunc}(\operatorname{lower}(Z + \alpha X Y),\epsilon)$. More...
 
void invert_hmatrix (phmatrix a, phmatrix work, pctruncmode tm, real eps)
 Invert an H-matrix, $A \gets \operatorname{succtrunc}(A^{-1},\epsilon)$. More...
 
phmatrix clone_lower_hmatrix (bool aunit, pchmatrix a)
 Create a duplicate of the lower triangular part of an H-matrix. More...
 
phmatrix clone_upper_hmatrix (bool aunit, pchmatrix a)
 Create a duplicate of the upper triangular part of an H-matrix. More...
 
void triangularinvmul_hmatrix_avector (bool alower, bool aunit, bool atrans, pchmatrix a, pavector xp)
 Solve a triangular system, $x \gets R^{-1} x$ or $x \gets R^{-*} x$. More...
 
void triangularsolve_hmatrix_avector (bool alower, bool aunit, bool atrans, pchmatrix a, pavector x)
 Solve a triangular system $x \gets A^{-1} x$ or $x \gets (A^*)^{-1} x$. More...
 
void triangularinvmul_hmatrix_amatrix (bool alower, bool aunit, bool atrans, pchmatrix a, bool xtrans, pamatrix xp)
 Solve a triangular system $X \gets A^{-1} X$ or $X \gets (A^*)^{-1} X$. More...
 
void triangularinvmul_hmatrix (bool alower, bool aunit, bool atrans, pchmatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix xp)
 Solve a triangular system, $X \gets \operatorname{succtrunc}(A^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(A^{-*} X, \epsilon)$. More...
 
void triangularinvmul_amatrix_hmatrix (bool alower, bool aunit, bool atrans, pcamatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix xp)
 Solve a triangular system, $X \gets \operatorname{succtrunc}(A^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(A^{-*} X, \epsilon)$. More...
 
void triangularmul_hmatrix_avector (bool alower, bool aunit, bool atrans, pchmatrix a, pavector xp)
 Evaluate a triangular system $x \gets A x$ or $x \gets A^* x$. More...
 
void triangulareval_hmatrix_avector (bool alower, bool aunit, bool atrans, pchmatrix a, pavector x)
 Evaluate a triangular system $x \gets A x$ or $x \gets A^* x$. More...
 
void triangularmul_hmatrix_amatrix (bool alower, bool aunit, bool atrans, pchmatrix a, bool xtrans, pamatrix xp)
 Evaluate a triangular system $X \gets A X$ or $X \gets A^* X$. More...
 
void triangularmul_hmatrix (bool alower, bool aunit, bool atrans, pchmatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix xp)
 Evaluate a triangular system, $X \gets \operatorname{succtrunc}(A X, \epsilon)$ or $X \gets \operatorname{succtrunc}(A^{*} X, \epsilon)$. More...
 
void lrdecomp_hmatrix (phmatrix a, pctruncmode tm, real eps)
 Compute the LR factorization, $A \approx L R$. More...
 
void lrsolve_n_hmatrix_avector (pchmatrix a, pavector x)
 Solve the linear systems $A x = b$ using the LR factorization provided by lrdecomp_hmatrix. More...
 
void lrsolve_t_hmatrix_avector (pchmatrix a, pavector x)
 Solve the linear systems $A^* x = b$ using the LR factorization provided by lrdecomp_hmatrix. More...
 
void lrsolve_hmatrix_avector (bool atrans, pchmatrix a, pavector x)
 Solve the linear systems $A x = b$ or $A^* x = b$ using the LR factorization provided by lrdecomp_hmatrix. More...
 
void lreval_n_hmatrix_avector (pchmatrix a, pavector x)
 Evaluate the linear system $x \gets A x$ using the LR factorization provided by lrdecomp_hmatrix. More...
 
void lreval_t_hmatrix_avector (pchmatrix a, pavector x)
 Evaluate the linear system $x \gets A^* x$ using the LR factorization provided by lrdecomp_hmatrix. More...
 
void lreval_hmatrix_avector (bool atrans, pchmatrix a, pavector x)
 Evaluate the linear system $x \gets A x$ or $x \gets A^* x$ using the LR factorization provided by lrdecomp_hmatrix. More...
 
void choldecomp_hmatrix (phmatrix a, pctruncmode tm, real eps)
 Compute the Cholesky factorization, $A \approx L L^*$. More...
 
void cholsolve_hmatrix_avector (pchmatrix a, pavector x)
 Solve the linear system $A x = b$ using the Cholesky factorization provided by choldecomp_hmatrix. More...
 
void choleval_hmatrix_avector (pchmatrix a, pavector x)
 Evaluate the linear system $x \gets A x$ using the Cholesky factorization provided by choldecomp_hmatrix. More...
 
phaccum new_haccum (phmatrix z, pctruncmode tm, real eps)
 Create a new H-matrix accumulator. More...
 
void del_haccum (phaccum ha)
 Delete an accumulator. More...
 
void addproduct_haccum (field alpha, bool xtrans, pchmatrix x, bool ytrans, pchmatrix y, phaccum za)
 Add a product $\alpha X Y$ to an accumulator. More...
 
phaccumsplit_haccum (phmatrix z, phaccum za)
 Create H-matrix accumulators for submatrices. More...
 
void flush_haccum (phaccum ha)
 Flush an H-matrix accumulator, i.e., add all accumulated products to the target matrix and clear the accumulator. More...
 
void add_amatrix_destructive_rkmatrix (field alpha, bool atrans, pamatrix a, pctruncmode tm, real eps, prkmatrix b)
 Truncated addition of a matrix to a low-rank matrices in rkmatrix representation, $B \gets \operatorname{trunc}(B+\alpha A,\epsilon)$. More...
 
void add_amatrix_destructive_hmatrix (field alpha, bool atrans, pamatrix a, pctruncmode tm, real eps, phmatrix b)
 Locally truncated addition of a matrix in amatrix representation to a hierarchical matrix in hmatrix representation, $B \gets \operatorname{blocktrunc}(B + \alpha A,\epsilon)$. More...
 
void addmul2_hmatrix (field alpha, bool xtrans, pchmatrix x, bool ytrans, pchmatrix y, pctruncmode tm, real eps, phmatrix z)
 Multiply two H-matrices using accumulators, $Z \gets \operatorname{succtrunc}(Z + \alpha X Y,\epsilon)$. More...
 
void lowersolve_haccum (bool aunit, bool atrans, pchmatrix a, bool xtrans, phaccum xa)
 Solve a lower triangular system using accumulators, $X \gets \operatorname{succtrunc}(L^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(L^{-*} X, \epsilon)$. More...
 
void lowersolve2_hmatrix_hmatrix (bool aunit, bool atrans, pchmatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix x)
 Solve a lower triangular system using accumulators, $X \gets \operatorname{succtrunc}(L^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(L^{-*} X, \epsilon)$. More...
 
void uppersolve_haccum (bool aunit, bool atrans, pchmatrix a, bool xtrans, phaccum xa)
 Solve an upper triangular system using accumulators, $X \gets \operatorname{succtrunc}(R^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(R^{-*} X, \epsilon)$. More...
 
void uppersolve2_hmatrix_hmatrix (bool aunit, bool atrans, pchmatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix x)
 Solve an upper triangular system using accumulators, $X \gets \operatorname{succtrunc}(R^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(R^{-*} X, \epsilon)$. More...
 
void lrdecomp_haccum (phaccum aa)
 Compute the LR factorization using accumulators, $A \approx L R$. More...
 
void lrdecomp2_hmatrix (phmatrix a, pctruncmode tm, real eps)
 Compute the LR factorization using accumulators, $A \approx L R$. More...
 
void choldecomp_haccum (phaccum aa)
 Compute the Cholesky factorization using accumulators, $A \approx L L^*$. More...
 
void choldecomp2_hmatrix (phmatrix a, pctruncmode tm, real eps)
 Compute the Cholesky factorization using accumulators, $A \approx L L^*$. More...
 

Detailed Description

Algebraic operations with hierarchical matrices.

Function Documentation

void add_amatrix_destructive_hmatrix ( field  alpha,
bool  atrans,
pamatrix  a,
pctruncmode  tm,
real  eps,
phmatrix  b 
)

Locally truncated addition of a matrix in amatrix representation to a hierarchical matrix in hmatrix representation, $B \gets \operatorname{blocktrunc}(B + \alpha A,\epsilon)$.

Attention
The matrix $A$ is overwritten by intermediate results.
Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aMatrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget hmatrix $B$.
void add_amatrix_destructive_rkmatrix ( field  alpha,
bool  atrans,
pamatrix  a,
pctruncmode  tm,
real  eps,
prkmatrix  b 
)

Truncated addition of a matrix to a low-rank matrices in rkmatrix representation, $B \gets \operatorname{trunc}(B+\alpha A,\epsilon)$.

Attention
The matrix $A$ is overwritten by intermediate results.
Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aMatrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget matrix $B$.
void add_amatrix_hmatrix ( field  alpha,
bool  atrans,
pcamatrix  a,
pctruncmode  tm,
real  eps,
phmatrix  b 
)

Truncated addition of two low-rank matrices in rkmatrix.

*Locally truncated addition of a matrix in amatrix representation to a hierarchical matrix in hmatrix representation, $B \gets \operatorname{blocktrunc}(B + \alpha A,\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aMatrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget hmatrix $B$.
void add_amatrix_rkmatrix ( field  alpha,
bool  atrans,
pcamatrix  a,
pctruncmode  tm,
real  eps,
prkmatrix  b 
)

Truncated addition of a matrix to a low-rank matrices in rkmatrix representation, $B \gets \operatorname{trunc}(B+\alpha A,\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aMatrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget matrix $B$.
void add_hmatrix ( field  alpha,
pchmatrix  a,
pctruncmode  tm,
real  eps,
phmatrix  b 
)

Locally truncated addition of a hierarchical matrix to a hierarchical matrix representation, $B \gets \operatorname{blocktrunc}(B + \alpha A,\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
ahmatrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget hmatrix $B$.
void add_hmatrix_amatrix ( field  alpha,
bool  atrans,
pchmatrix  a,
pamatrix  b 
)

Addition of a hierarchical matrix to a matrix, $B \gets B+\alpha A$.

Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aMatrix $A$.
bTarget matrix $B$.
void add_lower_amatrix_hmatrix ( field  alpha,
bool  atrans,
pcamatrix  a,
pctruncmode  tm,
real  eps,
phmatrix  b 
)

Locally truncated addition of a matrix in amatrix representation to the lower triangular part of a hierarchical matrix in hmatrix representation, $\operatorname{lower}(B) \gets \operatorname{blocktrunc}(\operatorname{lower}(B + \alpha A),\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aMatrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget hmatrix $B$.
void add_lower_rkmatrix_hmatrix ( field  alpha,
pcrkmatrix  a,
pctruncmode  tm,
real  eps,
phmatrix  b 
)

Locally truncated addition of a low-rank matrix in rkmatrix representation to the lower triangular part of a hierarchical matrix in hmatrix representation, $\operatorname{lower}(B) \gets \operatorname{blocktrunc}(\operatorname{lower}(B + \alpha A),\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
aLow-rank matrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget hmatrix $B$.
void add_rkmatrix ( field  alpha,
pcrkmatrix  src,
pctruncmode  tm,
real  eps,
prkmatrix  trg 
)

Truncated addition of two low-rank matrices in rkmatrix representation, $A \gets \operatorname{trunc}(A+\alpha B,\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
srcSource matrix $B$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
trgTarget matrix $A$.
void add_rkmatrix_amatrix ( field  alpha,
bool  atrans,
pcrkmatrix  a,
pamatrix  b 
)

Truncated addition of a full matrix to a low-rank matrix.

*Addition of a low-rank matrix to a matrix, $B \gets B+\alpha A$.

Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aMatrix $A$.
bTarget matrix $B$.
void add_rkmatrix_hmatrix ( field  alpha,
pcrkmatrix  a,
pctruncmode  tm,
real  eps,
phmatrix  b 
)

Locally truncated addition of a low-rank matrix in rkmatrix representation to a hierarchical matrix in hmatrix representation, $B \gets \operatorname{blocktrunc}(B + \alpha A,\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
aLow-rank matrix $A$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
bTarget hmatrix $B$.
void addmul2_hmatrix ( field  alpha,
bool  xtrans,
pchmatrix  x,
bool  ytrans,
pchmatrix  y,
pctruncmode  tm,
real  eps,
phmatrix  z 
)

Multiply two H-matrices using accumulators, $Z \gets \operatorname{succtrunc}(Z + \alpha X Y,\epsilon)$.

Compared to the standard implementation addmul_hmatrix, this algorithm uses H-matrix accumulators to reduce the number of singular value decompositions.

Parameters
alphaScaling factor $\alpha$.
xtransSet if $X^*$ is to be used instead of $X$.
xHierarchical matrix $X$.
ytransSet if $Y^*$ is to be used instead of $Y$.
yHierarchical matrix $Y$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
zTarget matrix $Z$.
void addmul_hmatrix ( field  alpha,
bool  xtrans,
pchmatrix  x,
bool  ytrans,
pchmatrix  y,
pctruncmode  tm,
real  eps,
phmatrix  z 
)

Multiply two H-matrices, $Z \gets \operatorname{succtrunc}(Z + \alpha X Y,\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
xtransSet if $X^*$ is to be used instead of $X$.
xHierarchical matrix $X$.
ytransSet if $Y^*$ is to be used instead of $Y$.
yHierarchical matrix $Y$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
zTarget matrix $Z$.
void addmul_hmatrix_amatrix_amatrix ( field  alpha,
bool  atrans,
pchmatrix  a,
bool  btrans,
pcamatrix  bp,
bool  ctrans,
pamatrix  cp 
)

Multiply a matrix by a hmatrix, $C \gets C + \alpha A B$.

Remarks
The matrices $C$ and $B$ are assumed to be ordered matching the row and column cluster trees of $A$, respectively.
Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aHierarchical matrix $A$.
btransSet if $B^*$ is to be used instead of $B$.
bpMatrix $B$.
ctransSet if $C^*$ is to be used instead of $C$.
cpTarget matrix $C$.
void addmul_lower_hmatrix ( field  alpha,
bool  xtrans,
pchmatrix  x,
bool  ytrans,
pchmatrix  y,
pctruncmode  tm,
real  eps,
phmatrix  z 
)

Multiply two H-matrices, computing only the lower triangular part of the result, $\operatorname{lower}(Z) \gets \operatorname{succtrunc}(\operatorname{lower}(Z + \alpha X Y),\epsilon)$.

Parameters
alphaScaling factor $\alpha$.
xtransSet if $X^*$ is to be used instead of $X$.
xHierarchical matrix $X$.
ytransSet if $Y^*$ is to be used instead of $Y$.
yHierarchical matrix $Y$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
zTarget matrix $Z$.
void addmul_rkmatrix_amatrix_amatrix ( field  alpha,
bool  atrans,
pcrkmatrix  a,
bool  btrans,
pcamatrix  b,
bool  ctrans,
pamatrix  c 
)

Multiply a matrix by an rkmatrix, $C \gets C + \alpha A B$.

Parameters
alphaScaling factor $\alpha$.
atransSet if $A^*$ is to be used instead of $A$.
aLow-rank matrix $A$.
btransSet if $B^*$ is to be used instead of $B$.
bMatrix $B$.
ctransSet if $C^*$ is to be used instead of $C$.
cTarget matrix $C$.
void addproduct_haccum ( field  alpha,
bool  xtrans,
pchmatrix  x,
bool  ytrans,
pchmatrix  y,
phaccum  za 
)

Add a product $\alpha X Y$ to an accumulator.

Parameters
alphaScaling factor.
xtransSet if $X^*$ is to be used instead of $X$.
xH-matrix $X$.
ytransSet if $Y^*$ is to be used instead of $Y$.
yH-matrix $Y$.
zaAccumulator for the result.
void choldecomp2_hmatrix ( phmatrix  a,
pctruncmode  tm,
real  eps 
)

Compute the Cholesky factorization using accumulators, $A \approx L L^*$.

The lower triangular part of $L$ is stored in the lower triangular part of the source matrix.

The strictly upper triangular part of the source matrix is not used.

Parameters
aSource matrix $A$, lower triangular part will be overwritten by $L$.
tmTruncation mode.
epsTruncation accuracy.
void choldecomp_haccum ( phaccum  aa)

Compute the Cholesky factorization using accumulators, $A \approx L L^*$.

The lower triangular part of $L$ is stored in the lower triangular part of the source matrix.

The strictly upper triangular part of the source matrix is not used.

Parameters
aaAccumulator representation of $A$, lower triangular part will be overwritten
void choldecomp_hmatrix ( phmatrix  a,
pctruncmode  tm,
real  eps 
)

Compute the Cholesky factorization, $A \approx L L^*$.

The lower triangular part of $L$ is stored in the lower triangular part of the source matrix.

The strictly upper triangular part of the source matrix is not used.

Parameters
aSource matrix $A$, lower triangular part will be overwritten by $L$.
tmTruncation mode.
epsTruncation accuracy.
void choleval_hmatrix_avector ( pchmatrix  a,
pavector  x 
)

Evaluate the linear system $x \gets A x$ using the Cholesky factorization provided by choldecomp_hmatrix.

Parameters
aMatrix containing the Cholesky factorization in the form returned by choldecomp_hmatrix.
xRight-hand side vector $x$, that gets overwritten by the result.
void cholsolve_hmatrix_avector ( pchmatrix  a,
pavector  x 
)

Solve the linear system $A x = b$ using the Cholesky factorization provided by choldecomp_hmatrix.

Parameters
aMatrix containing the Cholesky factorization in the form returned by choldecomp_hmatrix.
xRight-hand side vector $b$, will be overwritten by solution vector $x$.
phmatrix clone_lower_hmatrix ( bool  aunit,
pchmatrix  a 
)

Create a duplicate of the lower triangular part of an H-matrix.

Parameters
aunitSet if the duplicate has to have unit diagonal.
aSource matrix $A$.
Returns
Lower triangular matrix containing the lower triangular part of $A$, with unit diagonal if aunit==true.
phmatrix clone_upper_hmatrix ( bool  aunit,
pchmatrix  a 
)

Create a duplicate of the upper triangular part of an H-matrix.

Parameters
aunitSet if the duplicate has to have unit diagonal.
aSource matrix $A$.
Returns
Upper triangular matrix containing the upper triangular part of $A$, with unit diagonal if aunit==true.
void del_haccum ( phaccum  ha)

Delete an accumulator.

Parameters
haAccumulator to be deleted.
void flush_haccum ( phaccum  ha)

Flush an H-matrix accumulator, i.e., add all accumulated products to the target matrix and clear the accumulator.

Parameters
haH-matrix accumulator.
void invert_hmatrix ( phmatrix  a,
phmatrix  work,
pctruncmode  tm,
real  eps 
)

Invert an H-matrix, $A \gets \operatorname{succtrunc}(A^{-1},\epsilon)$.

Parameters
aMatrix $A$, gets overwritten with $A^{-1}$.
workAuxiliary storage, has to be zero and be of the same structure as a. Consider using clonestructure_hmatrix.
tmTruncation mode.
epsTruncation accuracy.
void lowersolve2_hmatrix_hmatrix ( bool  aunit,
bool  atrans,
pchmatrix  a,
pctruncmode  tm,
real  eps,
bool  xtrans,
phmatrix  x 
)

Solve a lower triangular system using accumulators, $X \gets \operatorname{succtrunc}(L^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(L^{-*} X, \epsilon)$.

Parameters
aunitSet if $L$ has unit diagonal.
atransSet if $L^*$ instead of $L$ is to be used.
aMatrix containing the lower triangular part of $L$. The remainder is not used.
tmTruncation mode.
epsTruncation accuracy.
xtransset if $X^*$ instead of $X$ is to be used.
xInput matrix $X$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void lowersolve_haccum ( bool  aunit,
bool  atrans,
pchmatrix  a,
bool  xtrans,
phaccum  xa 
)

Solve a lower triangular system using accumulators, $X \gets \operatorname{succtrunc}(L^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(L^{-*} X, \epsilon)$.

Parameters
aunitSet if $L$ has unit diagonal.
atransSet if $L^*$ instead of $L$ is to be used.
aMatrix containing the lower triangular part of $L$. The remainder is not used.
xtransset if $X^*$ instead of $X$ is to be used.
xaAccumulator representation of $X$, will be overwritten by the solution.
void lrdecomp2_hmatrix ( phmatrix  a,
pctruncmode  tm,
real  eps 
)

Compute the LR factorization using accumulators, $A \approx L R$.

The lower triangular part $L$ has unit diagonal, only its strict lower triangular part is stored in the strict lower triangular part of the source matrix.

The upper triangular part is stored in the upper triangular part of the source matrix.

Parameters
aSource matrix $A$, will be overwritten with $L$ and $R$.
tmTruncation mode.
epsTruncation accuracy.
void lrdecomp_haccum ( phaccum  aa)

Compute the LR factorization using accumulators, $A \approx L R$.

The lower triangular part $L$ has unit diagonal, only its strict lower triangular part is stored in the strict lower triangular part of the source matrix.

The upper triangular part is stored in the upper triangular part of the source matrix.

Parameters
aaAccumulator representation of $A$, will be overwritten with $L$ and $R$.
void lrdecomp_hmatrix ( phmatrix  a,
pctruncmode  tm,
real  eps 
)

Compute the LR factorization, $A \approx L R$.

The lower triangular part $L$ has unit diagonal, only its strict lower triangular part is stored in the strict lower triangular part of the source matrix.

The upper triangular part is stored in the upper triangular part of the source matrix.

Parameters
aSource matrix $A$, will be overwritten by $L$ and $R$.
tmTruncation mode.
epsTruncation accuracy.
void lreval_hmatrix_avector ( bool  atrans,
pchmatrix  a,
pavector  x 
)

Evaluate the linear system $x \gets A x$ or $x \gets A^* x$ using the LR factorization provided by lrdecomp_hmatrix.

Parameters
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the LR factorization in the form returned by lrdecomp_hmatrix.
xRight-hand side vector $x$, that gets overwritten by the result.
void lreval_n_hmatrix_avector ( pchmatrix  a,
pavector  x 
)

Evaluate the linear system $x \gets A x$ using the LR factorization provided by lrdecomp_hmatrix.

Parameters
aMatrix containing the LR factorization in the form returned by lrdecomp_hmatrix.
xRight-hand side vector $x$, that gets overwritten by the result.
void lreval_t_hmatrix_avector ( pchmatrix  a,
pavector  x 
)

Evaluate the linear system $x \gets A^* x$ using the LR factorization provided by lrdecomp_hmatrix.

Parameters
aMatrix containing the LR factorization in the form returned by lrdecomp_hmatrix.
xRight-hand side vector $x$, that gets overwritten by the result.
void lrsolve_hmatrix_avector ( bool  atrans,
pchmatrix  a,
pavector  x 
)

Solve the linear systems $A x = b$ or $A^* x = b$ using the LR factorization provided by lrdecomp_hmatrix.

Parameters
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the LR factorization in the form returned by lrdecomp_hmatrix.
xRight-hand side vector $b$, will be overwritten by solution vector $x$.
void lrsolve_n_hmatrix_avector ( pchmatrix  a,
pavector  x 
)

Solve the linear systems $A x = b$ using the LR factorization provided by lrdecomp_hmatrix.

Parameters
aMatrix containing the LR factorization in the form returned by lrdecomp_hmatrix.
xRight-hand side vector $b$, will be overwritten by solution vector $x$.
void lrsolve_t_hmatrix_avector ( pchmatrix  a,
pavector  x 
)

Solve the linear systems $A^* x = b$ using the LR factorization provided by lrdecomp_hmatrix.

Parameters
aMatrix containing the LR factorization in the form returned by lrdecomp_hmatrix.
xRight-hand side vector $b$, will be overwritten by solution vector $x$.
prkmatrix merge_hmatrix_rkmatrix ( pchmatrix  s,
pctruncmode  tm,
real  eps 
)

Merge a depth-one hmatrix containing only rkmatrix submatrices into a new rkmatrix, $A \gets \operatorname{succtrunc}\left(\begin{pmatrix} B_{11} & \ldots & B_{1m}\\ \vdots & \ddots & \vdots\\ B_{n1} & \ldots & B_{nm} \end{pmatrix},\epsilon\right)$.

Parameters
sSource matrix $B$
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
Returns
rkmatrix approximating the block matrix.
void merge_rkmatrix ( bool  colmerge,
pcrkmatrix  src,
pctruncmode  tm,
real  eps,
prkmatrix  trg 
)

Split a leaf hmatrix into submatrices.

** Split a leaf hmatrix into submatrices referencing Block-merge two rk matrices into a new rk matrix, $A \gets \operatorname{trunc}\left(\begin{pmatrix} A & B \end{pmatrix},\epsilon\right)$ or $A \gets \operatorname{trunc}\left(\begin{pmatrix} A\\ B \end{pmatrix},\epsilon\right)$.

If colmerge is not set, $A$ and $B$ have to have the same number of rows, otherwise the same number of columns.

Attention
The matrix $A$ is the last parameter, although it is placed in the first row or column.
Parameters
colmergeSet to merge into a column instead of a row.
srcSource matrix $B$.
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
trgTarget matrix $A$, will be overwritten by result.
phaccum new_haccum ( phmatrix  z,
pctruncmode  tm,
real  eps 
)

Create a new H-matrix accumulator.

Parameters
zTarget matrix.
tmTruncation mode.
epsTruncation accuracy.
Returns
Accumulator for z
phaccum* split_haccum ( phmatrix  z,
phaccum  za 
)

Create H-matrix accumulators for submatrices.

Parameters
zTarget matrix, row and column cluster have to conincide with ha->z->rc and ha->z->cc.
zaH-matrix accumulator.
Returns
Array containing accumulators for the submatrices of z, numbered as in z->son.
phmatrix split_rkmatrix ( pcrkmatrix  r,
pccluster  rc,
pccluster  cc,
bool  rsplit,
bool  csplit,
bool  copy 
)

Split an rkmatrix into submatrices.

Parameters
rOriginal matrix.
rcRow cluster.
ccColumn cluster.
rsplitSet to split in row direction.
csplitSet to split in column direction.
copySet to copy values into new matrix, otherwise it will be zero.
Returns
Block matrix of depth one.
phmatrix split_sub_amatrix ( pcamatrix  f,
pccluster  rc,
pccluster  cc,
bool  rsplit,
bool  csplit 
)

Split an amatrix into submatrices.

*Split an amatrix into submatrices referencing the original matrix.

Parameters
fOriginal matrix.
rcRow cluster.
ccColumn cluster.
rsplitSet to split in row direction.
csplitSet to split in column direction.
Returns
Block matrix of depth one.
phmatrix split_sub_rkmatrix ( pcrkmatrix  r,
pccluster  rc,
pccluster  cc,
bool  rsplit,
bool  csplit 
)

Split an rkmatrix into submatrices referencing the original matrix.

Attention
Since the blocks of the result contain only references to the original matrix, they cannot be resized, and changing them changes the other submatrices in the same row and column. It's probably best to treat the submatrices as constant.
Parameters
rOriginal matrix.
rcRow cluster.
ccColumn cluster.
rsplitSet to split in row direction.
csplitSet to split in column direction.
Returns
Block matrix of depth one.
void triangulareval_hmatrix_avector ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
pavector  x 
)

Evaluate a triangular system $x \gets A x$ or $x \gets A^* x$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to be upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the lower or upper triangular part of $A$. The remainder of the matrix is not used.
xInput vector $x$, will be overwritten by the solution.
void triangularinvmul_amatrix_hmatrix ( bool  alower,
bool  aunit,
bool  atrans,
pcamatrix  a,
pctruncmode  tm,
real  eps,
bool  xtrans,
phmatrix  xp 
)

Solve a triangular system, $X \gets \operatorname{succtrunc}(A^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(A^{-*} X, \epsilon)$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to by upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the lower or upper triangular part of $A$. The remainder of the matrix is not used.
tmTruncation mode.
epsTruncation accuracy.
xtransset if $X^*$ instead of $X$ is to be used.
xpInput matrix $X$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void triangularinvmul_hmatrix ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
pctruncmode  tm,
real  eps,
bool  xtrans,
phmatrix  xp 
)

Solve a triangular system, $X \gets \operatorname{succtrunc}(A^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(A^{-*} X, \epsilon)$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to be upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the triangular part of $A$. The remainder is not used.
tmTruncation mode.
epsTruncation accuracy.
xtransset if $X^*$ instead of $X$ is to be used.
xpInput matrix $X$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void triangularinvmul_hmatrix_amatrix ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
bool  xtrans,
pamatrix  xp 
)

Solve a triangular system $X \gets A^{-1} X$ or $X \gets (A^*)^{-1} X$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to by upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the lower or upper triangular part of $A$. The remainder of the matrix is not used.
xtransset if $X^*$ instead of $X$ is to be used.
xpInput matrix $X$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void triangularinvmul_hmatrix_avector ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
pavector  xp 
)

Solve a triangular system, $x \gets R^{-1} x$ or $x \gets R^{-*} x$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to be upper triangular.
aunitSet if $R$ has unit diagonal.
atransSet if $R^*$ instead of $R$ is to be used.
aMatrix containing the upper triangular part of $R$. The strictly lower triangular part is not used.
xpInput vector $x$, will be overwritten by the solution. The vector is assumed to be in cluster numbering according to a->cc respectively a->rc.
void triangularmul_hmatrix ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
pctruncmode  tm,
real  eps,
bool  xtrans,
phmatrix  xp 
)

Evaluate a triangular system, $X \gets \operatorname{succtrunc}(A X, \epsilon)$ or $X \gets \operatorname{succtrunc}(A^{*} X, \epsilon)$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to by upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the lower or upper triangular part of $A$. The remainder of the matrix is not used.
tmTruncation mode.
epsTruncation accuracy.
xtransset if $X^*$ instead of $X$ is to be used.
xpInput matrix $X$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void triangularmul_hmatrix_amatrix ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
bool  xtrans,
pamatrix  xp 
)

Evaluate a triangular system $X \gets A X$ or $X \gets A^* X$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to by upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the lower or upper triangular part of $A$. The remainder of the matrix is not used.
xtransset if $X^*$ instead of $X$ is to be used.
xpInput matrix $X$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void triangularmul_hmatrix_avector ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
pavector  xp 
)

Evaluate a triangular system $x \gets A x$ or $x \gets A^* x$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to be upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the lower or upper triangular part of $A$. The remainder of the matrix is not used.
xpInput vector $x$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void triangularsolve_hmatrix_avector ( bool  alower,
bool  aunit,
bool  atrans,
pchmatrix  a,
pavector  x 
)

Solve a triangular system $x \gets A^{-1} x$ or $x \gets (A^*)^{-1} x$.

Parameters
alowerSet if $A$ is lower triangular, otherwise it is assumed to be upper triangular.
aunitSet if $A$ has unit diagonal.
atransSet if $A^*$ instead of $A$ is to be used.
aMatrix containing the lower or upper triangular part of $A$. The remainder of the matrix is not used.
xInput vector $x$, will be overwritten by the solution.
void trunc_rkmatrix ( pctruncmode  tm,
real  eps,
prkmatrix  r 
)

Truncate an rkmatrix, $A \gets \operatorname{trunc}(A,\epsilon)$.

Parameters
tmTruncation mode.
epsTruncation accuracy $\epsilon$.
rSource matrix, will be overwritten by truncated matrix.
void uppersolve2_hmatrix_hmatrix ( bool  aunit,
bool  atrans,
pchmatrix  a,
pctruncmode  tm,
real  eps,
bool  xtrans,
phmatrix  x 
)

Solve an upper triangular system using accumulators, $X \gets \operatorname{succtrunc}(R^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(R^{-*} X, \epsilon)$.

Parameters
aunitSet if $R$ has unit diagonal.
atransSet if $R^*$ instead of $R$ is to be used.
aMatrix containing the lower triangular part of $R$. The remainder is not used.
tmTruncation mode.
epsTruncation accuracy.
xtransset if $X^*$ instead of $X$ is to be used.
xInput matrix $X$, will be overwritten by the solution. The rows are assumed to be in cluster numbering according to a->rc.
void uppersolve_haccum ( bool  aunit,
bool  atrans,
pchmatrix  a,
bool  xtrans,
phaccum  xa 
)

Solve an upper triangular system using accumulators, $X \gets \operatorname{succtrunc}(R^{-1} X, \epsilon)$ or $X \gets \operatorname{succtrunc}(R^{-*} X, \epsilon)$.

Parameters
aunitSet if $R$ has unit diagonal.
atransSet if $R^*$ instead of $R$ is to be used.
aMatrix containing the upper triangular part of $R$. The remainder is not used.
xtransset if $X^*$ instead of $X$ is to be used.
xaAccumulator representation of $X$, will be overwritten by the solution.