heat.linalg.svdtools

distributed hierarchical SVD

Module Contents

hsvd_rank(A: heat.core.dndarray.DNDarray, maxrank: int, compute_sv: bool = False, maxmergedim: int | None = None, safetyshift: int = 5, silent: bool = True) Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, float] | Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray] | heat.core.dndarray.DNDarray

Hierarchical SVD (hSVD) with prescribed truncation rank maxrank. If A = U diag(sigma) V^T is the true SVD of A, this routine computes an approximation for U[:,:maxrank] (and sigma[:maxrank], V[:,:maxrank]).

The accuracy of this approximation depends on the structure of A (“low-rank” is best) and appropriate choice of parameters.

One can expect a similar outcome from this routine as for sci-kit learn’s TruncatedSVD (with algorithm=’randomized’) although a different, determinstic algorithm is applied here. Hereby, the parameters n_components and n_oversamples (sci-kit learn) roughly correspond to maxrank and safetyshift (see below).

ADNDarray

2D-array (float32/64) of which the hSVD has to be computed.

maxrankint

truncation rank. (This parameter corresponds to n_components in sci-kit learn’s TruncatedSVD.)

compute_svbool, optional

compute_sv=True implies that also Sigma and V are computed and returned. The default is False.

maxmergedimint, optional

maximal size of the concatenation matrices during the merging procedure. The default is None and results in an appropriate choice depending on the size of the local slices of A and maxrank. Too small choices for this parameter will result in failure if the maximal size of the concatenation matrices does not allow to merge at least two matrices. Too large choices for this parameter can cause memory errors if the resulting merging problem becomes too large.

safetyshiftint, optional

Increases the actual truncation rank within the computations by a safety shift. The default is 5. (There is some similarity to n_oversamples in sci-kit learn’s TruncatedSVD.)

silentbool, optional

silent=False implies that some information on the computations are printed. The default is True.

(Union[ Tuple[DNDarray, DNDarray, DNDarray, float], Tuple[DNDarray, DNDarray, DNDarray], DNDarray])

if compute_sv=True: U, Sigma, V, a-posteriori error estimate for the reconstruction error ||A-U Sigma V^T ||_F / ||A||_F (computed according to [2] along the “true” merging tree). if compute_sv=False: U, a-posteriori error estimate

The size of the process local SVDs to be computed during merging is proportional to the non-split size of the input A and (maxrank + safetyshift). Therefore, conservative choice of maxrank and safetyshift is advised to avoid memory issues. Note that, as sci-kit learn’s randomized SVD, this routine is different from numpy.linalg.svd because not all singular values and vectors are computed and even those computed may be inaccurate if the input matrix exhibts a unfavorable structure.

See also

hsvd()

hsvd_rtol()

References ——- [1] Iwen, Ong. A distributed and incremental SVD algorithm for agglomerative data analysis on large networks. SIAM J. Matrix Anal. Appl., 37(4), 2016. [2] Himpe, Leibner, Rave. Hierarchical approximate proper orthogonal decomposition. SIAM J. Sci. Comput., 40 (5), 2018.

hsvd_rtol(A: heat.core.dndarray.DNDarray, rtol: float, compute_sv: bool = False, maxrank: int | None = None, maxmergedim: int | None = None, safetyshift: int = 5, no_of_merges: int | None = None, silent: bool = True) Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, float] | Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray] | heat.core.dndarray.DNDarray

Hierchical SVD (hSVD) with prescribed upper bound on the relative reconstruction error. If A = U diag(sigma) V^T is the true SVD of A, this routine computes an approximation for U[:,:r] (and sigma[:r], V[:,:r]) such that the rel. reconstruction error ||A-U[:,:r] diag(sigma[:r]) V[:,:r]^T ||_F / ||A||_F does not exceed rtol.

The accuracy of this approximation depends on the structure of A (“low-rank” is best) and appropriate choice of parameters. This routine is similar to hsvd_rank with the difference that truncation is not performed after a fixed number (namly maxrank many) singular values but after such a number of singular values that suffice to capture a prescribed fraction of the amount of information contained in the input data (rtol).

ADNDarray

2D-array (float32/64) of which the hSVD has to be computed.

rtolfloat

desired upper bound on the relative reconstruction error ||A-U Sigma V^T ||_F / ||A||_F. This upper bound is processed into ‘local’ tolerances during the actual computations assuming the worst case scenario of a binary “merging tree”; therefore, the a-posteriori error for the relative error using the true “merging tree” (see output) may be significantly smaller than rtol. Prescription of maxrank or maxmergedim (disabled in default) can result in loss of desired precision, but can help to avoid memory issues.

compute_svbool, optional

compute_sv=True implies that also Sigma and V are computed and returned. The default is False.

no_of_mergesint, optional

Maximum number of processes to be merged at each step. If no further arguments are provided (see below), this completely determines the “merging tree” and may cause memory issues. The default is None and results in a binary merging tree. Note that no_of_merges dominates maxrank and maxmergedim in the sense that at most no_of_merges processes are merged even if maxrank and maxmergedim would allow merging more processes.

