FORTRAN Library Interface¶
New in version 2.0.
The next enumerations and functions are declared in primme_svds_f77.h
.
sprimme_svds_f77¶
-
subroutine
sprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a real singular value problem using single precision.
All arrays should be hosted on CPU. The computations are performed on CPU (see
magma_sprimme_svds_f77()
for using GPUs).- Parameters
svals (*) [real] :: (output) 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 (*) [real] :: 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 (*) [real] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.
cprimme_svds_f77¶
-
subroutine
cprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a complex singular value problem using single precision.
All arrays should be hosted on CPU. The computations are performed on CPU (see
magma_cprimme_svds_f77()
for using GPUs).- Parameters
svals (*) [real] :: (output) 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 (*) [complex] :: 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 (*) [real] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.
dprimme_svds_f77¶
-
subroutine
dprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a real singular value problem using double precision.
All arrays should be hosted on CPU. The computations are performed on CPU (see
magma_dprimme_svds_f77()
for using GPUs).- Parameters
svals (*) [double precision] :: (output) 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 (*) [double precision] :: 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 (*) [double precision] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.
zprimme_svds_f77¶
-
subroutine
zprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a complex singular value problem using double precision.
All arrays should be hosted on CPU. The computations are performed on CPU (see
magma_zprimme_svds_f77()
for using GPUs).- Parameters
svals (*) [double precision] :: (output) 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 (*) [complex double precision] :: 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 (*) [double precision] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.
magma_sprimme_svds_f77¶
-
subroutine
magma_sprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a real singular value problem using single precision.
Most of the computations are performed on GPU (see
sprimme_svds_f77()
for using only the CPU).- Parameters
svals (*) [real] :: (output) 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 (*) [real] :: 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 (*) [real] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.New in version 3.0.
magma_cprimme_svds_f77¶
-
subroutine
magma_cprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a complex singular value problem using single precision.
Most of the computations are performed on GPU (see
cprimme_svds_f77()
for using only the CPU).- Parameters
svals (*) [real] :: (output) 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 (*) [complex] :: 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 (*) [real] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.New in version 3.0.
magma_dprimme_svds_f77¶
-
subroutine
magma_dprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a real singular value problem using double precision.
Most of the computations are performed on GPU (see
dprimme_svds_f77()
for using only the CPU).- Parameters
svals (*) [double precision] :: (output) 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 (*) [double precision] :: 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 (*) [double precision] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.New in version 3.0.
magma_zprimme_svds_f77¶
-
subroutine
magma_zprimme_svds_f77
(svals, svecs, resNorms, primme_svds, ierr)¶ Solve a complex singular value problem using double precision.
Most of the computations are performed on GPU (see
zprimme_svds_f77()
for using only the CPU).- Parameters
svals (*) [double precision] :: (output) 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 (*) [complex double precision] :: 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 (*) [double precision] :: 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 [ptr] :: parameters structure.
ierr [integer] :: (output) error indicator; see Error Codes.
On input,
svecs
should start with the content of thenumOrthoConst
left vectors, followed by theinitSize
left vectors, followed by thenumOrthoConst
right vectors and followed by theinitSize
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 setsinternalPrecision
otherwise.The type and precision of the callbacks depends on the type and precision of svecs. See details for
matrixMatvec
,applyPreconditioner
,globalSumReal
,broadcastReal
, andconvTestFun
.New in version 3.0.
primme_svds_initialize_f77¶
-
subroutine
primme_svds_initialize_f77
(primme_svds, ierr)¶ Set PRIMME SVDS parameters structure to the default values.
After calling
dprimme_svds_f77()
(or a variant), callprimme_svds_free_f77()
to release allocated resources by PRIMME.- Parameters
primme_svds [ptr] :: (output) parameters structure.
primme_svds_set_method_f77¶
-
subroutine
primme_svds_set_method_f77
(method, methodStage1, methodStage2, primme_svds, ierr)¶ Set PRIMME SVDS parameters to one of the preset configurations.
- Parameters
method [integer] ::
(input) preset configuration 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_augmented
if the required accuracy was not achieved.
methodStage1 [primme_preset_method] :: (input) preset method to compute the eigenpairs at the first stage; see available values at
primme_set_method_f77()
.methodStage2 [primme_preset_method] :: (input) preset method to compute the eigenpairs with the second stage of
PRIMME_SVDS_hybrid
; see available values atprimme_set_method_f77()
.primme_svds [ptr] :: (input/output) parameters structure.
ierr [integer] :: (output) if 0, successful; if negative, something went wrong.
primme_svds_display_params_f77¶
-
subroutine
primme_svds_display_params_f77
(primme_svds)¶ Display all printable settings of
primme_svds
into the file descriptoroutputFile
.- Parameters
primme_svds [ptr] :: (input) parameters structure.
primme_svds_free_f77¶
primme_svds_set_member_f77¶
-
subroutine
primme_svds_set_member_f77
(primme_svds, label, value, ierr)¶ Set a value in some field of the parameter structure.
- Parameters
primme_svds [ptr] :: (input) parameters structure.
label [integer] ::
field where to set value. One of:
PRIMME_SVDS_numOrthoConst
value :: (input) value to set.
Note
Don’t use this subroutine inside PRIMME SVDS’s callback functions, e.g.,
matrixMatvec
orapplyPreconditioner
, or in functions called by these functions.
primme_svdstop_get_member_f77¶
-
subroutine
primme_svdstop_get_member_f77
(primme_svds, label, value, ierr)¶ Get the value in some field of the parameter structure.
- Parameters
primme_svds [ptr] :: (input) parameters structure.
label [integer] :: (input) field where to get value. One of the detailed in subroutine
primmesvds_top_set_member_f77()
.value :: (output) value of the field.
Note
Don’t use this subroutine inside PRIMME SVDS’s callback functions, e.g.,
matrixMatvec
orapplyPreconditioner
, or in functions called by these functions. In those cases useprimme_svds_get_member_f77()
.Note
When
label
is one ofPRIMME_SVDS_matrixMatvec
,PRIMME_SVDS_applyPreconditioner
,PRIMME_SVDS_commInfo
,PRIMME_SVDS_intWork
,PRIMME_SVDS_realWork
,PRIMME_SVDS_matrix
andPRIMME_SVDS_preconditioner
, the returnedvalue
is a C pointer (void*
). Use Fortran pointer or other extensions to deal with it. For instance:use iso_c_binding MPI_Comm comm comm = MPI_COMM_WORLD call primme_svds_set_member_f77(primme_svds, PRIMME_SVDS_commInfo, comm) ... subroutine par_GlobalSumDouble(x,y,k,primme_svds) use iso_c_binding implicit none ... MPI_Comm, pointer :: comm type(c_ptr) :: pcomm call primme_svds_get_member_f77(primme_svds, PRIMME_SVDS_commInfo, pcomm) call c_f_pointer(pcomm, comm) call MPI_Allreduce(x,y,k,MPI_DOUBLE,MPI_SUM,comm,ierr)
Most users would not need to retrieve these pointers in their programs.
primme_svds_get_member_f77¶
-
subroutine
primme_svds_get_member_f77
(primme_svds, label, value, ierr)¶ Get the value in some field of the parameter structure.
- Parameters
primme_svds [ptr] :: (input) parameters structure.
label [integer] :: (input) field where to get value. One of the detailed in subroutine
primme_svdstop_set_member_f77()
.value :: (output) value of the field.
Note
Use this subroutine exclusively inside PRIMME SVDS’s callback functions, e.g.,
matrixMatvec
orapplyPreconditioner
, or in functions called by these functions. Otherwise, e.g., from the main program, use the subroutineprimme_svdstop_get_member_f77()
.Note
When
label
is one ofPRIMME_SVDS_matrixMatvec
,PRIMME_SVDS_applyPreconditioner
,PRIMME_SVDS_commInfo
,PRIMME_SVDS_intWork
,PRIMME_SVDS_realWork
,PRIMME_SVDS_matrix
andPRIMME_SVDS_preconditioner
, the returnedvalue
is a C pointer (void*
). Use Fortran pointer or other extensions to deal with it. For instance:use iso_c_binding MPI_Comm comm comm = MPI_COMM_WORLD call primme_svds_set_member_f77(primme_svds, PRIMME_SVDS_commInfo, comm) ... subroutine par_GlobalSumDouble(x,y,k,primme_svds) use iso_c_binding implicit none ... MPI_Comm, pointer :: comm type(c_ptr) :: pcomm call primme_svds_get_member_f77(primme_svds, PRIMME_SVDS_commInfo, pcomm) call c_f_pointer(pcomm, comm) call MPI_Allreduce(x,y,k,MPI_DOUBLE,MPI_SUM,comm,ierr)
Most users would not need to retrieve these pointers in their programs.