FORTRAN 90 Library Interface

New in version 3.0.

The next enumerations and functions are declared in primme_f90.inc.

subroutine primme_svds_matvec(x, ldx, y, ldy, blockSize, mode, primme_svds, ierr)

Abstract interface for the callbacks matrixMatvec and applyPreconditioner.

Parameters
  • x (ldx,*) [type(*),in] :: matrix with blockSize columns in column-major order with leading dimension ldx.

  • ldx [c_int64_t] :: the leading dimension of the array x.

  • y (ldy,*) [type(*),out] :: matrix with blockSize columns in column-major order with leading dimension ldy.

  • ldy [c_int64_t] :: the leading dimension of the array y.

  • blockSize [c_int,in] :: number of columns in x and y.

  • mode [c_int,in] :: a flag.

  • primme_svds [c_ptr,in] :: parameters structure created by primme_svds_params_create().

  • ierr [c_int,out] :: output error code; if it is set to non-zero, the current call to PRIMME will stop.

See more details about the precision and type and dimension for x and y, and the meaning of mode in the documentation of the callbacks.

primme_svds_params_create

function primme_svds_params_create()

Allocate and initialize a parameters structure to the default values.

After calling xprimme_svds() (or a variant), call primme_svds_params_destroy() to release allocated resources by PRIMME.

Return

primme_svds_params_create [c_ptr] :: pointer to a parameters structure.

primme_svds_set_method

function primme_svds_set_method(method, methodStage1, methodStage2, primme_svds)

Set PRIMME SVDS parameters to one of the preset configurations.

Parameters

xprimme_svds

function xprimme_svds(svals, svecs, resNorms, primme_svds)

Solve a real or complex singular value problem.

All arrays should be hosted on CPU. The computations are performed on CPU (see magma_xprimme_svds() for using GPUs).

Parameters
  • svals (*) [out] :: 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 (*) [out] :: 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 [c_ptr,in] :: parameters structure created by primme_params_create_svds().

Return

xprimme_svds [c_int] :: error indicator; see Error Codes.

The arrays svals, svecs, and resNorms should have the same kind.

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 - 1) * mLocal ). The i-th right singular vector starts at svecs(( numOrthoConst + initSize )* mLocal + ( numOrthoConst + i - 1)* nLocal ). The first vector has i=1.

All internal operations are performed at the same precision than svecs unless the user sets internalPrecision otherwise.

The type and precision of the callbacks depends on the type and precision of svecs. See details for matrixMatvec, applyPreconditioner, globalSumReal, broadcastReal, and convTestFun.

magma_xprimme_svds

function magma_xprimme_svds(svals, svecs, resNorms, primme_svds)

Solve a real or complex singular value problem.

Most of the computations are performed on GPU (see xprimme_svds() for using only the CPU).

Parameters
  • svals (*) [out] :: 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 (*) [out] :: 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 [c_ptr,in] :: parameters structure created by primme_params_create_svds().

Return

magma_xprimme_svds [c_int] :: error indicator; see Error Codes.

The arrays svals, svecs, and resNorms should have the same kind.

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 - 1) * mLocal ). The i-th right singular vector starts at svecs(( numOrthoConst + initSize )* mLocal + ( numOrthoConst + i - 1)* nLocal ). The first vector has i=1.

All internal operations are performed at the same precision than svecs unless the user sets internalPrecision otherwise.

The type and precision of the callbacks depends on the type and precision of svecs. See details for matrixMatvec, applyPreconditioner, globalSumReal, broadcastReal, and convTestFun.

primme_svds_params_destroy

function primme_svds_params_destroy(primme_svds)

Free memory allocated by PRIMME associated to a parameters structure created with primme_svds_params_create().

Parameters

primme_svds [c_ptr] :: parameters structure.

Return

primme_svds_params_destroy :: nonzero value if the call is not successful.

primme_svds_set_member

function primme_svds_set_member(primme_svds, label, value)

Set a value in some field of the parameter structure.

Parameters
Return

primme_svds_set_member [c_int] ::

nonzero value if the call is not successful.

Examples:

type(c_ptr) :: primme_svds
integer :: ierr
...

integer(c_int64_t) :: m               = 100
ierr = primme_svds_set_member(primme_svds, PRIMME_SVDS_m, m)
ierr = primme_svds_set_member(primme_svds, PRIMME_SVDS_n, m)

real(c_double) :: tol             = 1.0D-12
ierr = primme_svds_set_member(primme, PRIMME_SVDS_eps, tol)

integer(c_int64_t), parameter :: numTargetShifts = 2
real(c_double) :: TargetShifts(numTargetShifts) = (/3.0D0, 5.1D0/)
ierr = primme_svds_set_member(primme_svds, PRIMME_SVDS_numTargetShifts, numTargetShifts)
ierr = primme_svds_set_member(primme_svds, PRIMME_SVDS_targetShifts, TargetShifts)

ierr = primme_svds_set_member(primme_svds, PRIMME_SVDS_target, primme_svds_closest_abs)

procedure(primme_svds_matvec) :: MV, ApplyPrecon
ierr = primme_svds_set_member(primme_svds, PRIMME_SVDS_matrixMatvec, MV)

ierr = primme_svds_set_member(primme_svds, PRIMME_SVDS_applyPreconditioner,
                         c_funloc(ApplyPrecon))

type(c_ptr) :: primme
ierr = primme_svds_get_member(primme_svds, PRIMME_SVDS_primme, primme)
ierr = primme_set_member(primme, PRIMME_correctionParams_precondition,
                         1_c_int64_t)

primme_get_member

function primme_svds_get_member(primme, label, value)

Get the value in some field of the parameter structure.

Parameters
Return

primme_svds_get_member [c_int] :: nonzero value if the call is not successful.

Examples:

type(c_ptr) :: primme_svds
integer :: ierr
...

integer(c_int64_t) :: m
ierr = primme_svds_get_member(primme_svds, PRIMME_SVDS_m, m)

real(c_double) :: aNorm
ierr = primme_svds_get_member(primme_svds, PRIMME_SVDS_aNorm, aNorm)

type(c_ptr) :: primme
ierr = primme_svds_get_member(primme_svds, PRIMME_SVDS_primme, primme)
ierr = primme_set_member(primme, PRIMME_correctionParams_precondition,
                         1_c_int64_t)