maxrankint, optional

maximal truncation rank. The default is None. Setting at least one of maxrank and maxmergedim is recommended to avoid memory issues, but can result in loss of desired precision. Setting only maxrank (and not maxmergedim) results in an appropriate default choice for maxmergedim depending on the size of the local slices of A and the value of maxrank.

maxmergedimint, optional

maximal size of the concatenation matrices during the merging procedure. The default is None and results in an appropriate choice depending on the size of the local slices of A and maxrank. The default is None. Too small choices for this parameter will result in failure if the maximal size of the concatenation matrices does not allow to merge at least two matrices. Too large choices for this parameter can cause memory errors if the resulting merging problem becomes too large. Setting at least one of maxrank and maxmergedim is recommended to avoid memory issues, but can result in loss of desired precision. Setting only maxmergedim (and not maxrank) results in an appropriate default choice for maxrank.

safetyshiftint, optional

Increases the actual truncation rank within the computations by a safety shift. The default is 5.

silentbool, optional

silent=False implies that some information on the computations are printed. The default is True.

(Union[ Tuple[DNDarray, DNDarray, DNDarray, float], Tuple[DNDarray, DNDarray, DNDarray], DNDarray])

if compute_sv=True: U, Sigma, V, a-posteriori error estimate for the reconstruction error ||A-U Sigma V^T ||_F / ||A||_F (computed according to [2] along the “true” merging tree used in the computations). if compute_sv=False: U, a-posteriori error estimate

The maximum size of the process local SVDs to be computed during merging is proportional to the non-split size of the input A and (maxrank + safetyshift). Therefore, conservative choice of maxrank and safetyshift is advised to avoid memory issues. For similar reasons, prescribing only rtol and the number of processes to be merged in each step (without specifying maxrank or maxmergedim) may result in memory issues. Although prescribing maxrank is therefore strongly recommended to avoid memory issues, but may result in loss of desired precision (rtol). If this occures, a separate warning will be raised.

Note that this routine is different from numpy.linalg.svd because not all singular values and vectors are computed and even those computed may be inaccurate if the input matrix exhibts a unfavorable structure.

To avoid confusion, note that rtol in this routine does not have any similarity to tol in scikit learn’s TruncatedSVD.

See also

hsvd()

hsvd_rank()

References ——- [1] Iwen, Ong. A distributed and incremental SVD algorithm for agglomerative data analysis on large networks. SIAM J. Matrix Anal. Appl., 37(4), 2016. [2] Himpe, Leibner, Rave. Hierarchical approximate proper orthogonal decomposition. SIAM J. Sci. Comput., 40 (5), 2018.

hsvd(A: heat.core.dndarray.DNDarray, maxrank: int | None = None, maxmergedim: int | None = None, rtol: float | None = None, safetyshift: int = 0, no_of_merges: int | None = 2, compute_sv: bool = False, silent: bool = True, warnings_off: bool = False) Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, float] | Tuple[heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray, heat.core.dndarray.DNDarray] | heat.core.dndarray.DNDarray

This function computes an approximate truncated SVD of A utilizing a distributed hiearchical algorithm; see the references. The present function hsvd is a low-level routine, provides many options/parameters, but no default values, and is not recommended for usage by non-experts since conflicts arising from inappropriate parameter choice will not be catched. We strongly recommend to use the corresponding high-level functions hsvd_rank and hsvd_rtol instead.

Input

A: DNDarray

2D-array (float32/64) of which hSVD has to be computed

maxrank: int, optional

truncation rank of the SVD

maxmergedim: int, optional

maximal size of the concatenation matrices when “merging” the local SVDs

rtol: float, optional

upper bound on the relative reconstruction error ||A-U Sigma V^T ||_F / ||A||_F (may deteriorate due to other parameters)

safetyshift: int, optional

shift that increases the actual truncation rank of the local SVDs during the computations in order to increase accuracy

no_of_merges: int, optional

maximum number of local SVDs to be “merged” at one step

compute_sv: bool, optional

determines whether to compute U, Sigma, V (compute_sv=True) or not (then U only)

silent: bool, optional

determines whether to print infos on the computations performed (silent=False)

warnings_off: bool, optional

switch on and off warnings that are not intended for the high-level routines based on this function

returns:

if compute_sv=True: U, Sigma, V, a-posteriori error estimate for the reconstruction error ||A-U Sigma V^T ||_F / ||A||_F (computed according to [2] along the “true” merging tree used in the computations). if compute_sv=False: U, a-posteriori error estimate

rtype:

(Union[ Tuple[DNDarray, DNDarray, DNDarray, float], Tuple[DNDarray, DNDarray, DNDarray], DNDarray])

References

[1] Iwen, Ong. A distributed and incremental SVD algorithm for agglomerative data analysis on large networks. SIAM J. Matrix Anal. Appl., 37(4), 2016. [2] Himpe, Leibner, Rave. Hierarchical approximate proper orthogonal decomposition. SIAM J. Sci. Comput., 40 (5), 2018.