Discontinuous Galerkin Library
#include "dg/algorithm.h"
Linear and nonlinear solvers

Linear \( Ax = b\) and non-linear \( f(x) = b\). More...

Collaboration diagram for Linear and nonlinear solvers:

Classes

struct  dg::AndersonAcceleration< ContainerType >
 Anderson Acceleration of Fixed Point/Richardson Iteration for the nonlinear equation \( f(x) = b\). More...
 
class  dg::BICGSTABl< ContainerType >
 Preconditioned BICGSTAB(l) method to solve \( Ax=b\). More...
 
class  dg::ChebyshevIteration< ContainerType >
 Preconditioned Chebyshev iteration for solving \( PAx=Pb\). More...
 
struct  dg::ChebyshevPreconditioner< Matrix, ContainerType >
 Chebyshev Polynomial Preconditioner \( C( A)\). More...
 
struct  dg::ModifiedChebyshevPreconditioner< Matrix, ContainerType >
 Approximate inverse Chebyshev Polynomial Preconditioner \( A^{-1} = \frac{c_0}{2} I + \sum_{k=1}^{r}c_kT_k( Z)\). More...
 
struct  dg::LeastSquaresPreconditioner< Matrix, InnerPreconditioner, ContainerType >
 Least Squares Polynomial Preconditioner \( M^{-1} s( AM^{-1})\). More...
 
class  dg::EVE< ContainerType >
 The Eigen-Value-Estimator (EVE) finds largest Eigenvalue of \( M^{-1}Ax = \lambda_\max x\). More...
 
struct  dg::DefaultSolver< ContainerType >
 PCG Solver class for solving \( (y-\alpha\hat I(t,y)) = \rho\). More...
 
class  dg::LGMRES< ContainerType >
 Functor class for the right preconditioned LGMRES method to solve \( Ax=b\). More...
 
class  dg::PCG< ContainerType >
 Preconditioned conjugate gradient method to solve \( Ax=b\). More...
 

Typedefs

template<class ContainerType >
using dg::FixedPointIteration = AndersonAcceleration< ContainerType >
 If you are looking for fixed point iteration: it is a special case of Anderson Acceleration. More...
 

Functions

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...
 

Detailed Description

Linear \( Ax = b\) and non-linear \( f(x) = b\).

Typedef Documentation

◆ FixedPointIteration

template<class ContainerType >
using dg::FixedPointIteration = typedef AndersonAcceleration<ContainerType>

If you are looking for fixed point iteration: it is a special case of Anderson Acceleration.

Function Documentation

◆ 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
UnaryOpunary function operator
Parameters
opFunction or Functor
x_minleft boundary, contains new left boundary on execution
x_maxright boundary, contains new right boundary on execution
epsaccuracy 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
NoRoot1dif no root lies between the given boundaries
std::runtime_errorif 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,
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
ainput matrix in csr format (self-adjoint in weights), must be sorted by columns
soutput matrix S
ddiagonal output
weightsweights in which the input matrix a is self-adjoint
nnzmaxmaximum number of non-zeroes in a row in s
thresholdabsolute 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
Treal type