Visible to Intel only — GUID: GUID-56F0C6ED-0A49-48C4-AF7B-33C3851A39CA
Visible to Intel only — GUID: GUID-56F0C6ED-0A49-48C4-AF7B-33C3851A39CA
FGMRES Interface Description
All types in this documentation refer to the common Fortran types: INTEGER and DOUBLE PRECISION.
Routine Options
All of the RCI FGMRES routines have common parameters for passing various options to the routines (see FGMRES Common Parameters). The values for these parameters can be changed during computations.
User Data Arrays
Many of the RCI FGMRES routines take arrays of user data as input. For example, user arrays are passed to the routine dfgmresto compute the solution of a system of linear algebraic equations. To minimize storage requirements and improve overall run-time efficiency, the Intel® oneAPI Math Kernel Library (oneMKL) RCI FGMRES routines do not make copies of the user input arrays.
FGMRES Common Parameters
The default and initial values listed below are assigned to the parameters by calling the dfgmres_init routine.
- n
-
INTEGER, this parameter sets the size of the problem in the dfgmres_init routine. All the other routines use the ipar(1) parameter instead. Note that the coefficient matrix A is a square matrix of size n*n.
- x
-
DOUBLE PRECISION array, this parameter contains the current approximation to the solution vector. Before the first call to the dfgmres routine, it contains the initial approximation to the solution vector.
- b
-
DOUBLE PRECISION array, this parameter contains the right-hand side vector. Depending on user requests (see the parameter ipar(13)), it might contain the approximate solution after execution.
- RCI_request
-
INTEGER, this parameter gives information about the result of work of the RCI FGMRES routines. Negative values of the parameter indicate that the routine completed with errors or warnings. The 0 value indicates successful completion of the task. Positive values mean that you must perform specific actions:
RCI_request= 1
multiply the matrix by tmp(ipar(22)), put the result in tmp(ipar(23)), and return the control to the dfgmres routine;
RCI_request= 2
perform the stopping tests. If they fail, return the control to the dfgres routine. Otherwise, the solution can be updated by a subsequent call to dfgmres_get routine;
RCI_request= 3
apply the preconditioner to tmp(ipar(22)), put the result in tmp(ipar(23)), and return the control to the dfgmres routine.
RCI_request= 4
check if the norm of the current orthogonal vector is zero, within the rounding or computational errors. Return the control to the dfgmres routine if it is not zero, otherwise complete the solution process by calling dfgmres_get routine.
- ipar(128)[128]
-
INTEGER array, this parameter specifies the integer set of data for the RCI FGMRES computations:
ipar(1)
specifies the size of the problem. The dfgmres_init routine assigns ipar(1)=n. All the other routines uses this parameter instead of n. There is no default value for this parameter.
ipar(2)
specifies the type of output for error and warning messages that are generated by the RCI FGMRES routines. The default value 6 means that all messages are displayed on the screen. Otherwise the error and warning messages are written to the newly created file MKL_RCI_FGMRES_Log.txt. Note that if ipar(6) and ipar(7) parameters are set to 0, error and warning messages are not generated at all.
ipar(3)
contains the current stage of the RCI FGMRES computations. The initial value is 1.
warning:Avoid altering this variable during computations.
ipar(4)
contains the current iteration number. The initial value is 0.
ipar(5)
specifies the maximum number of iterations. The default value is min (150,n).
ipar(6)
if the value is not 0, the routines output error messages in accordance with the parameter ipar(2). If it is 0, the routines do not output error messages at all, but return a negative value of the parameter RCI_request. The default value is 1.
ipar(7)
if the value is not 0, the routines output warning messages in accordance with the parameter ipar(2). Otherwise, the routines do not output warning messages at all, but they return a negative value of the parameter RCI_request. The default value is 1.
ipar(8)
if the value is not equal to 0, the dfgmres routine performs the stopping test for the maximum number of iterations: ipar(4)≤ipar(5). If the value is 0, the dfgmres routine does not perform this stopping test. The default value is 1.
ipar(9)
if the value is not 0, the dfgmres routine performs the residual stopping test: dpar(5)≤dpar(4).If the value is 0, the dfgmres routine does not perform this stopping test. The default value is 0.
ipar(10)
if the value is not 0, the dfgmres routine indicates that the user-defined stopping test should be performed by setting RCI_request=2. If the value is 0, the dfgmres routine does not perform the user-defined stopping test. The default value is 1.
NOTE:At least one of the parameters ipar(8)-ipar(10) must be set to 1.
ipar(11)
if the value is 0, the dfgmres routine runs the non-preconditioned version of the FGMRES method. Otherwise, the routine runs the preconditioned version of the FGMRES method, and requests that you perform the preconditioning step by setting the output parameter RCI_request=3. The default value is 0.
ipar(12)
if the value is not equal to 0, the dfgmres routine performs the automatic test for zero norm of the currently generated vector: dpar(7)≤dpar(8), where dpar(8) contains the tolerance value. Otherwise, the routine indicates that you must perform this check by setting the output parameter RCI_request=4. The default value is 0.
ipar(13)
if the value is equal to 0, the dfgmres_get routine updates the solution to the vector x according to the computations done by the dfgmres routine. If the value is positive, the routine writes the solution to the right-hand side vector b. If the value is negative, the routine returns only the number of the current iteration, and does not update the solution. The default value is 0.
NOTE:It is possible to call the dfgmres_get routine at any place in the code, but you must pay special attention to the parameter ipar(13). The RCI FGMRES iterations can be continued after the call to dfgmres_get routine only if the parameter ipar(13) is not equal to zero. If ipar(13) is positive, then the updated solution overwrites the right-hand side in the vector b. If you want to run the restarted version of FGMRES with the same right-hand side, then it must be saved in a different memory location before the first call to the dfgmres_get routine with positive ipar(13).
ipar(14)
contains the internal iteration counter that counts the number of iterations before the restart takes place. The initial value is 0.
warning:Do not alter this variable during computations.
ipar(15)
specifies the number of the non-restarted FGMRES iterations. To run the restarted version of the FGMRES method, assign the number of iterations to ipar(15) before the restart. The default value is min(150, n), which means that by default the non-restarted version of FGMRES method is used.
ipar(16)
service variable specifying the location of the rotated Hessenberg matrix from which the matrix stored in the packed format (see Matrix Arguments in the Appendix B for details) is started in the tmp array.
ipar(17)
service variable specifying the location of the rotation cosines from which the vector of cosines is started in the tmp array.
ipar(18)
service variable specifying the location of the rotation sines from which the vector of sines is started in the tmp array.
ipar(19)
service variable specifying the location of the rotated residual vector from which the vector is started in the tmp array.
ipar(20)
service variable, specifies the location of the least squares solution vector from which the vector is started in the tmp array.
ipar(21)
service variable specifying the location of the set of preconditioned vectors from which the set is started in the tmp array. The memory locations in the tmp array starting from ipar(21) are used only for the preconditioned FGMRES method.
ipar(22)
specifies the memory location from which the first vector (source) used in operations requested via RCI_request is started in the tmp array.
ipar(23)
specifies the memory location from which the second vector (output) used in operations requested via RCI_request is started in the tmp array.
ipar(24:128)
are reserved and not used in the current RCI FGMRES routines.
NOTE:You must declare the array ipar with length 128. While defining the array in the code as INTEGERipar(23)works, there is no guarantee of future compatibility with Intel® oneAPI Math Kernel Library (oneMKL).
- dpar(128)
-
DOUBLE PRECISION array, this parameter specifies the double precision set of data for the RCI CG computations, specifically:
dpar(1)
specifies the relative tolerance. The default value is 1.0e-6.
dpar(2)
specifies the absolute tolerance. The default value is 0.0e-0.
dpar(3)
specifies the Euclidean norm of the initial residual (if it is computed in the dfgmres routine). The initial value is 0.0.
dpar(4)
service variable equal to dpar(1)*dpar(3)+dpar(2) (if it is computed in the dfgmres routine). The initial value is 0.0.
dpar(5)
specifies the Euclidean norm of the current residual. The initial value is 0.0.
dpar(6)
specifies the Euclidean norm of residual from the previous iteration step (if available). The initial value is 0.0.
dpar(7)
contains the norm of the generated vector. The initial value is 0.0.
NOTE:In terms of [Saad03] this parameter is the coefficient hk+1,k of the Hessenberg matrix.
dpar(8)
contains the tolerance for the zero norm of the currently generated vector. The default value is 1.0e-12.
dpar(9:128)
are reserved and not used in the current RCI FGMRES routines.
NOTE:You must declare the array dpar with length 128. While defining the array in the code as DOUBLE PRECISION dpar(8)works, there is no guarantee of future compatibility with Intel® oneAPI Math Kernel Library (oneMKL).
- tmp
-
DOUBLE PRECISION array of size ((2*ipar(15)+1)*n + ipar(15)*(ipar(15)+9)/2 + 1)) used to supply the double precision temporary space for the RCI FGMRES computations, specifically:
tmp(1:ipar(16)-1)
contains the sequence of vectors generated by the FGMRES method. The initial value is 0.0.
tmp(ipar(16):ipar(17)-1)
contains the rotated Hessenberg matrix generated by the FGMRES method; the matrix is stored in the packed format. There is no initial value for this part of tmp array.
tmp(ipar(17):ipar(18)-1)
contains the rotation cosines vector generated by the FGMRES method. There is no initial value for this part of tmp array.
tmp(ipar(18):ipar(19)-1)
contains the rotation sines vector generated by the FGMRES method. There is no initial value for this part of tmp array.
tmp(ipar(19):ipar(20)-1)
contains the rotated residual vector generated by the FGMRES method. There is no initial value for this part of tmp array.
tmp(ipar(20):ipar(21)-1)
contains the solution vector to the least squares problem generated by the FGMRES method. There is no initial value for this part of tmp array.
tmp(ipar(21):)
contains the set of preconditioned vectors generated for the FGMRES method by the user. This part of tmp array is not used if the non-preconditioned version of FGMRES method is called. There is no initial value for this part of tmp array.
NOTE:You can define this array in the code as DOUBLE PRECISION tmp((2*ipar(15)+1)*n + ipar(15)*(ipar(15)+9)/2 + 1)) if you run only non-preconditioned FGMRES iterations.
Scheme of Using the RCI FGMRES Routines
The following pseudocode shows the general scheme of using the RCI FGMRES routines.
...
generate matrix A
generate preconditioner C (optional)
call dfgmres_init(n, x, b, RCI_request, ipar, dpar, tmp)
change parameters in ipar, dpar if necessary
call dfgmres_check(n, x, b, RCI_request, ipar, dpar, tmp)
1 call dfgmres(n, x, b, RCI_request, ipar, dpar, tmp)
if (RCI_request.eq.1) then
multiply the matrix A by tmp(ipar(22)) and put the result in tmp(ipar(23))
It is possible to use MKL Sparse BLAS Level 2 subroutines for this operation
c proceed with FGMRES iterations
goto 1
endif
if (RCI_request.eq.2) then
do the stopping test
if (test not passed) then
c proceed with FGMRES iterations
go to 1
else
c stop FGMRES iterations
goto 2
endif
endif
if (RCI_request.eq.3) then (optional)
apply the preconditioner C inverse to tmp(ipar(22)) and put the result in tmp(ipar(23))
c proceed with FGMRES iterations
goto 1
endif
if (RCI_request.eq.4) then
check the norm of the next orthogonal vector, it is contained in dpar(7)
if (the norm is not zero up to rounding/computational errors) then
c proceed with FGMRES iterations
goto 1
else
c stop FGMRES iterations
goto 2
endif
endif
2 call dfgmres_get(n, x, b, RCI_request, ipar, dpar, tmp, itercount)
current iteration number is in itercount
the computed approximation is in the array x
For the FGMRES method, the array x initially contains the current approximation to the solution. It can be updated only by calling the routine dfgmres_get, which updates the solution in accordance with the computations performed by the routine dfgmres.
The above pseudocode demonstrates two main differences in the use of RCI FGMRES interface comparing with the CG Interface Description. The first difference relates to RCI_request=3: it uses different locations in the tmp array, which is two-dimensional for CG and one-dimensional for FGMRES. The second difference relates to RCI_request=4: the RCI CG interface never produces RCI_request=4.