UCFD_SPARSE  v1.1
Documentation
Loading...
Searching...
No Matches
krylov.h File Reference

Header file for Krylov subspace methods. More...

#include "mkl.h"
#include "ucfd_types.h"
#include "config.h"
Include dependency graph for krylov.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Macros

#define eps   2.22e-16
 

Functions

ucfd_status_t serial_gmres (sparse_matrix_t op, ucfd_precon_type_t precon_type, int bn, int block, int m, int *iter, double tol, int *row_ptr, int *col_ind, int *diag_ind, double *precon_nnz_data, double *x, double *b, double *H, double *V, double *g, double *y, double *w, double *r)
 Serial GMRES routine. More...
 
ucfd_status_t step_gmres (sparse_matrix_t op, ucfd_precon_solve psolve, const struct matrix_descr descr, int bn, int m, int *flag, int *row_ptr, int *col_ind, int *diag_ind, double *precon_nnz_data, double *x, double *b, double *H, double *V, double *g, double *y, double *w, double *r)
 Single GMRES iteration routine. More...
 
ucfd_status_t serial_bicgstab (sparse_matrix_t op, ucfd_precon_type_t precon_type, int bn, int *iter, double tol, int *row_ptr, int *col_ind, int *diag_ind, double *precon_nnz_data, double *x, double *b, double *r, double *p, double *v, double *s, double *t)
 Serial BiCGstab routine. More...
 

Detailed Description

Header file for Krylov subspace methods.

Definition in file krylov.h.

Macro Definition Documentation

◆ eps

#define eps   2.22e-16

Definition at line 12 of file krylov.h.

Function Documentation

◆ serial_bicgstab()

ucfd_status_t serial_bicgstab ( sparse_matrix_t  op,
ucfd_precon_type_t  precon_type,
int  bn,
int *  iter,
double  tol,
int *  row_ptr,
int *  col_ind,
int *  diag_ind,
double *  precon_nnz_data,
double *  x,
double *  b,
double *  r,
double *  p,
double *  v,
double *  s,
double *  t 
)

Serial BiCGstab routine.

Parameters
opIntel MKL sparse matrix handler (System matrix)
precon_typePreconditioner type
bnNumber of element cells (= Number of block rows)
blockBlock size of BSR matrix (= Size of blocks)
iterMaximum iteration number
tolTolerence of convergence
row_ptrRow-directional index pointer array of the matrix
col_indColumn index array of the matrix
diag_indDiagonal matrix index array based on row_ptr
precon_nnz_dataNon-zero value array of the preconditioner matrix
xSolution array
bRight-hand-side array
rResidual array
p

Overall process of BiCGstab routine. Outer iteration ends when solution is converged or reached in maximum iteration number. UCFD_STATUS_CONVERGED is returned when L-2 norm of the residual vector becomes smaller than tol, and UCFD_STATUS_NOT_CONVERGED is returned when maximum iteration finished.

Definition at line 286 of file krylov.c.

Here is the call graph for this function:

◆ serial_gmres()

ucfd_status_t serial_gmres ( sparse_matrix_t  op,
ucfd_precon_type_t  precon_type,
int  bn,
int  block,
int  m,
int *  iter,
double  tol,
int *  row_ptr,
int *  col_ind,
int *  diag_ind,
double *  precon_nnz_data,
double *  x,
double *  b,
double *  H,
double *  V,
double *  g,
double *  y,
double *  w,
double *  r 
)

Serial GMRES routine.

Parameters
opIntel MKL sparse matrix handler (System matrix)
precon_typePreconditioner type
bnNumber of element cells (= Number of block rows)
blockBlock size of BSR matrix (= Size of blocks)
mRestart number
iterMaximum iteration number
tolTolerence of convergence
row_ptrRow-directional index pointer array of the matrix
col_indColumn index array of the matrix
diag_indDiagonal matrix index array based on row_ptr
precon_nnz_dataNon-zero value array of the preconditioner matrix
xSolution array
bRight-hand-side array
HUpper Hessenberg matrix
VOrthogonal matrix
gGivens rotation prameters
yLeast squares array for argmin(||beta - Hy||)
wArnoldi iteration array
rResidual array

Overall process of GMRES routine. Outer iteration ends when solution is converged or reached in maximum iteration number. UCFD_STATUS_CONVERGED is returned when L-2 norm of the residual vector becomes smaller than tol, and UCFD_STATUS_NOT_CONVERGED is returned when maximum iteration finished.

Definition at line 46 of file krylov.c.

Here is the call graph for this function:

◆ step_gmres()

ucfd_status_t step_gmres ( sparse_matrix_t  op,
ucfd_precon_solve  psolve,
const struct matrix_descr  descr,
int  bn,
int  m,
int *  flag,
int *  row_ptr,
int *  col_ind,
int *  diag_ind,
double *  precon_nnz_data,
double *  x,
double *  b,
double *  H,
double *  V,
double *  g,
double *  y,
double *  w,
double *  r 
)

Single GMRES iteration routine.

Parameters
opIntel MKL sparse matrix handler (System matrix)
precon_typePreconditioner type
bnNumber of element cells (= Number of block rows)
blockBlock size of BSR matrix (= Size of blocks)
mRestart number
row_ptrRow-directional index pointer array of the matrix
col_indColumn index array of the matrix
diag_indDiagonal matrix index array based on row_ptr
precon_nnz_dataNon-zero value array of the preconditioner matrix
xSolution array
bRight-hand-side array
HUpper Hessenberg matrix
VOrthogonal matrix
gGivens rotation prameters
yLeast squares array for argmin(||beta - Hy||)
wArnoldi iteration array
rResidual array

Single iteration of GMRES.

Note
Residual array r must be initialized, r := b - A @ x.

Definition at line 185 of file krylov.c.