H2Lib  3.0
krylov.h
Go to the documentation of this file.
1 /* ------------------------------------------------------------
2  * This is the file "krylov.h" of the H2Lib package.
3  * All rights reserved, Steffen Boerm 2010
4  * ------------------------------------------------------------ */
5 
10 #ifndef KRYLOV_H
11 #define KRYLOV_H
12 
13 #ifndef NDEBUG
14 #include <stdio.h>
15 #endif
16 
17 #include "settings.h"
18 #include "avector.h"
19 #include "factorizations.h"
20 
41 typedef void (*addeval_t)(field alpha, void *matrix, pcavector x, pavector y);
42 
60 typedef void (*mvm_t)(field alpha, bool trans, void *matrix, pcavector x,
61  pavector y);
62 
75 typedef void (*prcd_t)(void *pdata, pavector r);
76 
77 /* ------------------------------------------------------------
78  * Vector iteration for norm2 and norm2diff estimation
79  * ------------------------------------------------------------ */
80 
93 norm2_matrix(mvm_t mvm, void *A, uint rows, uint cols);
94 
110 norm2diff_matrix(mvm_t mvmA, void *A, mvm_t mvmB, void *B, uint rows, uint cols);
111 
131 norm2diff_pre_matrix(mvm_t mvmA, void *A, prcd_t evalB, prcd_t evaltransB,
132  void *B, uint rows, uint cols);
133 
155 norm2diff_id_pre_matrix(mvm_t mvmA, void *A, prcd_t solveB, prcd_t solvetransB,
156  void *B, uint rows, uint cols);
157 
158 /* ------------------------------------------------------------
159  * Conjugate gradient method (CG)
160  * ------------------------------------------------------------ */
161 
176 HEADER_PREFIX void
177 init_cg(addeval_t addeval, void *matrix, pcavector b, /* Right-hand side */
178  pavector x, /* Approximate solution */
179  pavector r, /* Residual b-Ax */
180  pavector p, /* Search direction */
181  pavector a);
182 
193 HEADER_PREFIX void
194 step_cg(addeval_t addeval, void *matrix, pcavector b, /* Right-hand side */
195  pavector x, /* Approximate solution */
196  pavector r, /* Residual b-Ax */
197  pavector p, /* Search direction */
198  pavector a);
199 
210 evalfunctional_cg(addeval_t addeval, void *matrix, pcavector b, pcavector x,
211  pcavector r);
212 
213 /* ------------------------------------------------------------
214  * Preconditioned conjugated gradient method (PCG)
215  * ------------------------------------------------------------ */
216 
236 HEADER_PREFIX void
237 init_pcg(addeval_t addeval, void *matrix, prcd_t prcd, void *pdata, pcavector b, /* Right-hand side */
238  pavector x, /* Approximate solution */
239  pavector r, /* Residual b-Ax */
240  pavector q, /* Preconditioned residual */
241  pavector p, /* Search direction */
242  pavector a);
243 
257 HEADER_PREFIX void
258 step_pcg(addeval_t addeval, void *matrix, prcd_t prcd, void *pdata, pcavector b, /* Right-hand side */
259  pavector x, /* Approximate solution */
260  pavector r, /* Residual b-Ax */
261  pavector q, /* Preconditioned residual */
262  pavector p, /* Search direction */
263  pavector a);
264 
265 /* ------------------------------------------------------------
266  * Uzawa method
267  * ------------------------------------------------------------ */
268 
295 HEADER_PREFIX void
296 init_uzawa(prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21,
297  pcavector b1, pcavector b2, pavector x1, pavector x2, pavector r2,
298  pavector p2, pavector a1, pavector s2);
299 
316 HEADER_PREFIX void
317 step_uzawa(prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21,
318  pcavector b1, pcavector b2, pavector x1, pavector x2, pavector r2,
319  pavector p2, pavector a1, pavector s2);
320 
351 HEADER_PREFIX void
352 init_puzawa(prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21,
353  prcd_t prcd, void *pdata, pcavector b1, pcavector b2, pavector x1,
354  pavector x2, pavector r2, pavector q2, pavector p2, pavector a1,
355  pavector s2);
356 
376 HEADER_PREFIX void
377 step_puzawa(prcd_t solve_A11, void *matrix_A11, mvm_t mvm_A21, void *matrix_A21,
378  prcd_t prcd, void *pdata, pcavector b1, pcavector b2, pavector x1,
379  pavector x2, pavector r2, pavector q2, pavector p2, pavector a1,
380  pavector s2);
381 
396 HEADER_PREFIX void
397 init_bicg(addeval_t addeval, addeval_t addevaltrans, void *matrix, pcavector b, /* Right-hand side */
398  pavector x, /* Approximate solution */
399  pavector r, /* Residual b-Ax */
400  pavector rt, /* Adjoint residual */
401  pavector p, /* Search direction */
402  pavector pt, /* Adjoint search direction */
403  pavector a, pavector at);
404 
419 HEADER_PREFIX void
420 step_bicg(addeval_t addeval, addeval_t addevaltrans, void *matrix, pcavector b, /* Right-hand side */
421  pavector x, /* Approximate solution */
422  pavector r, /* Residual b-Ax */
423  pavector rt, /* Adjoint residual */
424  pavector p, /* Search direction */
425  pavector pt, /* Adjoint search direction */
426  pavector a, pavector at);
427 
428 /* ------------------------------------------------------------
429  * Stabilized biconjugated gradient method (BiCGstab)
430  * ------------------------------------------------------------ */
431 
444 HEADER_PREFIX void
445 init_bicgstab(addeval_t addeval, void *matrix, pcavector b, /* Right-hand side */
446  pavector x, /* Approximate solution */
447  pavector r, /* Residual b-Ax */
448  pavector rt, /* Adjoint residual */
449  pavector p, /* Search direction */
450  pavector a, pavector as);
451 
464 HEADER_PREFIX void
465 step_bicgstab(addeval_t addeval, void *matrix, pcavector b, /* Right-hand side */
466  pavector x, /* Approximate solution */
467  pavector r, /* Residual b-Ax */
468  pavector rt, /* Adjoint residual */
469  pavector p, /* Search direction */
470  pavector a, pavector as);
471 
472 /* ------------------------------------------------------------
473  * Generalized minimal residual method (GMRES)
474  * ------------------------------------------------------------ */
475 
500 HEADER_PREFIX void
501 init_gmres(addeval_t addeval, void *matrix, pcavector b, pavector x,
502  pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau);
503 
533 HEADER_PREFIX void
534 step_gmres(addeval_t addeval, void *matrix, pcavector b, pavector x,
535  pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau);
536 
560 HEADER_PREFIX void
561 finish_gmres(addeval_t addeval, void *matrix, pcavector b, pavector x,
562  pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau);
563 
572 
573 /* ------------------------------------------------------------
574  * Preconditioned generalized minimal residual method (GMRES)
575  * ------------------------------------------------------------ */
576 
604 HEADER_PREFIX void
605 init_pgmres(addeval_t addeval, void *matrix, prcd_t prcd, void *pdata,
606  pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr,
607  pavector tau);
608 
641 HEADER_PREFIX void
642 step_pgmres(addeval_t addeval, void *matrix, prcd_t prcd, void *pdata,
643  pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr,
644  pavector tau);
645 
672 HEADER_PREFIX void
673 finish_pgmres(addeval_t addeval, void *matrix, prcd_t prcd, void *pdata,
674  pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr,
675  pavector tau);
676 
686 
689 #endif
void init_gmres(addeval_t addeval, void *matrix, pcavector b, pavector x, pavector rhat, pavector q, uint *kk, pamatrix qr, pavector tau)
Initialize GMRES.
Definition: avector.h:39
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.