H2Lib  3.0
Typedefs | Functions
krylov.h File Reference
#include <stdio.h>
#include "settings.h"
#include "avector.h"
#include "factorizations.h"

Go to the source code of this file.

Typedefs

typedef void(* addeval_t) (field alpha, void *matrix, pcavector x, pavector y)
 Matrix callback. More...
 
typedef void(* mvm_t) (field alpha, bool trans, void *matrix, pcavector x, pavector y)
 Matrix callback. More...
 
typedef void(* prcd_t) (void *pdata, pavector r)
 Preconditioner callback. More...
 

Functions

real norm2_matrix (mvm_t mvm, void *A, uint rows, uint cols)
 Approximate the spectral norm $\|A\|_2$ of a matrix $A$. More...
 
real norm2diff_matrix (mvm_t mvmA, void *A, mvm_t mvmB, void *B, uint rows, uint cols)
 Approximate the spectral norm $\|A-B\|_2$ of the difference of two matrices $A$ and $B$. More...
 
real norm2diff_pre_matrix (mvm_t mvmA, void *A, prcd_t evalB, prcd_t evaltransB, void *B, uint rows, uint cols)
 Approximate the spectral norm $\|A-B\|_2$ of the difference of two matrices $A$ and $B$. The matrix $B$ is given as a factorization and can be applied to some vector. More...
 
real norm2diff_id_pre_matrix (mvm_t mvmA, void *A, prcd_t solveB, prcd_t solvetransB, void *B, uint rows, uint cols)
 Approximate the spectral norm $\|I - B^{-1}A\|_2$. The matrix $B$ is given as a factorization and can be applied to some vector. Usually $B$ is an approximation to $A$ and therefore this functions measures the quality of some preconditioner for $A$. More...
 
void init_cg (addeval_t addeval, void *matrix, pcavector b, pavector x, pavector r, pavector p, pavector a)
 Initialize a standard conjugate gradient method to solve $A x = b$. More...
 
void step_cg (addeval_t addeval, void *matrix, pcavector b, pavector x, pavector r, pavector p, pavector a)
 One step of a standard conjugate gradient method. More...
 
real evalfunctional_cg (addeval_t addeval, void *matrix, pcavector b, pcavector x, pcavector r)
 error information of a standard conjugate gradient method More...
 
void init_pcg (addeval_t addeval, void *matrix, prcd_t prcd, void *pdata, pcavector b, pavector x, pavector r, pavector q, pavector p, pavector a)
 Initialize a preconditioned conjugate gradient method to solve $N^{1/2} A N^{1/2} \widehat{x} = N^{1/2} b$ with $x = N^{1/2} \widehat{x}$. More...
 
void step_pcg (addeval_t addeval, void *matrix, prcd_t prcd, void *pdata, pcavector b, pavector x, pavector r, pavector q, pavector p, pavector a)
 One step of a preconditioned conjugate gradient method. More...
 
void init_uzawa (prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21, pcavector b1, pcavector b2, pavector x1, pavector x2, pavector r2, pavector p2, pavector a1, pavector s2)
 Initialize the standard Uzawa iteration. More...
 
void step_uzawa (prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21, pcavector b1, pcavector b2, pavector x1, pavector x2, pavector r2, pavector p2, pavector a1, pavector s2)
 Perform one step of the standard Uzawa iteration. More...
 
void init_puzawa (prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21, prcd_t prcd, void *pdata, pcavector b1, pcavector b2, pavector x1, pavector x2, pavector r2, pavector q2, pavector p2, pavector a1, pavector s2)
 Initialize a preconditioned Uzawa iteration. More...
 
void step_puzawa (prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21, prcd_t prcd, void *pdata, pcavector b1, pcavector b2, pavector x1, pavector x2, pavector r2, pavector q2, pavector p2, pavector a1, pavector s2)
 Perform one step of the preconditioned Uzawa iteration. More...
 
void init_bicg (addeval_t addeval, addeval_t addevaltrans, void *matrix, pcavector b, pavector x, pavector r, pavector rt, pavector p, pavector pt, pavector a, pavector at)
 Initialize a biconjugate gradient method. More...
 
void step_bicg (addeval_t addeval, addeval_t addevaltrans, void *matrix, pcavector b, pavector x, pavector r, pavector rt, pavector p, pavector pt, pavector a, pavector at)
 One step a biconjugate gradient method. More...
 
void init_bicgstab (addeval_t addeval, void *matrix, pcavector b, pavector x, pavector r, pavector rt, pavector p, pavector a, pavector as)
 Initialize a stabilized biconjugate gradient method. More...
 
void step_bicgstab (addeval_t addeval, void *matrix, pcavector b, pavector x, pavector r, pavector rt, pavector p, pavector a, pavector as)
 One step a stabilized biconjugate gradient method. More...
 
void init_gmres (addeval_t addeval, void *matrix, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
 Initialize GMRES. More...
 
void step_gmres (addeval_t addeval, void *matrix, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
 One step of the GMRES method. More...
 
void finish_gmres (addeval_t addeval, void *matrix, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
 Completes or restarts the GMRES method. More...
 
real residualnorm_gmres (pcavector rhat, uint k)
 Returns norm of current residual vector. More...
 
void init_pgmres (addeval_t addeval, void *matrix, prcd_t prcd, void *pdata, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
 Initialize preconditioned GMRES. More...
 
void step_pgmres (addeval_t addeval, void *matrix, prcd_t prcd, void *pdata, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
 One step of the preconditioned GMRES method. More...
 
void finish_pgmres (addeval_t addeval, void *matrix, prcd_t prcd, void *pdata, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
 Completes or restarts the preconditioned GMRES method. More...
 
real residualnorm_pgmres (pcavector rhat, uint k)
 Returns norm of current preconditioned residual vector. More...
 

Detailed Description

Author
Steffen Börm