next up previous contents index
Next: VECTOR(3LAS) Up: Manual Pages Previous: QMATRIX(3LAS)

  RTC(3LAS)

NAME

SetRTCAccuracy, SetRTCAuxProc, RTCResult, GetLastNoIter, GetLastAccuracy -- residual termination control of iterative solvers

SYNOPSIS

#include <laspack/rtc.h>

typedef enum { 
    JacobiIterId, 
    SORForwIterId, 
    SORBackwIterId, 
    SSORIterId, 

    ChebyshevIterId, 

    CGIterId, 
    CGNIterId, 
    GMRESIterId, 
    BiCGIterId, 
    QMRIterId, 
    CGSIterId, 
    BiCGSTABIterId, 

    MGIterId, 
    NestedMGIterId, 
    MGPCGIterId, 
    BPXPCGIterId 
} IterIdType;

void SetRTCAccuracy(double Eps); 
void SetRTCAuxProc(RTCAuxProcType AuxProc); 
Boolean RTCResult(int Iter, double rNorm, double bNorm, IterIdType IterId); 
int GetLastNoIter(void); 
double GetLastAccuracy(void);

DESCRIPTION

 
The procedure SetRTCAccuracy sets accuracy for iterative methods implemented in the modules ITERSOLV and MLSOLV. The default is 1e-8.

 
The procedure SetRTCAuxProc specifies a user defined routine which should be performed after each iteration step (exactly after each call of RTCResult). The prototype of such a procedure is
  void (*RTCAuxProcType)(int, double, double, IterIdType).
The parameters are identically with those of the function RTCResult (see below).

 
The procedure RTCResult is used internally by iterative solvers to check the termination condition
    || r ||_2 = <= eps || b ||_2,

where || r ||_2 = rNorm is the norm of the residuum, || b ||_2 = rNorm is the norm of the right hand side vector, and eps is the given accuracy. If this condition is fulfilled, the procedure returns True otherwise False. Furthermore, the number of performed iterations Iter and the identifier of the iteration method IterId (see synopsis above) are registered and passed to the auxiliary procedure defined by SetRTCAuxProc.

 
The procedure GetLastNoIter returns the number of performed iterations transmitted in the last call to the procedure RTCResult.

 
The procedure GetLastAccuracy returns the "accuracy" || r ||_2 / || b ||_2 transmitted in the last call to the procedure RTCResult.

REFERENCES

The stopping criterion applied is comprehensively discussed in:

R. Barrett, M. Berry, T. F. Chan, J. Demmel, J. Donato, J. Dongarra, V. Eijkhout, R. Pozo, Ch. Romine, H. van der Vorst: Templates for the Solution of Linear Systems: Building Blocks for Iterative Solvers, SIAM, Philadelphia, 1994. (HTML version is here).

FILES

rtc.h ... header file
rtc.c ... source file

EXAMPLES

A procedure which prints the norm of scaled residuals after each multigrid iteration could be defined in an application code as follows:

void PrintMGAccuracy(int Iter, double rNorm, double bNorm, IterIdType IterId);
    
int main(void)
{
    ...

    SetRTCAuxProc(PrintMGAccuracy);

    ...

    MGIter(...);

    ...
}

void PrintMGAccuracy(int Iter, double rNorm, double bNorm, IterIdType IterId)
/* put out accuracy reached after each multigrid iteration */
{
    if (IterId == MGIterId) {
        printf("%3d. iteration ... accuracy = ", Iter);
        if (!IsZero(bNorm))
            printf("%12.5e\n", rNorm / bNorm);
        else
            printf("    ---\n");
    }
}

SEE ALSO

itersolv(3LAS), mlsolv(3LAS)



next up previous contents index
Next: VECTOR(3LAS) Up: Manual Pages Previous: QMATRIX(3LAS)



Tomas Skalicky (skalicky@msmfs1.mw.tu-dresden.de)