Fabulous
fabulous_basic.h File Reference
#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 Documentation

◆ fabulous_arithmetic

◆ fabulous_callback_t

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.

Parameters
user_datasame pointer that was passed to fabulous_create()
NRHSnumber 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]Xeither NULL or a matrix of size ( 'dim' x NRHS ) corresponding to solution vector
ldxleading dimension of X
[in]Reither NULL or a matrix of size ( 'dim' x NRHS ) corresponding to residual vector
ldrleading dimension of R
handlecan be used by the user as a parameter to fabulous_set_iteration_user_data()
Note
this handle may be different from the solver handle but the user may not call fabulous_destroy() on this handle because it is internal to the library
dim is the parameter that was passed to fabulous_create()

◆ fabulous_dot_t

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.

Parameters
user_datasame pointer that was passed to fabulous_create()
Mnumber of vectors in A
Nnumber of vectors in B
[in]Amatrix of size ( 'dim' x M ), to be read by the user
ldaleading dimension of A
[in]Bmatrix of size ( 'dim' x N ), to be read by the user
ldbleading dimension of B
[out]Cmatrix of size ( M x N ) to be filled with the result by the user
ldcleading dimension of C
Returns
numbers of flops performed by the operation for (optionnal) perfomance analysis
Note
dim is the parameter that was passed to fabulous_create()

◆ fabulous_mvp_t

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

Parameters
user_datasame pointer that was passed to fabulous_create()
Nnumber of vector in the matrices X and B.
[in]Xmatrix of size 'dim' x N to be read. Values are stored in Column Major layout (as in Fortran.)
ldxleading dimension of X
[out]Bmatrix of size 'dim' x N to be written by the user. Values are stored in Column Major layout.
ldbleading dimension of B
Returns
numbers of flops performed by the operation for (optionnal) perfomance analysis
Note
dim is the parameter that was passed to fabulous_create()

◆ fabulous_orthoscheme

◆ fabulous_orthotype

◆ fabulous_rightprecond_t

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

Parameters
user_datasame pointer that was passed to fabulous_create()
Nnumber of vector inside the matrices B and X
[in]XMatrix of size 'dim' x N to be read. Values are stored in Column Major layout (as in Fortran)
ldxleading dimension of X
[out]Bmatrix of size 'dim' x N to be written by the user. Values are stored in Column Major layout.
ldbleading dimension of B
Returns
numbers of flops performed by the operation for (optionnal) perfomance analysis
Note
dim is the parameter that was passed to fabulous_create()

Enumeration Type Documentation

◆ fabulous_arithmetic

The different arithmetics available.

Enumerator
FABULOUS_REAL_FLOAT 

correspond to 'float' type

FABULOUS_REAL_DOUBLE 

correspond to 'double' type

FABULOUS_COMPLEX_FLOAT 

correspond to 'float complex' type

FABULOUS_COMPLEX_DOUBLE 

correspond to 'double complex' type

◆ fabulous_orthoscheme

The different available orthogonalization schemas.

Enumerator
FABULOUS_MGS 

Modified Gram-Schmidt

FABULOUS_IMGS 

Iterated Modified Gram-Schmidt

FABULOUS_CGS 

Classical Gram-Schmidt

FABULOUS_ICGS 

Iteratied Classical Gram-Schmidt

◆ fabulous_orthotype

The different available orthogonalization types.

Enumerator
FABULOUS_RUHE 

RUHE variant (vector-by-vector)

FABULOUS_BLOCK 

BLOCK variant (block-by-block)

Function Documentation

◆ fabulous_create()

fabulous_handle fabulous_create ( fabulous_arithmetic  ari,
int  dim,
void *  user_data 
)

Create a library handle.

Parameters
arithe arithmetic used (float/double) (real/complex)
dimLocal 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_dataa pointer that will be provided to each callback. The user may use this pointer to store extra informations needed by the callbacks (see examples.)
Returns
a fabulous_handle opaque object, it is needed to make calls to other functions of the library or NULL if an error occured

◆ fabulous_destroy()

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.

Parameters
handlethe fabulous opaque handle to be destroyed

◆ fabulous_error_type_str()

const char* fabulous_error_type_str ( int  errcode)

get the error number type description as a string

Returns
the string as a c-style NUL terminated string The user must not free the value returned by this function as it is stored statically. The returned string must not be modified by the user but it may be modified by subsequent call to fabulous_error_type_str()

◆ fabulous_get_logs()

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}.
Parameters
[out]sizefilled with the size of the returned array
handlethe fabulous opaque handle
Returns
the array containing the convergence history as described above.
Warning
The returned array is valid until fabulous_destroy() is called on the handle or another solve call is performed. If you want to access the data after calling fabulous_destroy() or a solve function, you need to copy the data

◆ fabulous_get_Uk_const_last_solve()

int fabulous_get_Uk_const_last_solve ( fabulous_handle  handle)

Get the constant property of the latest deflation space used for GCRO.

Returns
true if the latest deflation space used for GCRO was constant, false otherwise

◆ fabulous_last_error_str()

const char* fabulous_last_error_str ( )

get the error description as a string

Returns
the string as a c-style NUL terminated string The user must not free the value returned by this function as it is stored statically. The returned string must not be modified by the user but it may be modified by subsequent call to fabulous_last_error_str()

◆ fabulous_print_log_file()

int fabulous_print_log_file ( void *  file,
const char *  log_id,
fabulous_handle  handle 
)

Print the convergence history with advanced details.

Parameters
filemust be a valid pointer of type (FILE*) returned by fopen or similar
log_ida string to identify the log
handlethe fabulous opaque handle
Returns
0 on success or a negative value to indicate an error

