Fast Auxiliary Space Preconditioning 2.7.7 Aug/28/2022
Loading...
Searching...
No Matches
KryPvgmres.c File Reference

Krylov subspace methods – Preconditioned variable-restart GMRes. More...

#include "fasp.h"
#include "fasp_functs.h"
#include <math.h>
#include "KryUtil.inl"

Go to the source code of this file.

Functions

INT fasp_solver_dcsr_pvgmres (dCSRmat *A, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
 Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration. More...
 
INT fasp_solver_dbsr_pvgmres (dBSRmat *A, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
 Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration. More...
 
INT fasp_solver_dblc_pvgmres (dBLCmat *A, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
 Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration. More...
 
INT fasp_solver_dstr_pvgmres (dSTRmat *A, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
 Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration. More...
 
INT fasp_solver_pvgmres (mxv_matfree *mf, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, SHORT restart, const SHORT StopType, const SHORT PrtLvl)
 Solve "Ax=b" using PGMRES(right preconditioned) iterative method in which the restart parameter can be adaptively modified during iteration. More...
 

Detailed Description

Krylov subspace methods – Preconditioned variable-restart GMRes.

Note
This file contains Level-3 (Kry) functions. It requires: AuxArray.c, AuxMemory.c, AuxMessage.c, BlaArray.c, BlaSpmvBLC.c, BlaSpmvBSR.c, BlaSpmvCSR.c, and BlaSpmvSTR.c
See KrySPvgmres.c for a safer version

Reference: A.H. Baker, E.R. Jessup, and Tz.V. Kolev A Simple Strategy for Varying the Restart Parameter in GMRES(m) Journal of Computational and Applied Mathematics, 230 (2009) pp. 751-761. UCRL-JRNL-235266.


Copyright (C) 2010–Present by the FASP team. All rights reserved.

Released under the terms of the GNU Lesser General Public License 3.0 or later.

Definition in file KryPvgmres.c.

Function Documentation

◆ fasp_solver_dblc_pvgmres()

INT fasp_solver_dblc_pvgmres ( dBLCmat A,
dvector b,
dvector x,
precond pc,
const REAL  tol,
const REAL  abstol,
const INT  MaxIt,
const SHORT  restart,
const SHORT  StopType,
const SHORT  PrtLvl 
)

Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration.

Parameters
APointer to dCSRmat: coefficient matrix
bPointer to dvector: right hand side
xPointer to dvector: unknowns
pcPointer to precond: structure of precondition
tolTolerance for relative residual
abstolTolerance for absolute residual
MaxItMaximal number of iterations
restartRestarting steps
StopTypeStopping criteria type
PrtLvlHow much information to print out
Returns
Iteration number if converges; ERROR otherwise.
Author
Chensong Zhang
Date
04/05/2013

Definition at line 765 of file KryPvgmres.c.

◆ fasp_solver_dbsr_pvgmres()

INT fasp_solver_dbsr_pvgmres ( dBSRmat A,
dvector b,
dvector x,
precond pc,
const REAL  tol,
const REAL  abstol,
const INT  MaxIt,
const SHORT  restart,
const SHORT  StopType,
const SHORT  PrtLvl 
)

Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration.

Parameters
APointer to dCSRmat: coefficient matrix
bPointer to dvector: right hand side
xPointer to dvector: unknowns
pcPointer to precond: structure of precondition
tolTolerance for relative residual
abstolTolerance for absolute residual
MaxItMaximal number of iterations
restartRestarting steps
StopTypeStopping criteria type
PrtLvlHow much information to print out
Returns
Iteration number if converges; ERROR otherwise.
Author
Zhiyang Zhou
Date
12/21/2011

Modified by Chensong Zhang on 04/06/2013: Add stop type support

Definition at line 416 of file KryPvgmres.c.

◆ fasp_solver_dcsr_pvgmres()

INT fasp_solver_dcsr_pvgmres ( dCSRmat A,
dvector b,
dvector x,
precond pc,
const REAL  tol,
const REAL  abstol,
const INT  MaxIt,
const SHORT  restart,
const SHORT  StopType,
const SHORT  PrtLvl 
)

Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration.

Parameters
APointer to dCSRmat: coefficient matrix
bPointer to dvector: right hand side
xPointer to dvector: unknowns
pcPointer to precond: structure of precondition
tolTolerance for relative residual
abstolTolerance for absolute residual
MaxItMaximal number of iterations
restartRestarting steps
StopTypeStopping criteria type
PrtLvlHow much information to print out
Returns
Iteration number if converges; ERROR otherwise.
Author
Zhiyang Zhou
Date
2010/12/14

Modified by Chensong Zhang on 04/06/2013: Add stop type support Modified by Chunsheng Feng on 07/22/2013: Add adapt memory allocate

Definition at line 66 of file KryPvgmres.c.

◆ fasp_solver_dstr_pvgmres()

INT fasp_solver_dstr_pvgmres ( dSTRmat A,
dvector b,
dvector x,
precond pc,
const REAL  tol,
const REAL  abstol,
const INT  MaxIt,
const SHORT  restart,
const SHORT  StopType,
const SHORT  PrtLvl 
)

Right preconditioned GMRES method in which the restart parameter can be adaptively modified during iteration.

Parameters
APointer to dCSRmat: coefficient matrix
bPointer to dvector: right hand side
xPointer to dvector: unknowns
pcPointer to precond: structure of precondition
tolTolerance for relative residual
abstolTolerance for absolute residual
MaxItMaximal number of iterations
restartRestarting steps
StopTypeStopping criteria type
PrtLvlHow much information to print out
Returns
Iteration number if converges; ERROR otherwise.
Author
Zhiyang Zhou
Date
2010/12/14

Modified by Chensong Zhang on 04/06/2013: Add stop type support

Definition at line 1116 of file KryPvgmres.c.

◆ fasp_solver_pvgmres()

INT fasp_solver_pvgmres ( mxv_matfree mf,
dvector b,
dvector x,
precond pc,
const REAL  tol,
const REAL  abstol,
const INT  MaxIt,
SHORT  restart,
const SHORT  StopType,
const SHORT  PrtLvl 
)

Solve "Ax=b" using PGMRES(right preconditioned) iterative method in which the restart parameter can be adaptively modified during iteration.

Parameters
mfPointer to mxv_matfree: spmv operation
bPointer to dvector: right hand side
xPointer to dvector: unknowns
pcPointer to precond: structure of precondition
tolTolerance for relative residual
abstolTolerance for absolute residual
MaxItMaximal number of iterations
restartRestarting steps
StopTypeStopping criteria type – DOES not support this parameter
PrtLvlHow much information to print out
Returns
Iteration number if converges; ERROR otherwise.
Author
Zhiyang Zhou
Date
2010/12/14

Modified by Feiteng Huang on 09/26/2012: matrix free Modified by Chunsheng Feng on 07/22/2013: Add adapt memory allocate

Definition at line 1468 of file KryPvgmres.c.