.. highlight:: fortran FORTRAN Library Interface ------------------------- .. versionadded:: 2.0 The next enumerations and functions are declared in ``primme_svds_f77.h``. sprimme_svds_f77 """""""""""""""" .. f: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 :f:func:`magma_sprimme_svds_f77` for using GPUs). :param svals(*): (output) array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: real :param svecs(*): array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: real :param resNorms(*): array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: real :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. cprimme_svds_f77 """""""""""""""" .. f: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 :f:func:`magma_cprimme_svds_f77` for using GPUs). :param svals(*): (output) array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: real :param svecs(*): array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: complex :param resNorms(*): array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: real :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. dprimme_svds_f77 """""""""""""""" .. f: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 :f:func:`magma_dprimme_svds_f77` for using GPUs). :param svals(*): (output) array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: double precision :param svecs(*): array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: double precision :param resNorms(*): array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: double precision :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. zprimme_svds_f77 """""""""""""""" .. f: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 :f:func:`magma_zprimme_svds_f77` for using GPUs). :param svals(*): (output) array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: double precision :param svecs(*): array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: complex double precision :param resNorms(*): array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: double precision :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. magma_sprimme_svds_f77 """""""""""""""""""""" .. f: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 :f:func:`sprimme_svds_f77` for using only the CPU). :param svals(*): (output) CPU array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: real :param svecs(*): GPU array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: real :param resNorms(*): CPU array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: real :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. .. versionadded:: 3.0 magma_cprimme_svds_f77 """""""""""""""""""""" .. f: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 :f:func:`cprimme_svds_f77` for using only the CPU). :param svals(*): (output) CPU array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: real :param svecs(*): GPU array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: complex :param resNorms(*): CPU array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: real :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. .. versionadded:: 3.0 magma_dprimme_svds_f77 """""""""""""""""""""" .. f: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 :f:func:`dprimme_svds_f77` for using only the CPU). :param svals(*): (output) CPU array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: double precision :param svecs(*): GPU array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: double precision :param resNorms(*): CPU array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: double precision :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. .. versionadded:: 3.0 magma_zprimme_svds_f77 """""""""""""""""""""" .. f: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 :f:func:`zprimme_svds_f77` for using only the CPU). :param svals(*): (output) CPU array at least of size |SnumSvals| to store the computed singular values; all processes in a parallel run return this local array with the same values. :type svals: double precision :param svecs(*): GPU array at least of size (|SmLocal| + |SnLocal|) times (|SnumOrthoConst| + |SnumSvals|) to store column-wise the (local part for this process of the) computed left singular vectors and the right singular vectors. :type svecs: complex double precision :param resNorms(*): CPU array at least of size |SnumSvals| to store the residual norms of the computed triplets; all processes in parallel run return this local array with the same values. :type resNorms: double precision :param ptr primme_svds: parameters structure. :param integer ierr: (output) error indicator; see :ref:`error-codes-svds`. On input, ``svecs`` should start with the content of the |SnumOrthoConst| left vectors, followed by the |SinitSize| left vectors, followed by the |SnumOrthoConst| right vectors and followed by the |SinitSize| right vectors. On return, the i-th left singular vector starts at svecs(( |SnumOrthoConst| + i - 1) \* |SmLocal| ). The i-th right singular vector starts at svecs(( |SnumOrthoConst| + |SinitSize| )\* |SmLocal| + ( |SnumOrthoConst| + i - 1)\* |SnLocal| ). The first vector has i=1. All internal operations are performed at the same precision than ``svecs`` unless the user sets |SinternalPrecision| otherwise. The type and precision of the callbacks depends on the type and precision of `svecs`. See details for |SmatrixMatvec|, |SapplyPreconditioner|, |SglobalSumReal|, |SbroadcastReal|, and |SconvTestFun|. .. versionadded:: 3.0 primme_svds_initialize_f77 """""""""""""""""""""""""" .. f:subroutine:: primme_svds_initialize_f77(primme_svds, ierr) Set PRIMME SVDS parameters structure to the default values. After calling :f:func:`dprimme_svds_f77` (or a variant), call :f:func:`primme_svds_free_f77` to release allocated resources by PRIMME. :param ptr primme_svds: (output) parameters structure. primme_svds_set_method_f77 """""""""""""""""""""""""" .. f:subroutine:: primme_svds_set_method_f77(method, methodStage1, methodStage2, primme_svds, ierr) Set PRIMME SVDS parameters to one of the preset configurations. :param integer method: (input) preset configuration to compute the singular triplets; one of * |PRIMME_SVDS_default|, currently set as |PRIMME_SVDS_hybrid|. * |PRIMME_SVDS_normalequations|, compute the eigenvectors of :math:`A^*A` or :math:`A A^*`. * |PRIMME_SVDS_augmented|, compute the eigenvectors of the augmented matrix, :math:`\left(\begin{array}{cc} 0 & A^* \\ A & 0 \end{array}\right)`. * |PRIMME_SVDS_hybrid|, start with |PRIMME_SVDS_normalequations|; use the resulting approximate singular vectors as initial vectors for |PRIMME_SVDS_augmented| if the required accuracy was not achieved. :param primme_preset_method methodStage1: (input) preset method to compute the eigenpairs at the first stage; see available values at :f:func:`primme_set_method_f77`. :param primme_preset_method methodStage2: (input) preset method to compute the eigenpairs with the second stage of |PRIMME_SVDS_hybrid|; see available values at :f:func:`primme_set_method_f77`. :param ptr primme_svds: (input/output) parameters structure. :param integer ierr: (output) if 0, successful; if negative, something went wrong. primme_svds_display_params_f77 """""""""""""""""""""""""""""" .. f:subroutine:: primme_svds_display_params_f77(primme_svds) Display all printable settings of ``primme_svds`` into the file descriptor |SoutputFile|. :param ptr primme_svds: (input) parameters structure. primme_svds_free_f77 """""""""""""""""""" .. f:subroutine:: primme_svds_free_f77(primme_svds, ierr) Free memory allocated by PRIMME SVDS and delete all values set. :param ptr primme_svds: (input/output) parameters structure. primme_svds_set_member_f77 """""""""""""""""""""""""" .. f:subroutine:: primme_svds_set_member_f77(primme_svds, label, value, ierr) Set a value in some field of the parameter structure. :param ptr primme_svds: (input) parameters structure. :param integer label: field where to set value. One of: | :c:member:`PRIMME_SVDS_primme ` | :c:member:`PRIMME_SVDS_primmeStage2 ` | :c:member:`PRIMME_SVDS_m ` | :c:member:`PRIMME_SVDS_n ` | :c:member:`PRIMME_SVDS_matrixMatvec ` | :c:member:`PRIMME_SVDS_matrixMatvec_type ` | :c:member:`PRIMME_SVDS_applyPreconditioner ` | :c:member:`PRIMME_SVDS_applyPreconditioner_type ` | :c:member:`PRIMME_SVDS_numProcs ` | :c:member:`PRIMME_SVDS_procID ` | :c:member:`PRIMME_SVDS_mLocal ` | :c:member:`PRIMME_SVDS_nLocal ` | :c:member:`PRIMME_SVDS_commInfo ` | :c:member:`PRIMME_SVDS_globalSumReal ` | :c:member:`PRIMME_SVDS_globalSumReal_type ` | :c:member:`PRIMME_SVDS_broadcastReal ` | :c:member:`PRIMME_SVDS_broadcastReal_type ` | :c:member:`PRIMME_SVDS_numSvals ` | :c:member:`PRIMME_SVDS_target ` | :c:member:`PRIMME_SVDS_numTargetShifts ` | :c:member:`PRIMME_SVDS_targetShifts ` | :c:member:`PRIMME_SVDS_method ` | :c:member:`PRIMME_SVDS_methodStage2 ` | :c:member:`PRIMME_SVDS_matrix ` | :c:member:`PRIMME_SVDS_preconditioner ` | :c:member:`PRIMME_SVDS_locking ` | :c:member:`PRIMME_SVDS_numOrthoConst ` | :c:member:`PRIMME_SVDS_aNorm ` | :c:member:`PRIMME_SVDS_eps ` | :c:member:`PRIMME_SVDS_precondition ` | :c:member:`PRIMME_SVDS_initSize ` | :c:member:`PRIMME_SVDS_maxBasisSize ` | :c:member:`PRIMME_SVDS_maxBlockSize ` | :c:member:`PRIMME_SVDS_maxMatvecs ` | :c:member:`PRIMME_SVDS_iseed ` | :c:member:`PRIMME_SVDS_printLevel ` | :c:member:`PRIMME_SVDS_outputFile ` | :c:member:`PRIMME_SVDS_internalPrecision ` | :c:member:`PRIMME_SVDS_convTestFun ` | :c:member:`PRIMME_SVDS_convTestFun_type ` | :c:member:`PRIMME_SVDS_convtest ` | :c:member:`PRIMME_SVDS_monitorFun ` | :c:member:`PRIMME_SVDS_monitorFun_type ` | :c:member:`PRIMME_SVDS_monitor ` | :c:member:`PRIMME_SVDS_queue ` | :c:member:`PRIMME_SVDS_stats_numOuterIterations ` | :c:member:`PRIMME_SVDS_stats_numRestarts ` | :c:member:`PRIMME_SVDS_stats_numMatvecs ` | :c:member:`PRIMME_SVDS_stats_numPreconds ` | :c:member:`PRIMME_SVDS_stats_numGlobalSum ` | :c:member:`PRIMME_SVDS_stats_numBroadcast ` | :c:member:`PRIMME_SVDS_stats_volumeGlobalSum ` | :c:member:`PRIMME_SVDS_stats_volumeBroadcast ` | :c:member:`PRIMME_SVDS_stats_elapsedTime ` | :c:member:`PRIMME_SVDS_stats_timeMatvec ` | :c:member:`PRIMME_SVDS_stats_timePrecond ` | :c:member:`PRIMME_SVDS_stats_timeOrtho ` | :c:member:`PRIMME_SVDS_stats_timeGlobalSum ` | :c:member:`PRIMME_SVDS_stats_timeBroadcast ` | :c:member:`PRIMME_SVDS_stats_lockingIssue ` :param value: (input) value to set. .. note:: **Don't use** this subroutine inside PRIMME SVDS's callback functions, e.g., |SmatrixMatvec| or |SapplyPreconditioner|, or in functions called by these functions. primme_svdstop_get_member_f77 """"""""""""""""""""""""""""" .. f:subroutine:: primme_svdstop_get_member_f77(primme_svds, label, value, ierr) Get the value in some field of the parameter structure. :param ptr primme_svds: (input) parameters structure. :param integer label: (input) field where to get value. One of the detailed in subroutine :f:func:`primmesvds_top_set_member_f77`. :param value: (output) value of the field. .. note:: **Don't use** this subroutine inside PRIMME SVDS's callback functions, e.g., |SmatrixMatvec| or |SapplyPreconditioner|, or in functions called by these functions. In those cases use :f:func:`primme_svds_get_member_f77`. .. note:: When ``label`` is one of ``PRIMME_SVDS_matrixMatvec``, ``PRIMME_SVDS_applyPreconditioner``, ``PRIMME_SVDS_commInfo``, ``PRIMME_SVDS_intWork``, ``PRIMME_SVDS_realWork``, ``PRIMME_SVDS_matrix`` and ``PRIMME_SVDS_preconditioner``, the returned ``value`` 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 """""""""""""""""""""""""" .. f:subroutine:: primme_svds_get_member_f77(primme_svds, label, value, ierr) Get the value in some field of the parameter structure. :param ptr primme_svds: (input) parameters structure. :param integer label: (input) field where to get value. One of the detailed in subroutine :f:func:`primme_svdstop_set_member_f77`. :param value: (output) value of the field. .. note:: Use this subroutine exclusively inside PRIMME SVDS's callback functions, e.g., |SmatrixMatvec| or |SapplyPreconditioner|, or in functions called by these functions. Otherwise, e.g., from the main program, use the subroutine :f:func:`primme_svdstop_get_member_f77`. .. note:: When ``label`` is one of ``PRIMME_SVDS_matrixMatvec``, ``PRIMME_SVDS_applyPreconditioner``, ``PRIMME_SVDS_commInfo``, ``PRIMME_SVDS_intWork``, ``PRIMME_SVDS_realWork``, ``PRIMME_SVDS_matrix`` and ``PRIMME_SVDS_preconditioner``, the returned ``value`` 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. .. include:: epilog.inc