bslv.h File Reference
#include "bsr.h"
#include "prec.h"
#include <stdbool.h>
Include dependency graph for bslv.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  bslv_memory
 Linear solver memory. More...
 

Typedefs

typedef struct bslv_memory bslv_memory
 Linear solver memory. More...
 

Functions

bslv_memorybslv_alloc ()
 Create empty solver memory object. More...
 
void bslv_free (bslv_memory *mem)
 Delete solver memroy object. More...
 
void bslv_info (bslv_memory *mem, int count)
 Display convergence information. More...
 
void bslv_init (bslv_memory *mem, double tol, int max_iter, bsr_matrix const *A, bool use_dilu)
 Initialize solver memory object. More...
 
int bslv_pbicgstab3m (bslv_memory *mem, bsr_matrix *A, const double *b, double *x)
 Preconditioned bicgstab in mixed-precision. More...
 
int bslv_pbicgstab3d (bslv_memory *mem, bsr_matrix *A, const double *b, double *x)
 Preconditioned bicgstab in double-precision. More...
 

Typedef Documentation

◆ bslv_memory

typedef struct bslv_memory bslv_memory

Linear solver memory.

Function Documentation

◆ bslv_alloc()

bslv_memory * bslv_alloc ( )

Create empty solver memory object.

Returns
Pointer to solver memory object.

Referenced by Dune::MixedSolver< X, M >::MixedSolver().

◆ bslv_free()

void bslv_free ( bslv_memory mem)

Delete solver memroy object.

Parameters
memPointer to solver memory object.

Referenced by Dune::MixedSolver< X, M >::~MixedSolver().

◆ bslv_info()

void bslv_info ( bslv_memory mem,
int  count 
)

Display convergence information.

Parameters
memPointer to solver memory object.
countNumber of linear iterations (returned by linear solver)

◆ bslv_init()

void bslv_init ( bslv_memory mem,
double  tol,
int  max_iter,
bsr_matrix const *  A,
bool  use_dilu 
)

Initialize solver memory object.

Parameters
memPointer to solver memory object.
tolLinear solver tolerance. @apram max_iter Maximum number of linear iterations.
APointer to coeffient matrix in bsr format.
use_diluSelect DILU preconditioner if true. Otherwise use ILU9.

Referenced by Dune::MixedSolver< X, M >::MixedSolver().

◆ bslv_pbicgstab3d()

int bslv_pbicgstab3d ( bslv_memory mem,
bsr_matrix A,
const double *  b,
double *  x 
)

Preconditioned bicgstab in double-precision.

Note
Preconditioner is either ILU0 or DILU based on value of mem->use_dilu
Parameters
memPointer to solver memory object.
APointer to coeffient matrix in bsr format
bPointer to right-hand side vector @apram x Pointer to solution vector
Returns
Number of linear iterations.

◆ bslv_pbicgstab3m()

int bslv_pbicgstab3m ( bslv_memory mem,
bsr_matrix A,
const double *  b,
double *  x 
)

Preconditioned bicgstab in mixed-precision.

Note
Preconditioner is either ILU0 or DILU based on value of mem->use_dilu
Parameters
memPointer to solver memory object.
APointer to coeffient matrix in bsr format
bPointer to right-hand side vector @apram x Pointer to solution vector
Returns
Number of linear iterations.

Referenced by Dune::MixedSolver< X, M >::apply().