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:

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)

There are more versions for matrix problems working in other precisions:

Precision

Real

Complex

half

hprimme_svds() hsprimme_svds()

kprimme_svds() ksprimme_svds()

single

sprimme_svds()

cprimme_svds()

double

dprimme_svds()

zprimme_svds()

Other useful functions:

void primme_svds_initialize (primme_svds_params *primme_svds)
int primme_svds_set_method (primme_svds_preset_method method,
   primme_preset_method methodStage1,
   primme_preset_method methodStage2, primme_svds_params *primme_svds)
void primme_svds_display_params (primme_svds_params primme_svds)
void primme_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.

  1. Include:

    #include "primme.h"   /* header file is required to run primme */
    
  2. Initialize a PRIMME SVDS parameters structure for default settings:

    primme_svds_params primme_svds;
    primme_svds_initialize (&primme_svds);
  3. 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);
    ...
  4. 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() and cprimme_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.

  5. 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:

Basic
PRIMME_INT m, number of rows of the matrix.
PRIMME_INT n, number of columns of the matrix.
void (* matrixMatvec )(...), matrix-vector product.
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.

For parallel programs
int numProcs, number of processes
int procID, rank of this process
PRIMME_INT mLocal, number of rows stored in this process
PRIMME_INT nLocal, number of columns stored in this process
void (* globalSumReal )(...), sum reduction among processes

Accelerate the convergence
void (* applyPreconditioner )(...), preconditioner-vector product.
int initSize, initial vectors as approximate solutions.
int minRestartSize

User data
void * commInfo
void * matrix
void * convtest
void * monitor

Advanced options
int numTargetShifts, for targeting interior singular values.
double * targetShifts
int numOrthoConst, orthogonal constrains to the singular vectors.
PRIMME_INT maxMatvecs
PRIMME_INT iseed [4]
double aNorm
FILE * outputFile
primme_svds_operator method
primme_svds_operator methodStage2
void (* convTestFun )(...), custom convergence criterion.
void (* monitorFun )(...), custom convergence history.
primme_op_datatype matrixMatvec_type
primme_op_datatype applyPreconditioner_type
primme_op_datatype globalSumReal_type
primme_op_datatype broadcastReal_type
primme_op_datatype internalPrecision

PRIMME 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 numSvals to 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 numSvals to 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, svecs should start with the content of the numOrthoConst left vectors, followed by the initSize left vectors, followed by the numOrthoConst right vectors, and followed by the initSize right 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 svecs unless the user sets internalPrecision otherwise. The functions hsprimme_svds() and 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 for matrixMatvec, applyPreconditioner, globalSumReal, broadcastReal, and convTestFun.

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 numSvals to 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 numSvals to 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, svecs should start with the content of the numOrthoConst left vectors, followed by the initSize left vectors, followed by the numOrthoConst right vectors, and followed by the initSize right 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 svecs unless the user sets internalPrecision otherwise. The functions magma_hsprimme_svds() and magma_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 for matrixMatvec, applyPreconditioner, globalSumReal, broadcastReal, and convTestFun.

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), call primme_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), call primme_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

See also Preset Methods.

primme_svds_display_params

void primme_svds_display_params(primme_svds_params primme_svds)

Display all printable settings of primme_svds into the file descriptor outputFile.

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.