◆ fabulous_print_log_filename()

int fabulous_print_log_filename ( const char *  filename,
const char *  log_id,
fabulous_handle  handle 
)

Print the convergence history with advanced details.

Parameters
filenamethe output filename
log_ida string to identify the log
handlethe fabulous opaque handle
Returns
0 on success or a negative value to indicate an error

◆ fabulous_set_advanced_parameters()

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

Parameters
max_kept_directionmaximum 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_residualIn GMRES algorithms, compute X and R at each iteration such that the user can access them in CallBack
logger_user_data_sizenumber of slots for user data when calling fabulous_set_iteration_user_data()
quietboolean value to limit fabulous verbosity
handlethe fabulous library opaque handle
Returns
0 if the call was successful or a negative value to indicate an error
Warning
Using real_residual parameter set to true may slow down BGMRES a lot! it is designed to be used for debugging purposes

◆ fabulous_set_callback()

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.

Note
Setting a user callback is NOT mandatory for the fabulous library to work. If this function is never called, this is not a problem.
Parameters
callbackuser callback (function pointer)
handlethe fabulous opaque handle

◆ fabulous_set_dot_product()

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.

Parameters
user_dotuser dot product callback (function pointer)
handlethe fabulous opaque handle

◆ fabulous_set_iteration_user_data()

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

Parameters
idNumber 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.
fpfloating point user data
integinteger user data
handlethe last parameter of the callback of type fabulous_callback_t
Returns
0 on success or a negative value to indicate an error

◆ fabulous_set_mvp()

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.

Parameters
user_mvpmatrix vector product callback (function pointer), refer to documentation to get the prototype
handlethe fabulous opaque handle
Note
MVP stand for Matrix Vector Product

◆ fabulous_set_ortho_process()

int fabulous_set_ortho_process ( fabulous_orthoscheme  scheme,
fabulous_orthotype  type,
int  nb_iter,
fabulous_handle  handle 
)

Set the orthogonalization schema.

Parameters
schemeone 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
typeone of the following:
FABULOUS_RUHE  -> RUHE variant (vector-by-vector)
FABULOUS_BLOCK -> BLOCK variant (block-by-block)
nb_iternumber of iteration for Iterated Schemas (IMGS and ICGS); Must be positive integer >= 2
handlethe 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.

Note
if no call to this function is made, FABULOUS_MGS will be used as the default value for the orthoproc parameter
Returns
0 if the call was successful or a negative value to indicate an error

◆ fabulous_set_output_colored()

void fabulous_set_output_colored ( int  enable_color)

enable or disable color in fabulous output

Parameters
enable_colorboolean value to indicate whether colors should be enabled
Note
the effect is global and does not depend on any handle

◆ fabulous_set_parameters()

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

Parameters
max_mvpMaximum number of matrix vector product to be done
max_krylov_space_sizeMaximum size allowed for Krylov Search Space.
[in]tolerancean array of scalar containing the tolerance threshold used to consider convergence is reached (normalized backward error)
  • the type of value in the array must match the correct type to represent the module of the Arithetic used: for instance if the handle was created with FABULOUS_REAL_DOUBLE or FABULOUS_COMPLEX_DOUBLE, then tolerance must be an array of double
  • otherwise
    tolerance must be an array of float
nb_tolerancesize of the tolerance array
  • if nb_tolerance == 1
    then tolerance[0] will be used as the threshold for all vectors.
  • otherwise
    nb_tolerance must match the parameter nrhs (number of right hand sides) passed in the fabulous_solve_XX() family of function
  • if nb_tolerance is neither equal to 1 nor nrhs
    this is an error and the behaviour is undefined.
handlethe fabulous library opaque handle
Returns
0 if the call was successful or a negative value to indicate an error

◆ fabulous_set_rightprecond()

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.

Parameters
user_rpccallback(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.
handlethe fabulous library opaque handle
Note
If this function is never called, it is assumed no preconditioner is used.

◆ fabulous_solve()

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.

Parameters
nrhsnumber of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_DR()

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)

Parameters
nrhs(p) number of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
knumber of eigen pair used in Deflated Restarting
targettargeted 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()
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_GCR()

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.

Parameters
nrhsnumber of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_GCRO()

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.

Parameters
[in]BRight Hand Side (Column major layout)
[in,out]Xinitial solution guess (Column major layout)
[in,out]Ukdeflation space (Column major layout)
[in,out]ldudeflation space leading dimension
[in,out]kudeflation space true size
knumber of eigen pair used in Deflated Restarting (size of Uk)
targettargeted 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_constIndicate 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)
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_IB()

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)

Parameters
nrhsnumber of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout).
ldxleading dimension of X
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_IBDR()

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)

Parameters
nrhs(p) number of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
knumber of eigen pair used in Deflated Restarting
targettargeted 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()
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_QR()

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.

Parameters
nrhsnumber of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_QRDR()

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)

Parameters
nrhs(p) number of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
knumber of eigen pair used in Deflated Restarting
targettargeted 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()
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_QRIB()

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)

Parameters
nrhs(p) number of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

◆ fabulous_solve_QRIBDR()

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)

Parameters
nrhs(p) number of right hand side
[in]BRight Hand Side (Column major layout)
ldbleading dimension of B
[in,out]Xinitial solution guess (Column major layout)
ldxleading dimension of X
knumber of eigen pair used in Deflated Restarting
targettargeted 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()
handlethe fabulous opaque handle
Returns
number of matrix vector product performed or a negative value to indicate an error

Variable Documentation

◆ fabulous_handle

FABULOUS_BEGIN_C_DECL typedef void* fabulous_handle

opaque handle for the fabulous library

Examples
example_api.c, and example_api2.c.