void init_gmres(addeval_t addeval, void *matrix, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
Initialize GMRES.
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.
unsigned uint
Unsigned integer type.
Definition: settings.h:70
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.
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.
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.
double _Complex field
Field type.
Definition: settings.h:171
void(* mvm_t)(field alpha, bool trans, void *matrix, pcavector x, pavector y)
Matrix callback.
Definition: krylov.h:60
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.
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.
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.
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.
real residualnorm_gmres(pcavector rhat, uint k)
Returns norm of current residual vector.
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.
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 of the difference of two matrices and . The matrix is given as a fac...
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.
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 with .
#define HEADER_PREFIX
Prefix for function declarations.
Definition: settings.h:43
double real
real floating point type.
Definition: settings.h:97
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 .
real evalfunctional_cg(addeval_t addeval, void *matrix, pcavector b, pcavector x, pcavector r)
error information of a standard conjugate gradient method
real norm2diff_matrix(mvm_t mvmA, void *A, mvm_t mvmB, void *B, uint rows, uint cols)
Approximate the spectral norm of the difference of two matrices and .
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.
real residualnorm_pgmres(pcavector rhat, uint k)
Returns norm of current preconditioned residual vector.
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 . The matrix is given as a factorization and can be applied to some ve...
void(* prcd_t)(void *pdata, pavector r)
Preconditioner callback.
Definition: krylov.h:75
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.
void(* addeval_t)(field alpha, void *matrix, pcavector x, pavector y)
Matrix callback.
Definition: krylov.h:41
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.
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.
Representation of a matrix as an array in column-major order.
Definition: amatrix.h:43
real norm2_matrix(mvm_t mvm, void *A, uint rows, uint cols)
Approximate the spectral norm of a matrix .
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.