H2Lib
3.0
|
Iterative solvers of Krylov type. More...
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 of a matrix . More... | |
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 . 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 of the difference of two matrices and . The matrix 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 . The matrix is given as a factorization and can be applied to some vector. Usually is an approximation to and therefore this functions measures the quality of some preconditioner for . 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 . 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 with . 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... | |
Iterative solvers of Krylov type.
Matrix callback.
Used to evaluate the system matrix or its adjoint, i.e., to perform .
Functions like addeval_amatrix_avector, addevaltrans_amatrix_avector, addeval_hmatrix_avector, addevaltrans_hmatrix_avector, addeval_h2matrix_avector, addevaltrans_h2matrix_avector or addeval_sparsematrix_avector can be cast to addeval_t
.
alpha | Scaling factor . |
matrix | Matrix data describing . |
x | Source vector . |
y | Target vector . |
Matrix callback.
Used to evaluate the system matrix or its adjoint, i.e., to perform .
Compared to addeval_t, callbacks of this type allow us to evaluate both the matrix and its adjoint.
Functions like mvm_amatrix_avector, mvm_hmatrix_avector, mvm_h2matrix_avector, or mvm_sparsematrix_avector can be cast to mvm_t
.
alpha | Scaling factor . |
trans | Set if is to be used instead of . |
matrix | Matrix data describing . |
x | Source vector . |
y | Target vector . |
typedef void(* prcd_t) (void *pdata, pavector r) |
Preconditioner callback.
Used to apply a precondtioner to a vector, i.e., to perform .
Functions like lrsolve_amatrix_avector, cholsolve_amatrix_avector, ldltsolve_amatrix_avector, qrsolve_amatrix_avector, lrsolve_hmatrix_avector or cholsolve_hmatrix_avector can be cast to prcd_t
.
pdata | Preconditioner data describing . |
r | Source vector, will be overwritten by result. |
error information of a standard conjugate gradient method
addeval | Callback function name |
matrix | untyped pointer to matrix data (A) |
b | Right-hand side. |
x | Approximate solution. |
r | Residual |
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.
Solves the least-squares problem and performs the update .
The function calls init_gmres to reset the iteration and prepare for a restart.
addeval | Callback function representing the matrix . |
matrix | Data for the addeval callback. |
b | Right-hand side vector . |
x | Initial guess for the solution , will eventually be replaced by an improved approximation. |
rhat | Transformed residual. The absolute value of rhat[*kk] is the Euclidean norm of the residual. |
q | Next vector of the Krylov basis, constructed by Householder's elementary reflectors. |
kk | Pointer to current dimension of the Krylov space. |
qr | Representation of the Arnoldi basis and the transformed matrix . |
tau | Scaling factors of elementary reflectors, provided by qrdecomp_amatrix. |
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.
Solves the least-squares problem and performs the update .
The function calls init_pgmres to reset the iteration and prepare for a restart.
addeval | Callback function representing the matrix . |
matrix | Data for the addeval callback. |
prcd | Callback function representing the preconditioner . |
pdata | Data for the prcd callback. |
b | Right-hand side vector . |
x | Initial guess for the solution , will eventually be replaced by an improved approximation. |
rhat | Transformed residual. The absolute value of rhat[*kk] is the Euclidean norm of the preconditioned residual. |
q | Next vector of the Krylov basis, constructed by Householder's elementary reflectors. |
kk | Pointer to current dimension of the Krylov space. |
qr | Representation of the Arnoldi basis and the transformed matrix . |
tau | Scaling factors of elementary reflectors, provided by qrdecomp_amatrix. |
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.
addeval | Callback function name |
addevaltrans | Callback function name |
matrix | untyped pointer to matrix data |
b | Right-hand side. |
x | Approximate solution. |
r | Residual b-Ax |
rt | Adjoint residual |
p | Search direction |
pt | Adjoint search direction |
a | auxiliary vector |
at | auxiliary vector |
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.
addeval | Callback function name |
matrix | untyped pointer to matrix data |
b | Right-hand side. |
x | Approximate solution. |
r | Residual b-Ax |
rt | Adjoint residual |
p | Search direction |
a | auxiliary vector |
as | auxiliary vector |
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 .
The matrix has to be self-adjoint and positive definite (or at least positive semidefinite).
addeval | Callback function name |
matrix | untyped pointer to matrix data |
b | Right-hand side. |
x | Approximate solution. |
r | Residual . |
p | Search direction. |
a | Auxiliary vector. |
void init_gmres | ( | addeval_t | addeval, |
void * | matrix, | ||
pcavector | b, | ||
pavector | x, | ||
pavector | rhat, | ||
pavector | q, | ||
uint * | kk, | ||
pamatrix | qr, | ||
pavector | tau | ||
) |
Initialize GMRES.
The parameters are prepared for solving with the GMRES method.
The maximal dimension of the Krylov subspace is determined by the number of columns of qr
: for a k
-dimensional subspace, qr->cols==k+1
is required.
addeval | Callback function representing the matrix . |
matrix | Data for the addeval callback. |
b | Right-hand side vector . |
x | Initial guess for the solution , will eventually be replaced by an improved approximation. |
rhat | Transformed residual. The absolute value of rhat[*kk] is the Euclidean norm of the residual. |
q | Next vector of the Krylov basis, constructed by Householder's elementary reflectors. |
kk | Pointer to current dimension of the Krylov space. |
qr | Representation of the Arnoldi basis and the transformed matrix . |
tau | Scaling factors of elementary reflectors, provided by qrdecomp_amatrix. |
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 .
The matrix has to be self-adjoint and positive definite (or at least positive semidefinite). The matrix has to be self-adjoint and positive definite.
addeval | Callback function name |
matrix | untyped pointer to matrix data |
prcd | Callback function name |
pdata | untyped pointer to matrix data |
b | Right-hand side. |
x | Approximate solution. |
r | Residual b-Ax |
q | Preconditioned residual |
p | Search direction |
a | auxiliary vector |
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.
The parameters are prepared for solving with the GMRES method.
The maximal dimension of the Krylov subspace is determined by the number of columns of qr
: for a k
-dimensional subspace, qr->cols==k+1
is required.
addeval | Callback function representing the matrix . |
matrix | Data for the addeval callback. |
prcd | Callback function representing the preconditioner . |
pdata | Data for the prcd callback. |
b | Right-hand side vector . |
x | Initial guess for the solution , will eventually be replaced by an improved approximation. |
rhat | Transformed residual. The absolute value of rhat[*kk] is the Euclidean norm of the preconditioned residual. |
q | Next vector of the Krylov basis, constructed by Householder's elementary reflectors. |
kk | Pointer to current dimension of the Krylov space. |
qr | Representation of the Arnoldi basis and the transformed matrix . |
tau | Scaling factors of elementary reflectors, provided by qrdecomp_amatrix. |
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.
The preconditioned Uzawa iteration solves saddle point systems of the form
by applying the preconditioned conjugate gradient method to the Schur complement. The matrix has to be symmetric and positive definite. If is not surjective, will only be determined up to an element of the null space of .
solve_A11 | Callback for solving |
matrix_A11 | Data for solve_A11 callback |
mvm_A21 | Callback for evaluating or |
matrix_A21 | Data for mvm_A21 callback |
prcd | Callback function name |
pdata | untyped pointer to matrix data |
b1 | Right-hand side vector |
b2 | Right-hand side vector |
x1 | Approximate solution vector |
x2 | Approximate solution vector |
r2 | Residual of the Schur complement system, may be useful for a stopping criterion |
p2 | Search direction for Schur complement system |
q2 | Auxiliary variable of the same dimension as |
a1 | Auxiliary variable of the same dimension as |
s2 | Auxiliary variable of the same dimension as |
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.
The Uzawa iteration solves saddle point systems of the form
by applying the conjugate gradient method to the Schur complement. The matrix has to be symmetric and positive definite. If is not surjective, will only be determined up to an element of the null space of .
solve_A11 | Callback for solving |
matrix_A11 | Data for solve_A11 callback |
mvm_A21 | Callback for evaluating or |
matrix_A21 | Data for mvm_A21 callback |
b1 | Right-hand side vector |
b2 | Right-hand side vector |
x1 | Approximate solution vector |
x2 | Approximate solution vector |
r2 | Residual of the Schur complement system, may be useful for a stopping criterion |
p2 | Search direction for Schur complement system |
a1 | Auxiliary variable of the same dimension as |
s2 | Auxiliary variable of the same dimension as |
Approximate the spectral norm of a matrix .
The spectral norm is approximated by applying a few steps of the power iteration to the matrix and computing the square root of the resulting eigenvalue approximation.
mvm | Callback function for evaluating A and its adjoint. |
A | Matrix . |
rows | Number of rows of A . |
cols | Number of columns of A . |
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 vector. Usually is an approximation to and therefore this functions measures the quality of some preconditioner for .
The spectral norm is approximated by applying a few steps of the power iteration to the matrix and computing the square root of the resulting eigenvalue approximation. is usually given as a LR decomposition or cholesky factorization of .
mvmA | Callback function for evaluating A and its adjoint. |
A | Matrix . |
solveB | Callback function for solving a linear system with B . |
solvetransB | Callback function for solving a linear system with the adjoint of B . |
B | Factorized matrix . |
rows | Number of rows of either A and B . |
cols | Number of columns of either A and B . |
Approximate the spectral norm of the difference of two matrices and .
The spectral norm is approximated by applying a few steps of the power iteration to the matrix and computing the square root of the resulting eigenvalue approximation.
mvmA | Callback function for evaluating A and its adjoint. |
A | Matrix . |
mvmB | Callback function for evaluating B and its adjoint. |
B | Matrix . |
rows | Number of rows of either A and B . |
cols | Number of columns of either A and B . |
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 factorization and can be applied to some vector.
The spectral norm is approximated by applying a few steps of the power iteration to the matrix and computing the square root of the resulting eigenvalue approximation. is usually given as a LR decomposition or cholesky factorization of .
mvmA | Callback function for evaluating A and its adjoint. |
A | Matrix . |
evalB | Callback function for evaluating B . |
evaltransB | Callback function for evaluating the adjoint of B . |
B | Factorized matrix . |
rows | Number of rows of either A and B . |
cols | Number of columns of either A and B . |
Returns norm of current residual vector.
rhat | Transformed residual. The absolute value of rhat[k] is the Euclidean norm of the residual. |
k | Current dimension of the Krylov space. |
Returns norm of current preconditioned residual vector.
rhat | Transformed residual. The absolute value of rhat[k] is the Euclidean norm of the preconditioned residual. |
k | Current dimension of the Krylov space. |
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.
addeval | Callback function name |
addevaltrans | Callback function name |
matrix | untyped pointer to matrix data |
b | Right-hand side. |
x | Approximate solution. |
r | Residual b-Ax |
rt | Adjoint residual |
p | Search direction |
pt | Adjoint search direction |
a | auxiliary vector |
at | auxiliary vector |
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.
addeval | Callback function name |
matrix | untyped pointer to matrix data |
b | Right-hand side. |
x | Approximate solution. |
r | Residual b-Ax |
rt | Adjoint residual |
p | Search direction |
a | auxiliary vector |
as | auxiliary vector |
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.
addeval | Callback function name |
matrix | untyped pointer to matrix data (A) |
b | Right-hand side. |
x | Approximate solution. |
r | Residual b-Ax |
p | Search direction |
a | auxiliary vector |
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.
If *kk+1 >= qr->cols
, there is no room for the next Arnoldi basis vector and the function returns immediately. It can be restarted using finish_gmres.
Otherwise a new vector is added to the Arnoldi basis and the matrix is updated.
b
and does not update x
. The current residual can be tracked via rhat[*kk]
. Once it is sufficiently small, the improved solution x
can be obtained by using finish_gmres.addeval | Callback function representing the matrix . |
matrix | Data for the addeval callback. |
b | Right-hand side vector . |
x | Initial guess for the solution , will eventually be replaced by an improved approximation. |
rhat | Transformed residual. The absolute value of rhat[*kk] is the Euclidean norm of the residual. |
q | Next vector of the Krylov basis, constructed by Householder's elementary reflectors. |
kk | Pointer to current dimension of the Krylov space. |
qr | Representation of the Arnoldi basis and the transformed matrix . |
tau | Scaling factors of elementary reflectors, provided by qrdecomp_amatrix. |
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.
addeval | Callback function name |
matrix | untyped pointer to matrix data |
prcd | Callback function name |
pdata | untyped pointer to matrix data |
b | Right-hand side. |
x | Approximate solution. |
r | Residual b-Ax |
q | Preconditioned residual |
p | Search direction |
a | auxiliary vector |
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.
If *kk+1 >= qr->cols
, there is no room for the next Arnoldi basis vector and the function returns immediately. It can be restarted using finish_gmres.
Otherwise a new vector is added to the Arnoldi basis and the matrix is updated.
b
and does not update x
. The current preconditioned residual can be tracked via rhat[*kk]
. Once it is sufficiently small, the improved solution x
can be obtained by using finish_gmres.addeval | Callback function representing the matrix . |
matrix | Data for the addeval callback. |
prcd | Callback function representing the preconditioner . |
pdata | Data for the prcd callback. |
b | Right-hand side vector . |
x | Initial guess for the solution , will eventually be replaced by an improved approximation. |
rhat | Transformed residual. The absolute value of rhat[*kk] is the Euclidean norm of the preconditioned residual. |
q | Next vector of the Krylov basis, constructed by Householder's elementary reflectors. |
kk | Pointer to current dimension of the Krylov space. |
qr | Representation of the Arnoldi basis and the transformed matrix . |
tau | Scaling factors of elementary reflectors, provided by qrdecomp_amatrix. |
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.
solve_A11 | Callback for solving |
matrix_A11 | Data for solve_A11 callback |
mvm_A21 | Callback for evaluating or |
matrix_A21 | Data for mvm_A21 callback |
prcd | Callback function name |
pdata | untyped pointer to matrix data |
b1 | Right-hand side vector |
b2 | Right-hand side vector |
x1 | Approximate solution vector |
x2 | Approximate solution vector |
r2 | Residual of the Schur complement system, may be useful for a stopping criterion |
q2 | Auxiliary variable of the same dimension as |
p2 | Search direction for Schur complement system |
a1 | Auxiliary variable of the same dimension as |
s2 | Auxiliary variable of the same dimension as |
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.
solve_A11 | Callback for solving |
matrix_A11 | Data for solve_A11 callback |
mvm_A21 | Callback for evaluating or |
matrix_A21 | Data for mvm_A21 callback |
b1 | Right-hand side vector |
b2 | Right-hand side vector |
x1 | Approximate solution vector |
x2 | Approximate solution vector |
r2 | Residual of the Schur complement system, may be useful for a stopping criterion |
p2 | Search direction for Schur complement system |
a1 | Auxiliary variable of the same dimension as |
s2 | Auxiliary variable of the same dimension as |