C Library Interface¶
New in version 2.0.
The PRIMME SVDS interface is composed of the following functions. To solve real and complex singular value problems call respectively:
intdprimme_svds(double *svals, double *svecs, double *resNorms, primme_svds_params *primme_svds) intzprimme_svds(double *svals,PRIMME_COMPLEX_DOUBLE*svecs, double *resNorms, primme_svds_params *primme_svds)
There are more versions for matrix problems working in other precisions:
Precision |
Real |
Complex |
|---|---|---|
half |
||
single |
||
double |
Other useful functions:
voidprimme_svds_initialize(primme_svds_params *primme_svds) intprimme_svds_set_method(primme_svds_preset_method method, primme_preset_method methodStage1, primme_preset_method methodStage2, primme_svds_params *primme_svds) voidprimme_svds_display_params(primme_svds_params primme_svds) voidprimme_svds_free(primme_svds_params *primme_svds)
PRIMME SVDS stores its data on the structure primme_svds_params.
See Parameters Guide for an introduction about its fields.
Running¶
To use PRIMME SVDS, follow these basic steps.
Include:
#include "primme.h" /* header file is required to run primme */
Initialize a PRIMME SVDS parameters structure for default settings:
primme_svds_paramsprimme_svds;primme_svds_initialize(&primme_svds);Set problem parameters (see also Parameters Guide), and, optionally, set one of the
preset methods:primme_svds.
matrixMatvec= matrixMatvec; /* MV product */ primme_svds.m= 1000; /* set the matrix dimensions */ primme_svds.n= 100; primme_svds.numSvals= 10; /* Number of singular values */primmesvds_set_method(primme_svds_hybrid, PRIMME_DEFAULT_METHOD, PRIMME_DEFAULT_METHOD, &primme_svds); ...Then to solve a real singular value problem call:
ret =
dprimme_svds(svals, svecs, resNorms, &primme_svds);The previous is the double precision call. There is available calls for complex double, single and complex single; check
zprimme_svds(),sprimme_svds()andcprimme_svds().To solve complex singular value problems call:
ret =
zprimme_svds(svals, svecs, resNorms, &primme_svds);The call arguments are:
svals, array to return the found singular values;
svecs, array to return the found left and right singular vectors;
resNorms, array to return the residual norms of the found triplets; and
ret, returned error code.
To free the work arrays in PRIMME SVDS:
primme_svds_free(&primme_svds);
Parameters Guide¶
PRIMME SVDS stores the data on the structure primme_svds_params, which has the next fields:
PRIMME_INT m, number of rows of the matrix.PRIMME_INT n, number of columns of the matrix.int numSvals, how many singular triplets to find.primme_svds_target target, which singular values to find.double eps, tolerance of the residual norm of converged triplets.int numProcs, number of processesint procID, rank of this processPRIMME_INT mLocal, number of rows stored in this processPRIMME_INT nLocal, number of columns stored in this processint initSize, initial vectors as approximate solutions.int maxBasisSizeint minRestartSizeint maxBlockSizevoid * commInfovoid * matrixvoid * preconditionervoid * convtestvoid * monitorint numTargetShifts, for targeting interior singular values.double * targetShiftsint numOrthoConst, orthogonal constrains to the singular vectors.int lockingPRIMME_INT maxMatvecsdouble aNormint printLevelFILE * outputFileprimme_svds_operator methodprimme_svds_operator methodStage2primme_op_datatype matrixMatvec_typeprimme_op_datatype applyPreconditioner_typeprimme_op_datatype globalSumReal_typeprimme_op_datatype broadcastReal_typeprimme_op_datatype internalPrecisionPRIMME SVDS requires the user to set at least the matrix dimensions (m x n) and
the matrix-vector product (matrixMatvec), as they define the problem to be solved.
For parallel programs, mLocal, nLocal, procID and globalSumReal are also required.
In addition, most users would want to specify how many singular triplets to find, and provide a preconditioner (if available).
It is useful to have set all these before calling primme_svds_set_method().
Also, if users have a preference on maxBasisSize, maxBlockSize, etc, they
should also provide them into primme_svds_params prior to the
primme_svds_set_method() call. This helps primme_svds_set_method() make
the right choice on other parameters. It is sometimes useful to check the actual
parameters that PRIMME SVDS is going to use (before calling it) or used (on return)
by printing them with primme_svds_display_params().
Interface Description¶
The next enumerations and functions are declared in primme.h.
?primme_svds¶
-
int
hprimme_svds(PRIMME_HALF *svals, PRIMME_HALF *svecs, PRIMME_HALF *resNorms, primme_svds_params *primme_svds)¶
-
int
hsprimme_svds(float *svals, PRIMME_HALF *svecs, float *resNorms, primme_svds_params *primme_svds)¶
-
int
kprimme_svds(PRIMME_HALF *svals, PRIMME_COMPLEX_HALF *svecs, PRIMME_HALF *resNorms, primme_svds_params *primme_svds)¶
-
int
ksprimme_svds(float *svals, PRIMME_COMPLEX_HALF *svecs, float *resNorms, primme_svds_params *primme_svds)¶ New in version 3.0.
-
int
sprimme_svds(float *svals, float *svecs, float *resNorms, primme_svds_params *primme_svds)¶
-
int
cprimme_svds(float *svals, PRIMME_COMPLEX_FLOAT *svecs, float *resNorms, primme_svds_params *primme_svds)¶ New in version 2.0.
-
int
dprimme_svds(double *svals, double *svecs, double *resNorms, primme_svds_params *primme_svds)¶
-
int
zprimme_svds(double *svals, PRIMME_COMPLEX_DOUBLE *svecs, double *resNorms, primme_svds_params *primme_svds)¶ Solve a real singular value problem.
All arrays should be hosted on CPU. The computations are performed on CPU (see
magma_dprimme_svds()for using GPUs).- Parameters
svals – array at least of size
numSvalsto store the computed singular values; all processes in a parallel run return this local array with the same values.svecs – array at least of size (
mLocal+nLocal) times (numOrthoConst+numSvals) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors.resNorms – array at least of size
numSvalsto store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values.primme_svds – parameters structure.
- Returns
error indicator; see Error Codes.
On input,
svecsshould start with the content of thenumOrthoConstleft vectors, followed by theinitSizeleft vectors, followed by thenumOrthoConstright vectors, and followed by theinitSizeright vectors.On return, the i-th left singular vector starts at svecs[(
numOrthoConst+i)*mLocal]. The i-th right singular vector starts at svecs[(numOrthoConst+initSize)*mLocal+ (numOrthoConst+i)*nLocal]. The first vector has i=0.All internal operations are performed at the same precision than
svecsunless the user setsinternalPrecisionotherwise. The functionshsprimme_svds()andksprimme_svds()perform all computations in half precision by default and report the eigenvalues and the residual norms in single precision. These functions may help in applications that may be not built with a compiler supporting half precision.The type and precision of the callbacks depends on the type and precision of
svecs. Although this can be changed. See details formatrixMatvec,applyPreconditioner,globalSumReal,broadcastReal, andconvTestFun.
magma_?primme_svds¶
-
int
magma_hprimme_svds(PRIMME_HALF *svals, PRIMME_HALF *svecs, PRIMME_HALF *resNorms, primme_svds_params *primme_svds)¶
-
int
magma_hsprimme_svds(float *svals, PRIMME_HALF *svecs, float *resNorms, primme_svds_params *primme_svds)¶
-
int
magma_kprimme_svds(PRIMME_HALF *svals, PRIMME_COMPLEX_HALF *svecs, PRIMME_HALF *resNorms, primme_svds_params *primme_svds)¶
-
int
magma_ksprimme_svds(float *svals, PRIMME_COMPLEX_HALF *svecs, float *resNorms, primme_svds_params *primme_svds)¶
-
int
magma_sprimme_svds(float *svals, float *svecs, float *resNorms, primme_svds_params *primme_svds)¶
-
int
magma_cprimme_svds(float *svals, PRIMME_COMPLEX_FLOAT *svecs, float *resNorms, primme_svds_params *primme_svds)¶
-
int
magma_dprimme_svds(double *svals, double *svecs, double *resNorms, primme_svds_params *primme_svds)¶
-
int
magma_zprimme_svds(double *svals, PRIMME_COMPLEX_DOUBLE *svecs, double *resNorms, primme_svds_params *primme_svds)¶ Solve a real singular value problem.
Most of the computations are performed on GPU (see
dprimme_svds()for using only the CPU).- Parameters
svals – CPU array at least of size
numSvalsto store the computed singular values; all processes in a parallel run return this local array with the same values.svecs – GPU array at least of size (
mLocal+nLocal) times (numOrthoConst+numSvals) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors.resNorms – CPU array at least of size
numSvalsto store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values.primme_svds – parameters structure.
- Returns
error indicator; see Error Codes.
On input,
svecsshould start with the content of thenumOrthoConstleft vectors, followed by theinitSizeleft vectors, followed by thenumOrthoConstright vectors, and followed by theinitSizeright vectors.On return, the i-th left singular vector starts at svecs[(
numOrthoConst+i)*mLocal]. The i-th right singular vector starts at svecs[(numOrthoConst+initSize)*mLocal+ (numOrthoConst+i)*nLocal]. The first vector has i=0.All internal operations are performed at the same precision than
svecsunless the user setsinternalPrecisionotherwise. The functionsmagma_hsprimme_svds()andmagma_ksprimme_svds()perform all computations in half precision by default and report the eigenvalues and the residual norms in single precision. These functions may help in applications that may be not built with a compiler supporting half precision.The type and precision of the callbacks depends on the type and precision of
svecs. Although this can be changed. See details formatrixMatvec,applyPreconditioner,globalSumReal,broadcastReal, andconvTestFun.New in version 3.0.
primme_svds_initialize¶
-
void
primme_svds_initialize(primme_svds_params *primme_svds)¶ Initialize PRIMME SVDS parameters structure to the default values.
After calling
dprimme_svds()(or a variant), callprimme_svds_free()to release allocated resources by PRIMME.- Parameters
primme_svds – parameters structure.
Example:
primme_svds_params primme_svds; primme_svds_initialize(&primme_svds); primme_svds.n = 100; ... dprimme_svds(svals, svecs, rnorms, &primme_svds); ... primme_svds_free(&primme_svds);
See the alternative function
primme_svds_params_create()that also allocates the structure.
primme_svds_create¶
-
primme_svds_params *
primme_svds_params_create(void)¶ Allocate and initialize a parameters structure to the default values.
After calling
dprimme_svds()(or a variant), callprimme_svds_params_destroy()to release allocated resources by PRIMME.- Parameters
primme_sv – parameters structure.
Example:
primme_svds_params *primme_svds = primme_svds_params_create(); primme_svds->n = 100; ... dprimme_svds(svals, svecs, rnorms, primme_svds); ... primme_svds_params_destroy(primme_svds);
See the alternative function
primme_svds_initialize()that only initializes the structure.New in version 3.0.
primme_svds_set_method¶
-
int
primme_svds_set_method(primme_svds_preset_method method, primme_preset_method methodStage1, primme_preset_method methodStage2, primme_svds_params *primme_svds)¶ Set PRIMME SVDS parameters to one of the preset configurations.
- Parameters
method –
preset method to compute the singular triplets; one of
primme_svds_default, currently set asprimme_svds_hybrid.primme_svds_normalequations, compute the eigenvectors of \(A^*A\) or \(A A^*\).primme_svds_augmented, compute the eigenvectors of the augmented matrix, \(\left(\begin{array}{cc} 0 & A^* \\ A & 0 \end{array}\right)\).primme_svds_hybrid, start withprimme_svds_normalequations; use the resulting approximate singular vectors as initial vectors forprimme_svds_augmentedif the required accuracy was not achieved.
methodStage1 – preset method to compute the eigenpairs at the first stage; see available values at
primme_set_method().methodStage2 – preset method to compute the eigenpairs with the second stage of
primme_svds_hybrid; see available values atprimme_set_method().primme_svds – parameters structure.
See also Preset Methods.
primme_svds_display_params¶
-
void
primme_svds_display_params(primme_svds_params primme_svds)¶ Display all printable settings of
primme_svdsinto the file descriptoroutputFile.- Parameters
primme_svds – parameters structure.
primme_svds_free¶
-
void
primme_svds_free(primme_svds_params *primme_svds)¶ Free memory allocated by PRIMME SVDS.
- Parameters
primme_svds – parameters structure.
primme_svds_params_destroy¶
-
int
primme_svds_params_destroy(primme_svds_params *primme)¶ Free memory allocated by PRIMME associated to a parameters structure created with
primme_svds_params_create().- Parameters
primme_svds – parameters structure.
- Returns
nonzero value if the call is not successful.
New in version 3.0.