DUNE PDELab (2.8)

Dune::PDELab::StationaryLinearProblemSolver< GO, LS, V > Class Template Reference

Solve linear problems using a residual formulation. More...

#include <dune/pdelab/stationary/linearproblem.hh>

Public Member Functions

 StationaryLinearProblemSolver (const GO &go, LS &ls, V &x, const ParameterTree &params)
 Construct a StationaryLinearProblemSolver for the given objects and read parameters from a ParameterTree. More...
 
 StationaryLinearProblemSolver (const GO &go, LS &ls, const ParameterTree &params)
 Construct a StationaryLinearProblemSolver for the given objects and read parameters from a ParameterTree. More...
 
void setHangingNodeModifications (bool b)
 Set whether the solver should apply the necessary transformations for calculations on hanging nodes.
 
bool hangingNodeModifications () const
 Return whether the solver performs the necessary transformations for calculations on hanging nodes.
 
void setKeepMatrix (bool b)
 Set whether the jacobian matrix should be kept across calls to apply().
 
bool keepMatrix () const
 Return whether the jacobian matrix is kept across calls to apply().
 
const Result & result () const
 Return result object.
 
void apply (V &x, bool reuse_matrix=false)
 Solve linear problem with the provided initial guess.
 
void apply (bool reuse_matrix=false)
 Solve linear problem (use initial guess that was passed at construction)
 
void discardMatrix ()
 Discard the stored Jacobian matrix.
 

Detailed Description

template<typename GO, typename LS, typename V>
class Dune::PDELab::StationaryLinearProblemSolver< GO, LS, V >

Solve linear problems using a residual formulation.

This work for matrix based and matrix free solvers. It uses a residual formulation solving \(r(u,v)=0\) instead of solving \(a(u,v)=l(v)\). In the matrix based case this means doing the following:

  1. Calculate residual: r = A z_0-b
  2. Solve: A z_u = r
  3. Update: z = z_0 - z_u
Examples
recipe-geometry-grid.cc, and recipe-linear-system-solution-pdelab.cc.

Constructor & Destructor Documentation

◆ StationaryLinearProblemSolver() [1/2]

template<typename GO , typename LS , typename V >
Dune::PDELab::StationaryLinearProblemSolver< GO, LS, V >::StationaryLinearProblemSolver ( const GO &  go,
LS &  ls,
V &  x,
const ParameterTree params 
)
inline

Construct a StationaryLinearProblemSolver for the given objects and read parameters from a ParameterTree.

This constructor reads the parameter controlling its operation from a passed-in ParameterTree instead of requiring the user to specify all of them as individual constructor parameters. Currently the following parameters are read:

Name Default Value Explanation
reduction Required relative defect reduction
min_defect 1e-99 minimum absolute defect at which to stop
hanging_node_modifications false perform required transformations for hanging nodes
keep_matrix true keep matrix between calls to apply() (but reassemble values every time)
verbosity 1 control amount of debug output

Apart from reduction, all parameters have a default value and are optional. The actual reduction for a call to apply() is calculated as r = max(reduction,min_defect/start_defect), where start defect is the norm of the residual of x.

◆ StationaryLinearProblemSolver() [2/2]

template<typename GO , typename LS , typename V >
Dune::PDELab::StationaryLinearProblemSolver< GO, LS, V >::StationaryLinearProblemSolver ( const GO &  go,
LS &  ls,
const ParameterTree params 
)
inline

Construct a StationaryLinearProblemSolver for the given objects and read parameters from a ParameterTree.

This constructor reads the parameter controlling its operation from a passed-in ParameterTree instead of requiring the user to specify all of them as individual constructor parameters. Currently the following parameters are read:

Name Default Value Explanation
reduction Required relative defect reduction
min_defect 1e-99 minimum absolute defect at which to stop
hanging_node_modifications false perform required transformations for hanging nodes
keep_matrix true keep matrix between calls to apply() (but reassemble values every time)
verbosity 1 control amount of debug output

Apart from reduction, all parameters have a default value and are optional. The actual reduction for a call to apply() is calculated as r = max(reduction,min_defect/start_defect), where start defect is the norm of the residual of x.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)