template<class T = real_t>
class gismo::gsGradientMethod< T >
The gradient method.
Iterates x^{new} = x - damping * precond * ( mat * x - rhs ),
where by default damping = 1 and precond = identity.
|
bool | adaptiveStepSize () const |
| Returns true iff adaptive step sizes are activated.
|
|
virtual std::string | detail () const |
| Prints the object as a string with extended details.
|
|
T | error () const |
| The relative residual error of the current iterate.
|
|
virtual void | finalizeIteration (VectorType &) |
| Some post-processing might be required.
|
|
template<typename OperatorType > |
| gsGradientMethod (const OperatorType &mat, const LinOpPtr &precond, T step_size) |
| Constructor using a matrix (operator) and optionally a preconditionner.
|
|
template<typename OperatorType > |
| gsGradientMethod (const OperatorType &mat, const LinOpPtr &precond=LinOpPtr()) |
| Constructor using a matrix (operator) and optionally a preconditionner.
|
|
bool | initIteration (const VectorType &rhs, VectorType &x) |
| Init the iteration.
|
|
index_t | iterations () const |
| The number of iterations needed to reach the error criteria.
|
|
LinOpPtr | preconditioner () const |
| Get the preconditioner.
|
|
std::ostream & | print (std::ostream &os) const |
| Prints the object as a string.
|
|
void | setAdaptiveStepSize () |
| Activate adaptive step size. Then in each step, the step size is chosen such that the norm of the residual is minimized.
|
|
void | setMaxIterations (index_t max_iters) |
| Set the maximum number of iterations (default: 1000)
|
|
gsGradientMethod & | setOptions (const gsOptionList &opt) |
| Set the options based on a gsOptionList.
|
|
void | setPreconditioner (const LinOpPtr &precond) |
| Set the preconditionner.
|
|
void | setStepSize (T step_size) |
| Set the step size.
|
|
void | setTolerance (T tol) |
| Set the tolerance for the error criteria on the relative residual error (default: 1e-10)
|
|
index_t | size () const |
| Returns the size of the linear system.
|
|
void | solve (const VectorType &rhs, VectorType &x) |
| Solves the linear system and stores the solution in x.
|
|
void | solveDetailed (const VectorType &rhs, VectorType &x, VectorType &error_history) |
| Solves the linear system and stores the solution in x and records the error histroy.
|
|
bool | step (VectorType &x) |
| Perform one step, requires initIteration.
|
|
T | stepSize () const |
| Returns the chosen step size.
|
|
T | tolerance () const |
| The chosen tolerance for the error criteria on the relative residual error.
|
|
LinOpPtr | underlying () const |
| Get the underlying matrix/operator to be solved for.
|
|
template<class T = real_t>
template<typename OperatorType >
static uPtr make |
( |
const OperatorType & |
mat, |
|
|
const LinOpPtr & |
precond = LinOpPtr() |
|
) |
| |
|
inlinestatic |
Constructor using a matrix (operator) and optionally a preconditionner.
- Parameters
-
mat | The operator to be solved for, see gsIterativeSolver for details |
precond | The preconditioner, defaulted to the identity |
In each step, the step width is chosen such that the residual is minimized (adaptive step sizes)