Fabulous
|
#include <stdint.h>
Go to the source code of this file.
Typedefs | |
typedef enum fabulous_arithmetic | fabulous_arithmetic |
typedef enum fabulous_orthoscheme | fabulous_orthoscheme |
typedef enum fabulous_orthotype | fabulous_orthotype |
typedef int64_t(* | fabulous_mvp_t) (void *user_data, int N, const void *alpha, const void *X, int ldx, const void *beta, void *B, int ldb) |
User defined function to compute product between user's Matrix (A) and Matrix given in argument (X) More... | |
typedef int64_t(* | fabulous_rightprecond_t) (void *user_data, int N, const void *X, int ldx, void *B, int ldb) |
User defined function to compute product between the user's right-preconditioner M and a matrix X. More... | |
typedef int64_t(* | fabulous_dot_t) (void *user_data, int M, int N, const void *A, int lda, const void *B, int ldb, void *C, int ldc) |
Function to perform multiple dot products between two sets of vectors C = A^{H} * B. More... | |
typedef void(* | fabulous_callback_t) (void *user_data, int iter, int NRHS, const void *X, int ldx, const void *R, int ldr, fabulous_handle handle) |
Function to be called once at the end of each iteration. More... | |
Enumerations | |
enum | fabulous_arithmetic { FABULOUS_REAL_FLOAT = 0, FABULOUS_REAL_DOUBLE = 1, FABULOUS_COMPLEX_FLOAT = 2, FABULOUS_COMPLEX_DOUBLE = 3 } |
The different arithmetics available. More... | |
enum | fabulous_orthoscheme { FABULOUS_MGS = 0, FABULOUS_IMGS = 1, FABULOUS_CGS = 2, FABULOUS_ICGS = 3 } |
The different available orthogonalization schemas. More... | |
enum | fabulous_orthotype { FABULOUS_RUHE = 0, FABULOUS_BLOCK = 1 } |
The different available orthogonalization types. More... | |
Functions | |
fabulous_handle | fabulous_create (fabulous_arithmetic ari, int dim, void *user_data) |
Create a library handle. More... | |
void | fabulous_destroy (fabulous_handle handle) |
free internal ressources held by the opaque fabulous_handle. After calling this function on a handle, no other library function may be called. More... | |
void | fabulous_set_mvp (fabulous_mvp_t user_mvp, fabulous_handle handle) |
Set the callback to compute product between user's Matrix and a block of vector. More... | |
void | fabulous_set_rightprecond (fabulous_rightprecond_t user_rpc, fabulous_handle handle) |
Set the callback to compute the product between user's right-preconditionner and a block of vector. More... | |
int | fabulous_set_parameters (int max_mvp, int max_krylov_space_size, void *tolerance, int nb_tolerance, fabulous_handle handle) |
set the solver's parameters More... | |
int | fabulous_set_advanced_parameters (int max_kept_direction, int real_residual, int logger_user_data_size, int quiet, fabulous_handle handle) |
set the solver's advanced parameters More... | |
int | fabulous_set_ortho_process (fabulous_orthoscheme scheme, fabulous_orthotype type, int nb_iter, fabulous_handle handle) |
Set the orthogonalization schema. More... | |
void | fabulous_set_dot_product (fabulous_dot_t user_dot, fabulous_handle handle) |
Set the user dot product. More... | |
void | fabulous_set_callback (fabulous_callback_t callback, fabulous_handle handle) |
Set the user callback. More... | |
int | fabulous_solve (int nrhs, void *B, int ldb, void *X, int ldx, fabulous_handle handle) |
Solve the system AX=B starting with X as initial solution guess using Block General Minimum Residual (BGMRES) algorithm. More... | |
int | fabulous_solve_GCR (int nrhs, void *B, int ldb, void *X, int ldx, fabulous_handle handle) |
Solve the system AX=B starting with X as initial solution guess using Block General Conjugate Residual (BGCR) algorithm. More... | |
int | fabulous_solve_GCRO (int nrhs, void *B, int ldb, void *X, int ldx, void **Uk, int *ldu, int *ku, int k, void *target, int Uk_const, fabulous_handle handle) |
Solve the system AX=B starting with X as initial solution guess using Block General Conjugate Residual O... (BGCRO) algorithm, starting with an empty deflation space Uk. More... | |
int | fabulous_get_Uk_const_last_solve (fabulous_handle handle) |
Get the constant property of the latest deflation space used for GCRO. More... | |
int | fabulous_solve_IB (int nrhs, void *B, int ldb, void *X, int ldx, fabulous_handle handle) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with detection of Inexact Breakdowns (IB-BGMRES) More... | |
int | fabulous_solve_DR (int nrhs, void *B, int ldb, void *X, int ldx, int k, void *target, fabulous_handle handle) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Deflated restarting (BGMRES-DR) More... | |
int | fabulous_solve_IBDR (int nrhs, void *B, int ldb, void *X, int ldx, int k, void *target, fabulous_handle handle) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Inexact Breakdown and Deflated restarting (IB-BGMRES-DR) More... | |
int | fabulous_solve_QR (int nrhs, void *B, int ldb, void *X, int ldx, fabulous_handle handle) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm (BGMRES) with Incremental QR factorization. More... | |
int | fabulous_solve_QRIB (int nrhs, void *B, int ldb, void *X, int ldx, fabulous_handle handle) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Inexact Breakdown and Incremental QR factorization (QR-IB-BGMRES) More... | |
int | fabulous_solve_QRDR (int nrhs, void *B, int ldb, void *X, int ldx, int k, void *target, fabulous_handle handle) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Deflated restarting and Incremental QR factorization (BGMRES-DR) More... | |
int | fabulous_solve_QRIBDR (int nrhs, void *B, int ldb, void *X, int ldx, int k, void *target, fabulous_handle handle) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Deflated restarting, Inexact Breakdown and Incremental QR factorization (QR-IB-BGMRES-DR) More... | |
const double * | fabulous_get_logs (int *size, fabulous_handle handle) |
Get the convergence history. More... | |
int | fabulous_set_iteration_user_data (unsigned id, double fp, int64_t integ, fabulous_handle handle) |
Allow the user to add its own data into fabulous internal logger. More... | |
int | fabulous_print_log_file (void *file, const char *log_id, fabulous_handle handle) |
Print the convergence history with advanced details. More... | |
int | fabulous_print_log_filename (const char *filename, const char *log_id, fabulous_handle handle) |
Print the convergence history with advanced details. More... | |
const char * | fabulous_last_error_str () |
get the error description as a string More... | |
const char * | fabulous_error_type_str (int errcode) |
get the error number type description as a string More... | |
void | fabulous_set_output_colored (int enable_color) |
enable or disable color in fabulous output More... | |
Variables | |
FABULOUS_BEGIN_C_DECL typedef void * | fabulous_handle |
opaque handle for the fabulous library More... | |
typedef enum fabulous_arithmetic fabulous_arithmetic |
typedef void(* fabulous_callback_t) (void *user_data, int iter, int NRHS, const void *X, int ldx, const void *R, int ldr, fabulous_handle handle) |
Function to be called once at the end of each iteration.
user_data | same pointer that was passed to fabulous_create() | |
NRHS | number of vectors in X and R In some case it may be 0 (For instance, in GMRES, the actual residual and solution are not computed on each iteration by default, because there no short term recurrence to compute them. This behaviour can be changed, see fabulous_set_advanced_parameters() (parameter 'real_residual')) | |
[in] | X | either NULL or a matrix of size ( 'dim' x NRHS ) corresponding to solution vector |
ldx | leading dimension of X | |
[in] | R | either NULL or a matrix of size ( 'dim' x NRHS ) corresponding to residual vector |
ldr | leading dimension of R | |
handle | can be used by the user as a parameter to fabulous_set_iteration_user_data() |
typedef int64_t(* fabulous_dot_t) (void *user_data, int M, int N, const void *A, int lda, const void *B, int ldb, void *C, int ldc) |
Function to perform multiple dot products between two sets of vectors C = A^{H} * B.
user_data | same pointer that was passed to fabulous_create() | |
M | number of vectors in A | |
N | number of vectors in B | |
[in] | A | matrix of size ( 'dim' x M ), to be read by the user |
lda | leading dimension of A | |
[in] | B | matrix of size ( 'dim' x N ), to be read by the user |
ldb | leading dimension of B | |
[out] | C | matrix of size ( M x N ) to be filled with the result by the user |
ldc | leading dimension of C |
typedef int64_t(* fabulous_mvp_t) (void *user_data, int N, const void *alpha, const void *X, int ldx, const void *beta, void *B, int ldb) |
User defined function to compute product between user's Matrix (A) and Matrix given in argument (X)
Operation to be performed: B := A * X
user_data | same pointer that was passed to fabulous_create() | |
N | number of vector in the matrices X and B. | |
[in] | X | matrix of size 'dim' x N to be read. Values are stored in Column Major layout (as in Fortran.) |
ldx | leading dimension of X | |
[out] | B | matrix of size 'dim' x N to be written by the user. Values are stored in Column Major layout. |
ldb | leading dimension of B |
typedef enum fabulous_orthoscheme fabulous_orthoscheme |
typedef enum fabulous_orthotype fabulous_orthotype |
typedef int64_t(* fabulous_rightprecond_t) (void *user_data, int N, const void *X, int ldx, void *B, int ldb) |
User defined function to compute product between the user's right-preconditioner M and a matrix X.
Operation to be performed: B := M*X
user_data | same pointer that was passed to fabulous_create() | |
N | number of vector inside the matrices B and X | |
[in] | X | Matrix of size 'dim' x N to be read. Values are stored in Column Major layout (as in Fortran) |
ldx | leading dimension of X | |
[out] | B | matrix of size 'dim' x N to be written by the user. Values are stored in Column Major layout. |
ldb | leading dimension of B |
enum fabulous_arithmetic |
enum fabulous_orthoscheme |
enum fabulous_orthotype |
fabulous_handle fabulous_create | ( | fabulous_arithmetic | ari, |
int | dim, | ||
void * | user_data | ||
) |
Create a library handle.
ari | the arithmetic used (float/double) (real/complex) |
dim | Local dimension of the problem to be solved. This parameter is very important as it will be the assumed number of lines for all input matrices in the user callbacks. It is possible that 'dim' does not match the full size of the matrix, and the library can work in a distributed fashion. To achieve the distributed version of the library the user must implements the callback to perform the necessary communications. |
user_data | a pointer that will be provided to each callback. The user may use this pointer to store extra informations needed by the callbacks (see examples.) |
void fabulous_destroy | ( | fabulous_handle | handle | ) |
free internal ressources held by the opaque fabulous_handle. After calling this function on a handle, no other library function may be called.
handle | the fabulous opaque handle to be destroyed |
const char* fabulous_error_type_str | ( | int | errcode | ) |
get the error number type description as a string
const double* fabulous_get_logs | ( | int * | size, |
fabulous_handle | handle | ||
) |
Get the convergence history.
This function return the convergence history (in terms of backward error) as an array containing 3 double for each iteration step:
{MIN over the vectors; MAX over the vectors; TIMESTAMP}.
[out] | size | filled with the size of the returned array |
handle | the fabulous opaque handle |
int fabulous_get_Uk_const_last_solve | ( | fabulous_handle | handle | ) |
Get the constant property of the latest deflation space used for GCRO.
const char* fabulous_last_error_str | ( | ) |
get the error description as a string
int fabulous_print_log_file | ( | void * | file, |
const char * | log_id, | ||
fabulous_handle | handle | ||
) |
Print the convergence history with advanced details.
file | must be a valid pointer of type (FILE*) returned by fopen or similar |
log_id | a string to identify the log |
handle | the fabulous opaque handle |
int fabulous_print_log_filename | ( | const char * | filename, |
const char * | log_id, | ||
fabulous_handle | handle | ||
) |
Print the convergence history with advanced details.
filename | the output filename |
log_id | a string to identify the log |
handle | the fabulous opaque handle |
int fabulous_set_advanced_parameters | ( | int | max_kept_direction, |
int | real_residual, | ||
int | logger_user_data_size, | ||
int | quiet, | ||
fabulous_handle | handle | ||
) |
set the solver's advanced parameters
max_kept_direction | maximum number of kept direction per iteration This parameter is meaningful only for IB variants, otherwise it is not taken into account. You can set this parameter to -1 if you want it to match the number of right hand sides |
real_residual | In GMRES algorithms, compute X and R at each iteration such that the user can access them in CallBack |
logger_user_data_size | number of slots for user data when calling fabulous_set_iteration_user_data() |
quiet | boolean value to limit fabulous verbosity |
handle | the fabulous library opaque handle |
void fabulous_set_callback | ( | fabulous_callback_t | callback, |
fabulous_handle | handle | ||
) |
Set the user callback.
Set the user callback to be called once in each iteration.
callback | user callback (function pointer) |
handle | the fabulous opaque handle |
void fabulous_set_dot_product | ( | fabulous_dot_t | user_dot, |
fabulous_handle | handle | ||
) |
Set the user dot product.
It is needed for computation of norms and orthogonalization coefficents in the orthogonalization schemas.
user_dot | user dot product callback (function pointer) |
handle | the fabulous opaque handle |
int fabulous_set_iteration_user_data | ( | unsigned | id, |
double | fp, | ||
int64_t | integ, | ||
fabulous_handle | handle | ||
) |
Allow the user to add its own data into fabulous internal logger.
This function can be called only inside the user callback of type fabulous_callback_t which was set with the fabulous_set_user_callback() function
id | Number identifying the user data. Must be an integer greater or equal than 0 and less than the maximum number of user data slots availables (The default is 4). This maximum number of slots can be changed with the fabulous_set_advanced_parameters() function. |
fp | floating point user data |
integ | integer user data |
handle | the last parameter of the callback of type fabulous_callback_t |
void fabulous_set_mvp | ( | fabulous_mvp_t | user_mvp, |
fabulous_handle | handle | ||
) |
Set the callback to compute product between user's Matrix and a block of vector.
user_mvp | matrix vector product callback (function pointer), refer to documentation to get the prototype |
handle | the fabulous opaque handle |
int fabulous_set_ortho_process | ( | fabulous_orthoscheme | scheme, |
fabulous_orthotype | type, | ||
int | nb_iter, | ||
fabulous_handle | handle | ||
) |
Set the orthogonalization schema.
scheme | one of the following: FABULOUS_MGS -> Modified Gram Schmidt (default value) FABULOUS_CGS -> Classical Gram Schmidt FABULOUS_IMGS -> Iterated Modified Gram Schmidt FABULOUS_ICGS -> Iterated Classical Gram Schmidt |
type | one of the following: FABULOUS_RUHE -> RUHE variant (vector-by-vector) FABULOUS_BLOCK -> BLOCK variant (block-by-block) |
nb_iter | number of iteration for Iterated Schemas (IMGS and ICGS); Must be positive integer >= 2 |
handle | the fabulous library opaque handle |
The choice of the schema will have influence on the numerical quality of the orthogonalization coefficients. A better orthogonalization schema may improve the convergence speed in terms of number of iteration but it may also slow the speed of the algorithm if you use the library in a distributed fashion since it requires a greater amount of synchronization.
less synchronizations -------------> more synchronizations CLASSICAL Gram-Schidt -------------> MODIFIED Gram-Schmidt
The iterated version of these algorithms denote the facts than a schema may be repeated several times (the default is two time) to improve the numerical quality. Of the course the number of synchronization is also multiplied by the amount of times the schema is repeated.
void fabulous_set_output_colored | ( | int | enable_color | ) |
enable or disable color in fabulous output
enable_color | boolean value to indicate whether colors should be enabled |
int fabulous_set_parameters | ( | int | max_mvp, |
int | max_krylov_space_size, | ||
void * | tolerance, | ||
int | nb_tolerance, | ||
fabulous_handle | handle | ||
) |
set the solver's parameters
max_mvp | Maximum number of matrix vector product to be done | |
max_krylov_space_size | Maximum size allowed for Krylov Search Space. | |
[in] | tolerance | an array of scalar containing the tolerance threshold used to consider convergence is reached (normalized backward error)
|
nb_tolerance | size of the tolerance array
| |
handle | the fabulous library opaque handle |
void fabulous_set_rightprecond | ( | fabulous_rightprecond_t | user_rpc, |
fabulous_handle | handle | ||
) |
Set the callback to compute the product between user's right-preconditionner and a block of vector.
user_rpc | callback(function pointer), refer to documentation to get the prototype To disable the use of a right preconditionner, the user may pass NULL as the first argument instead of a valid function pointer. |
handle | the fabulous library opaque handle |
int fabulous_solve | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X as initial solution guess using Block General Minimum Residual (BGMRES) algorithm.
nrhs | number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
handle | the fabulous opaque handle |
int fabulous_solve_DR | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
int | k, | ||
void * | target, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Deflated restarting (BGMRES-DR)
nrhs | (p) number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
k | number of eigen pair used in Deflated Restarting | |
target | targeted eigen value: the eigen values used will be the k eigen values whose are the closest to target. The type of target must match the arithmetic chosen in fabulous_create() | |
handle | the fabulous opaque handle |
int fabulous_solve_GCR | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X as initial solution guess using Block General Conjugate Residual (BGCR) algorithm.
nrhs | number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
handle | the fabulous opaque handle |
int fabulous_solve_GCRO | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
void ** | Uk, | ||
int * | ldu, | ||
int * | ku, | ||
int | k, | ||
void * | target, | ||
int | Uk_const, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X as initial solution guess using Block General Conjugate Residual O... (BGCRO) algorithm, starting with an empty deflation space Uk.
[in] | B | Right Hand Side (Column major layout) |
[in,out] | X | initial solution guess (Column major layout) |
[in,out] | Uk | deflation space (Column major layout) |
[in,out] | ldu | deflation space leading dimension |
[in,out] | ku | deflation space true size |
k | number of eigen pair used in Deflated Restarting (size of Uk) | |
target | targeted eigen value: the eigen values used will be the k eigen values whose are the closest to target. The type of target must match the arithmetic chosen in fabulous_g_create() | |
Uk_const | Indicate if the deflation space considered a good enough approximation and thus does not need to be updated (if set to true, the deflation space will not be updated) | |
handle | the fabulous opaque handle |
int fabulous_solve_IB | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with detection of Inexact Breakdowns (IB-BGMRES)
nrhs | number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout). |
ldx | leading dimension of X | |
handle | the fabulous opaque handle |
int fabulous_solve_IBDR | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
int | k, | ||
void * | target, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Inexact Breakdown and Deflated restarting (IB-BGMRES-DR)
nrhs | (p) number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
k | number of eigen pair used in Deflated Restarting | |
target | targeted eigen value: the eigen values used will be the k eigen values whose are the closest to target. The type of target must match the arithmetic chosen in fabulous_create() | |
handle | the fabulous opaque handle |
int fabulous_solve_QR | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm (BGMRES) with Incremental QR factorization.
nrhs | number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
handle | the fabulous opaque handle |
int fabulous_solve_QRDR | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
int | k, | ||
void * | target, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Deflated restarting and Incremental QR factorization (BGMRES-DR)
nrhs | (p) number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
k | number of eigen pair used in Deflated Restarting | |
target | targeted eigen value: the eigen values used will be the k eigen values whose are the closest to target. The type of target must match the arithmetic chosen in fabulous_create() | |
handle | the fabulous opaque handle |
int fabulous_solve_QRIB | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Inexact Breakdown and Incremental QR factorization (QR-IB-BGMRES)
nrhs | (p) number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
handle | the fabulous opaque handle |
int fabulous_solve_QRIBDR | ( | int | nrhs, |
void * | B, | ||
int | ldb, | ||
void * | X, | ||
int | ldx, | ||
int | k, | ||
void * | target, | ||
fabulous_handle | handle | ||
) |
Solve the system AX=B starting with X0 using Block General Minimum Residual algorithm with Deflated restarting, Inexact Breakdown and Incremental QR factorization (QR-IB-BGMRES-DR)
nrhs | (p) number of right hand side | |
[in] | B | Right Hand Side (Column major layout) |
ldb | leading dimension of B | |
[in,out] | X | initial solution guess (Column major layout) |
ldx | leading dimension of X | |
k | number of eigen pair used in Deflated Restarting | |
target | targeted eigen value: the eigen values used will be the k eigen values whose are the closest to target. The type of target must match the arithmetic chosen in fabulous_create() | |
handle | the fabulous opaque handle |
FABULOUS_BEGIN_C_DECL typedef void* fabulous_handle |
opaque handle for the fabulous library