Linear \( Ax = b\) and non-linear \( f(x) = b\).
More...
|
template<typename UnaryOp > |
int | dg::bisection1d (UnaryOp &op, double &x_min, double &x_max, const double eps) |
| Find a root of a 1d function \( f(x) = 0\). More...
|
|
template<class T > |
void | dg::create::sainv_precond (const cusp::csr_matrix< int, T, cusp::host_memory > &a, cusp::csr_matrix< int, T, cusp::host_memory > &s, thrust::host_vector< T > &d, const thrust::host_vector< T > &weights, unsigned nnzmax, T threshold) |
| Left looking sparse inverse preconditioner for self-adjoint positive definit matrices. More...
|
|
Linear \( Ax = b\) and non-linear \( f(x) = b\).
◆ FixedPointIteration
template<class ContainerType >
If you are looking for fixed point iteration: it is a special case of Anderson Acceleration.
◆ bisection1d()
template<typename UnaryOp >
int dg::bisection1d |
( |
UnaryOp & |
op, |
|
|
double & |
x_min, |
|
|
double & |
x_max, |
|
|
const double |
eps |
|
) |
| |
Find a root of a 1d function \( f(x) = 0\).
in given boundaries using bisection
It is assumed that a sign change occurs at the root. Function jumps closer to the root by checking the sign.
- Template Parameters
-
UnaryOp | unary function operator |
- Parameters
-
op | Function or Functor |
x_min | left boundary, contains new left boundary on execution |
x_max | right boundary, contains new right boundary on execution |
eps | accuracy of the root finding. Algorithm successful if \( |x_{\max} - x_{\min}| < \varepsilon |x_{\max}| + \varepsilon \) |
- Returns
- number of used steps to reach the desired accuracy
- Exceptions
-
NoRoot1d | if no root lies between the given boundaries |
std::runtime_error | if after 60 steps the accuracy wasn't reached |
bisection1d(f, x_min, x_max, eps);
- Note
- If the root is found exactly the x_min = x_max
◆ sainv_precond()
template<class T >
void dg::create::sainv_precond |
( |
const cusp::csr_matrix< int, T, cusp::host_memory > & |
a, |
|
|
cusp::csr_matrix< int, T, cusp::host_memory > & |
s, |
|
|
thrust::host_vector< T > & |
d, |
|
|
const thrust::host_vector< T > & |
weights, |
|
|
unsigned |
nnzmax, |
|
|
T |
threshold |
|
) |
| |
Left looking sparse inverse preconditioner for self-adjoint positive definit matrices.
\[ A^{-1} = S^T D^{-1} S W \]
where S is a lower triangular matrix, D is a diagonal matrix and W are the weights
- Note
- Implements Bertaccini and Filippone "Sparse approximate inverse preconditioners on high performance GPU platforms" (2016). S is the transpose of Z with respect to the original algorithm.
- Parameters
-
a | input matrix in csr format (self-adjoint in weights), must be sorted by columns |
s | output matrix S |
d | diagonal output |
weights | weights in which the input matrix a is self-adjoint |
nnzmax | maximum number of non-zeroes in a row in s |
threshold | absolute limit under which entries are dropped from s and entries in d are strictly greater than |
- Note
- Sparse inverse preconditioners can be applied directly like any other sparse matrix and do not need a linear system solve like in sparse LU factorization type methods
- Attention
- Use
dg::nested_iterations
if you have a geometry available, which is the faster method. Tests on the dg::Elliptic
matrix show no significant speedup that could justify the additional costs to build and apply the preconditioner. (So we did not bother to implement an MPI version)
- Template Parameters
-