H2Lib  3.0
harith.h
Go to the documentation of this file.
1 /* ------------------------------------------------------------
2  * This is the file "harith.h" of the H2Lib package.
3  * All rights reserved, Steffen Boerm 2014
4  * ------------------------------------------------------------ */
5 
10 #ifndef HARITH_H
11 #define HARITH_H
12 
13 #include "hmatrix.h"
14 #include "truncation.h"
15 
20 /* ------------------------------------------------------------
21  Truncation of an rkmatrix
22  ------------------------------------------------------------ */
23 
30 HEADER_PREFIX void
32 
33 /* ------------------------------------------------------------
34  * Matrix addition
35  * ------------------------------------------------------------ */
36 
47 HEADER_PREFIX void
48 add_amatrix_rkmatrix(field alpha, bool atrans, pcamatrix a, pctruncmode tm,
49  real eps, prkmatrix b);
50 
52 // * in @ref rkmatrix representation,
53 // * @f$B \gets \operatorname{trunc}(B+\alpha A,\epsilon)@f$.
54 // *
55 // * This version uses a rank-revealing QR factorization instead
56 // * of the standard singular value decomposition. This method
57 // * may lead to larger ranks, but may also be faster in some
58 // * applications.
59 // *
60 // * @param alpha Scaling factor @f$\alpha@f$.
61 // * @param a Source matrix @f$A@f$.
62 // * @param atrans Set if @f$A^*@f$ is to be used instead of @f$A@f$.
63 // * @param tm Truncation mode.
64 // * @param eps Truncation accuracy @f$\epsilon@f$.
65 // * @param b Target matrix @f$B@f$. */
66 //HEADER_PREFIX void
67 //add2_amatrix_rkmatrix(field alpha, bool atrans, pcamatrix a,
68 // pctruncmode tm, real eps, prkmatrix b);
69 
77 HEADER_PREFIX void
78 add_rkmatrix_amatrix(field alpha, bool atrans, pcrkmatrix a, pamatrix b);
79 
87 HEADER_PREFIX void
88 add_hmatrix_amatrix(field alpha, bool atrans, pchmatrix a, pamatrix b);
89 
99 HEADER_PREFIX void
100 add_rkmatrix(field alpha, pcrkmatrix src, pctruncmode tm, real eps,
101  prkmatrix trg);
102 
104 // * representation,
105 // * @f$B \gets \operatorname{trunc}(B+\alpha A,\epsilon)@f$.
106 // *
107 // * This version uses a rank-revealing QR factorization instead
108 // * of the standard singular value decomposition. This method
109 // * may lead to larger ranks, but may also be faster in some
110 // * applications.
111 // *
112 // * @param alpha Scaling factor @f$\alpha@f$.
113 // * @param a Source matrix @f$A@f$.
114 // * @param tm Truncation mode.
115 // * @param eps Truncation accuracy @f$\epsilon@f$.
116 // * @param b Target matrix @f$B@f$. */
117 //HEADER_PREFIX void
118 //add2_rkmatrix(field alpha, pcrkmatrix a,
119 // pctruncmode tm, real eps, prkmatrix b);
120 
131 HEADER_PREFIX void
132 add_amatrix_hmatrix(field alpha, bool atrans, pcamatrix a, pctruncmode tm,
133  real eps, phmatrix b);
134 
146 HEADER_PREFIX void
147 add_lower_amatrix_hmatrix(field alpha, bool atrans, pcamatrix a, pctruncmode tm,
148  real eps, phmatrix b);
149 
159 HEADER_PREFIX void
161  phmatrix b);
162 
173 HEADER_PREFIX void
175  phmatrix b);
176 
186 HEADER_PREFIX void
187 add_hmatrix(field alpha, pchmatrix a, pctruncmode tm, real eps, phmatrix b);
188 
189 /* ------------------------------------------------------------
190  * Splitting and merging H-matrices
191  * ------------------------------------------------------------ */
192 
194 // *
195 // * @param f Original matrix.
196 // * @param rc Row cluster.
197 // * @param cc Column cluster.
198 // * @param rsplit Set to split in row direction.
199 // * @param csplit Set to split in column direction.
200 // * @param copy Set to copy values into new matrix, otherwise it
201 // * will be zero.
202 // * @returns Block matrix of depth one. */
203 //HEADER_PREFIX phmatrix
204 //split_amatrix(pcamatrix f, pccluster rc, pccluster cc,
205 // bool rsplit, bool csplit, bool copy);
216 split_sub_amatrix(pcamatrix f, pccluster rc, pccluster cc, bool rsplit,
217  bool csplit);
218 
230 split_rkmatrix(pcrkmatrix r, pccluster rc, pccluster cc, bool rsplit,
231  bool csplit, bool copy);
232 
248 split_sub_rkmatrix(pcrkmatrix r, pccluster rc, pccluster cc, bool rsplit,
249  bool csplit);
250 
252 // *
253 // * @param hm Leaf @ref hmatrix.
254 // * @param rsplit Set to split in row direction.
255 // * @param csplit Set to split in column direction.
256 // * @param copy Set to copy values into new matrix, otherwise it
257 // * will be zero.
258 // * @returns Block matrix of depth one or zero according to
259 // * <tt>rsplit</tt> and <tt>csplit</tt>. */
260 //HEADER_PREFIX phmatrix
261 //split_hmatrix(phmatrix hm, bool rsplit, bool csplit, bool copy);
262 //
264 // * the original matrix.
265 // *
266 // * @attention Since the blocks of the result contain only references
267 // * to the original matrix, they cannot be resized, and changing
268 // * them changes the original matrix.
269 // *
270 // * @param hm Leaf @ref hmatrix.
271 // * @param rsplit Set to split in row direction.
272 // * @param csplit Set to split in column direction.
273 // * @returns Block matrix of depth one or zero according to
274 // * <tt>rsplit</tt> and <tt>csplit</tt>. */
275 //HEADER_PREFIX phmatrix
276 //split_sub_hmatrix(phmatrix hm,
277 // bool rsplit, bool csplit);
278 
294 HEADER_PREFIX void
295 merge_rkmatrix(bool colmerge, pcrkmatrix src, pctruncmode tm, real eps,
296  prkmatrix trg);
297 
310 
311 /* ------------------------------------------------------------
312  * Matrix multiplication
313  * ------------------------------------------------------------ */
314 
325 HEADER_PREFIX void
326 addmul_rkmatrix_amatrix_amatrix(field alpha, bool atrans, pcrkmatrix a,
327  bool btrans, pcamatrix b, bool ctrans, pamatrix c);
328 
343 HEADER_PREFIX void
344 addmul_hmatrix_amatrix_amatrix(field alpha, bool atrans, pchmatrix a,
345  bool btrans, pcamatrix bp, bool ctrans, pamatrix cp);
346 
358 HEADER_PREFIX void
359 addmul_hmatrix(field alpha, bool xtrans, pchmatrix x, bool ytrans, pchmatrix y,
360  pctruncmode tm, real eps, phmatrix z);
361 
374 HEADER_PREFIX void
375 addmul_lower_hmatrix(field alpha, bool xtrans, pchmatrix x, bool ytrans,
376  pchmatrix y, pctruncmode tm, real eps, phmatrix z);
377 
378 /* ------------------------------------------------------------
379  * Matrix inversion
380  * ------------------------------------------------------------ */
381 
391 HEADER_PREFIX void
393 
394 /* ------------------------------------------------------------
395  * Triangular matrices
396  * ------------------------------------------------------------ */
397 
406 clone_lower_hmatrix(bool aunit, pchmatrix a);
407 
416 clone_upper_hmatrix(bool aunit, pchmatrix a);
417 
418 /* ------------------------------------------------------------
419  * Solve triangular systems
420  * ------------------------------------------------------------ */
421 
434 HEADER_PREFIX void
435 triangularinvmul_hmatrix_avector(bool alower, bool aunit, bool atrans,
436  pchmatrix a, pavector xp);
437 
448 HEADER_PREFIX void
449 triangularsolve_hmatrix_avector(bool alower, bool aunit, bool atrans,
450  pchmatrix a, pavector x);
451 
465 HEADER_PREFIX void
466 triangularinvmul_hmatrix_amatrix(bool alower, bool aunit, bool atrans,
467  pchmatrix a, bool xtrans, pamatrix xp);
468 
485 HEADER_PREFIX void
486 triangularinvmul_hmatrix(bool alower, bool aunit, bool atrans, pchmatrix a,
487  pctruncmode tm, real eps, bool xtrans, phmatrix xp);
488 
505 HEADER_PREFIX void
506 triangularinvmul_amatrix_hmatrix(bool alower, bool aunit, bool atrans,
507  pcamatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix xp);
508 
509 /* ------------------------------------------------------------
510  * Evaluate triangular matrices
511  * ------------------------------------------------------------ */
512 
525 HEADER_PREFIX void
526 triangularmul_hmatrix_avector(bool alower, bool aunit, bool atrans, pchmatrix a,
527  pavector xp);
528 
539 HEADER_PREFIX void
540 triangulareval_hmatrix_avector(bool alower, bool aunit, bool atrans,
541  pchmatrix a, pavector x);
542 
556 HEADER_PREFIX void
557 triangularmul_hmatrix_amatrix(bool alower, bool aunit, bool atrans, pchmatrix a,
558  bool xtrans, pamatrix xp);
559 
576 HEADER_PREFIX void
577 triangularmul_hmatrix(bool alower, bool aunit, bool atrans, pchmatrix a,
578  pctruncmode tm, real eps, bool xtrans, phmatrix xp);
579 
580 /* ------------------------------------------------------------
581  * Triangular factorizations
582  * ------------------------------------------------------------ */
583 
598 HEADER_PREFIX void
600 
608 HEADER_PREFIX void
610 
618 HEADER_PREFIX void
620 
629 HEADER_PREFIX void
630 lrsolve_hmatrix_avector(bool atrans, pchmatrix a, pavector x);
631 
640 HEADER_PREFIX void
642 
651 HEADER_PREFIX void
653 
663 HEADER_PREFIX void
664 lreval_hmatrix_avector(bool atrans, pchmatrix a, pavector x);
665 
679 HEADER_PREFIX void
681 
689 HEADER_PREFIX void
691 
700 HEADER_PREFIX void
702 
705 #endif
void choleval_hmatrix_avector(pchmatrix a, pavector x)
Evaluate the linear system using the Cholesky factorization provided by choldecomp_hmatrix.
void triangularmul_hmatrix_avector(bool alower, bool aunit, bool atrans, pchmatrix a, pavector xp)
Evaluate a triangular system or .
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.
phmatrix split_sub_rkmatrix(pcrkmatrix r, pccluster rc, pccluster cc, bool rsplit, bool csplit)
Split an rkmatrix into submatrices referencing the original matrix.
void choldecomp_hmatrix(phmatrix a, pctruncmode tm, real eps)
Compute the Cholesky factorization, .
Definition: avector.h:39
void triangularinvmul_hmatrix_amatrix(bool alower, bool aunit, bool atrans, pchmatrix a, bool xtrans, pamatrix xp)
Solve a triangular system or .
void lreval_n_hmatrix_avector(pchmatrix a, pavector x)
Evaluate the linear system using the LR factorization provided by lrdecomp_hmatrix.
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 pa...
void lrsolve_t_hmatrix_avector(pchmatrix a, pavector x)
Solve the linear systems using the LR factorization provided by lrdecomp_hmatrix.
void add_hmatrix_amatrix(field alpha, bool atrans, pchmatrix a, pamatrix b)
Addition of a hierarchical matrix to a matrix, .
void cholsolve_hmatrix_avector(pchmatrix a, pavector x)
Solve the linear system using the Cholesky factorization provided by choldecomp_hmatrix.
void triangularinvmul_hmatrix_avector(bool alower, bool aunit, bool atrans, pchmatrix a, pavector xp)
Solve a triangular system, or .
void add_rkmatrix_amatrix(field alpha, bool atrans, pcrkmatrix a, pamatrix b)
Truncated addition of a full matrix to a low-rank matrix.
Representation of -matrices.
Definition: hmatrix.h:49
void lreval_t_hmatrix_avector(pchmatrix a, pavector x)
Evaluate the linear system using the LR factorization provided by lrdecomp_hmatrix.
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, .
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 hi...
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, .
void lrsolve_hmatrix_avector(bool atrans, pchmatrix a, pavector x)
Solve the linear systems or using the LR factorization provided by lrdecomp_hmatrix.
void triangularsolve_hmatrix_avector(bool alower, bool aunit, bool atrans, pchmatrix a, pavector x)
Solve a triangular system or .
phmatrix clone_upper_hmatrix(bool aunit, pchmatrix a)
Create a duplicate of the upper triangular part of an H-matrix.
void invert_hmatrix(phmatrix a, phmatrix work, pctruncmode tm, real eps)
Invert an H-matrix, .
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, ...
double _Complex field
Field type.
Definition: settings.h:171
void lrdecomp_hmatrix(phmatrix a, pctruncmode tm, real eps)
Compute the LR factorization, .
void triangularmul_hmatrix(bool alower, bool aunit, bool atrans, pchmatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix xp)
Evaluate a triangular system, or .
void lreval_hmatrix_avector(bool atrans, pchmatrix a, pavector x)
Evaluate the linear system or using the LR factorization provided by lrdecomp_hmatrix.
phmatrix split_rkmatrix(pcrkmatrix r, pccluster rc, pccluster cc, bool rsplit, bool csplit, bool copy)
Split an rkmatrix into submatrices.
void add_rkmatrix(field alpha, pcrkmatrix src, pctruncmode tm, real eps, prkmatrix trg)
Truncated addition of two low-rank matrices in rkmatrix representation, .
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, .
phmatrix split_sub_amatrix(pcamatrix f, pccluster rc, pccluster cc, bool rsplit, bool csplit)
Split an amatrix into submatrices.
prkmatrix merge_hmatrix_rkmatrix(pchmatrix s, pctruncmode tm, real eps)
Merge a depth-one hmatrix containing only rkmatrix submatrices into a new rkmatrix, .
Define different strategies used by various truncation and compression algorithms for hmatrices and h...
Definition: truncation.h:43
void addmul_hmatrix(field alpha, bool xtrans, pchmatrix x, bool ytrans, pchmatrix y, pctruncmode tm, real eps, phmatrix z)
Multiply two H-matrices, .
#define HEADER_PREFIX
Prefix for function declarations.
Definition: settings.h:43
double real
real floating point type.
Definition: settings.h:97
void triangularmul_hmatrix_amatrix(bool alower, bool aunit, bool atrans, pchmatrix a, bool xtrans, pamatrix xp)
Evaluate a triangular system or .
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, .
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 i...
void triangulareval_hmatrix_avector(bool alower, bool aunit, bool atrans, pchmatrix a, pavector x)
Evaluate a triangular system or .
void merge_rkmatrix(bool colmerge, pcrkmatrix src, pctruncmode tm, real eps, prkmatrix trg)
Split a leaf hmatrix into submatrices.
Representation of cluster trees.
Definition: cluster.h:40
Representation of a low-rank matrix in factorized form .
Definition: rkmatrix.h:36
phmatrix clone_lower_hmatrix(bool aunit, pchmatrix a)
Create a duplicate of the lower triangular part of an H-matrix.
Representation of a matrix as an array in column-major order.
Definition: amatrix.h:43
void trunc_rkmatrix(pctruncmode tm, real eps, prkmatrix r)
Truncate an rkmatrix, .
void triangularinvmul_amatrix_hmatrix(bool alower, bool aunit, bool atrans, pcamatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix xp)
Solve a triangular system, or .
void triangularinvmul_hmatrix(bool alower, bool aunit, bool atrans, pchmatrix a, pctruncmode tm, real eps, bool xtrans, phmatrix xp)
Solve a triangular system, or .
void lrsolve_n_hmatrix_avector(pchmatrix a, pavector x)
Solve the linear systems using the LR factorization provided by lrdecomp_hmatrix.