dfgmres

Makes the FGMRES iterations.

Syntax

dfgmres(n, x, b, RCI_request, ipar, dpar, tmp)

Include Files

The Fortran interface is specified in the mkl_rci.fi include file and the C interface is specified in the mkl_rci.h include file.

Description

The routine dfgmres performs the FGMRES iterations [Saad03]. The routine dfgmres uses the value that was in the vector x before the first call as an initial approximation to the solution. To update the current approximation to the solution, the dfgmres_get routine must be called. The RCI FGMRES iterations can be continued after the call to the dfgmres_get routine only if the value of the parameter ipar(13) is not equal to 0 (default value). Note that the updated solution overwrites the right hand side in the vector b if the parameter ipar(13) is positive, and the restarted version of the FGMRES method can not be run. If the right hand side should be kept, it must be saved in a different memory location before the first call to the dfgmres_get routine with a positive ipar(13).

The parameter RCI_request informs the user about the task completion and requests results of certain operations that are needed to the solver.

Note that lengths of all the vectors are assumed to have been defined in a previous call to the dfgmres_init routine.

Input Parameters

n

INTEGER. Contains the size of the problem, and the sizes of arrays x and b.

x

DOUBLE PRECISION array of size n. Contains the initial approximation to the solution vector.

b

DOUBLE PRECISION array of size n. Contains the right-hand side vector.

tmp

DOUBLE PRECISION array of size ((2*ipar(15)+1)*n+ipar(15)*ipar(15)+9)/2 + 1). Refer to the FGMRES Common Parameters.

Output Parameters

RCI_request

INTEGER. Informs about result of work of the routine.

ipar

INTEGER array of size 128. Refer to the FGMRES Common Parameters.

dpar

DOUBLE PRECISION array of size 128. Refer to the FGMRES Common Parameters.

tmp

DOUBLE PRECISION array of size ((2*ipar(15)+1)*n+ipar(15)*ipar(15)+9)/2 + 1). Refer to the FGMRES Common Parameters.

Return Values

RCI_request= 0

The routine completes task normally, the solution is found. This occurs only if the stopping tests are fully automatic. For the user defined stopping tests, see the comments to the RCI_request= 2 or 4.

RCI_request= -1

The routine is interrupted because the maximal number of iterations is reached, but the relative stopping criterion is not met. This situation occurs only if both tests are requested by the user.

RCI_request= -10

The routine is interrupted because attempt to divide by zero occurs. Normally it happens if the matrix is (almost) degenerate. However it may happen if the parameter dpar is altered by mistake, or if the method is not stopped when the solution is found.

RCI_request= -11

The routine is interrupted because it enters the infinite cycle. Probably, the values ipar(8), ipar(9), ipar(10) were altered outside of the routine, or the dfgmres_check routine was not called.

RCI_request= -12

The routine is interrupted because errors are found in the method parameters. Normally this happens if the parameters ipar and dpar are altered by mistake outside the routine.

RCI_request= 1

Requests the user to multiply the matrix by tmp(ipar(22)), put the result in the tmp(ipar(23)), and return the control back to the routine dfgmres.

RCI_request= 2

Requests the user to perform the stopping tests. If they fail, the control back must be returned to the dfgmres routine. Otherwise, the FGMRES solution is found, and the fgmres_get routine can be run to update the computed solution in the vector x.

RCI_request= 3

Requests the user to apply the inverse preconditioner to ipar(22), put the result in the ipar(23), and return the control back to the routine dfgmres.

RCI_request= 4

Requests the user to check the norm of the currently generated vector. If it is not zero within the computational/rounding errors, the control must be returned back to the dfgmres routine. Otherwise, the FGMRES solution is found, and the dfgmres_get routine can be run to update the computed solution in the vector x.


Submit feedback on this help topic

Copyright © 1994 - 2010, Intel Corporation. All rights reserved.