H2Lib  3.0
krylovsolvers.h
Go to the documentation of this file.
1 /* ------------------------------------------------------------
2  * This is the file "krylovsolvers.h" of the H2Lib package.
3  * All rights reserved, Steffen Boerm 2016
4  * ------------------------------------------------------------ */
5 
9 #ifndef KRYLOVSOLVERS_H
10 #define KRYLOVSOLVERS_H
11 
12 #include "amatrix.h"
13 #include "sparsematrix.h"
14 #include "hmatrix.h"
15 #include "h2matrix.h"
16 #include "dh2matrix.h"
17 #include "krylov.h"
18 
38 solve_cg_avector(void* A, addeval_t addeval_A, pcavector b, pavector x,
39  real eps, uint maxiter);
40 
54  uint maxiter);
55 
69  real eps, uint maxiter);
70 
84  uint maxiter);
85 
99  uint maxiter);
100 
114  uint maxiter);
115 
132 solve_pcg_avector(void *A, addeval_t addeval_A, prcd_t prcd, void *pdata,
133  pcavector b, pavector x, real eps, uint maxiter);
134 
149 solve_pcg_amatrix_avector(pcamatrix A, prcd_t prcd, void *pdata, pcavector b,
150  pavector x, real eps, uint maxiter);
151 
167  pcavector b, pavector x, real eps, uint maxiter);
168 
183 solve_pcg_hmatrix_avector(pchmatrix A, prcd_t prcd, void *pdata, pcavector b,
184  pavector x, real eps, uint maxiter);
185 
200 solve_pcg_h2matrix_avector(pch2matrix A, prcd_t prcd, void *pdata, pcavector b,
201  pavector x, real eps, uint maxiter);
202 
217 solve_pcg_dh2matrix_avector(pcdh2matrix A, prcd_t prcd, void *pdata,
218  pcavector b, pavector x, real eps, uint maxiter);
219 
235 solve_gmres_avector(void* A, addeval_t addeval_A, pcavector b, pavector x,
236  real eps, uint maxiter, uint kmax);
237 
252  uint maxiter, uint kmax);
253 
268  real eps, uint maxiter, uint kmax);
269 
284  uint maxiter, uint kmax);
285 
300  uint maxiter, uint kmax);
301 
316  uint maxiter, uint kmax);
317 
334 uint
335 solve_pgmres_avector(void* A, addeval_t addeval_A, prcd_t prcd, void *pdata,
336  pcavector b, pavector x, real eps, uint maxiter, uint kmax);
337 
352 uint
354  pavector x, real eps, uint maxiter, uint kmax);
355 
370 uint
372  pcavector b, pavector x, real eps, uint maxiter, uint kmax);
373 
388 uint
390  pavector x, real eps, uint maxiter, uint kmax);
391 
406 uint
407 solve_pgmres_h2matrix_avector(pch2matrix A, prcd_t prcd, void *pdata,
408  pcavector b, pavector x, real eps, uint maxiter, uint kmax);
409 
424 uint
426  pcavector b, pavector x, real eps, uint maxiter, uint kmax);
427 
430 #endif
uint solve_cg_avector(void *A, addeval_t addeval_A, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the conjugate gradient method and a general matri...
uint solve_pcg_h2matrix_avector(pch2matrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the preconditioned conjugate gradient method...
uint solve_gmres_sparsematrix_avector(pcsparsematrix A, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the generalized minimal residual method.
uint solve_cg_h2matrix_avector(pch2matrix A, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the conjugate gradient method.
uint solve_cg_hmatrix_avector(pchmatrix A, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the conjugate gradient method.
uint solve_pgmres_hmatrix_avector(pchmatrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the preconditioned generalized minimal residual method.
uint solve_pcg_hmatrix_avector(pchmatrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the preconditioned conjugate gradient method...
Definition: avector.h:39
uint solve_gmres_avector(void *A, addeval_t addeval_A, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the generalized minimal residual method.
uint solve_pgmres_amatrix_avector(pcamatrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the preconditioned generalized minimal residual method.
uint solve_gmres_hmatrix_avector(pchmatrix A, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the generalized minimal residual method.
Representation of -matrices.
Definition: h2matrix.h:48
uint solve_gmres_amatrix_avector(pcamatrix A, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the generalized minimal residual method.
Representation of -matrices.
Definition: hmatrix.h:49
uint solve_cg_dh2matrix_avector(pcdh2matrix A, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the conjugate gradient method.
unsigned uint
Unsigned integer type.
Definition: settings.h:70
Tree structure representing a -matrix.
Definition: dh2matrix.h:50
uint solve_pgmres_sparsematrix_avector(pcsparsematrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the preconditioned generalized minimal residual method.
uint solve_pcg_sparsematrix_avector(pcsparsematrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the preconditioned conjugate gradient method...
uint solve_pgmres_h2matrix_avector(pch2matrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the preconditioned generalized minimal residual method.
uint solve_gmres_dh2matrix_avector(pcdh2matrix A, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the generalized minimal residual method.
uint solve_pcg_avector(void *A, addeval_t addeval_A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the preconditioned conjugate gradient method...
uint solve_cg_sparsematrix_avector(pcsparsematrix A, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the conjugate gradient method.
uint solve_pcg_dh2matrix_avector(pcdh2matrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the preconditioned conjugate gradient method...
#define HEADER_PREFIX
Prefix for function declarations.
Definition: settings.h:43
double real
real floating point type.
Definition: settings.h:97
uint solve_gmres_h2matrix_avector(pch2matrix A, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the generalized minimal residual method.
uint solve_pcg_amatrix_avector(pcamatrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the preconditioned conjugate gradient method...
void(* prcd_t)(void *pdata, pavector r)
Preconditioner callback.
Definition: krylov.h:75
void(* addeval_t)(field alpha, void *matrix, pcavector x, pavector y)
Matrix callback.
Definition: krylov.h:41
Representation of a sparse matrix in compressed row format.
Definition: sparsematrix.h:42
Representation of a matrix as an array in column-major order.
Definition: amatrix.h:43
uint solve_pgmres_dh2matrix_avector(pcdh2matrix A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the preconditioned generalized minimal residual method.
uint solve_pgmres_avector(void *A, addeval_t addeval_A, prcd_t prcd, void *pdata, pcavector b, pavector x, real eps, uint maxiter, uint kmax)
Solve a linear system with the preconditioned generalized minimal residual method.
uint solve_cg_amatrix_avector(pcamatrix A, pcavector b, pavector x, real eps, uint maxiter)
Solve a self-adjoint positive definite system with the conjugate gradient method.