KASKADE 7 development version
|
A base class for multigrid solver for QPs. More...
#include <qpmg.hh>
A base class for multigrid solver for QPs.
This base class implements a multigrid solver for problems of the form
\[ \min_x \frac{1}{2} x^T Ax + c^T x \quad \text{s.t. } Bx \le b, \]
where \( A \) is positive semidefinite defined on a hierarchy of grids, i.e. for a given prolongation stack defining the (geometric) relation between a sequence of grids.
To the QP on fine grids we apply a block Gauss-Seidel method, which treats DOFs coupled by the constraints
\[ Bx \le b \]
jointly. These local QPs are then solved by the smoother on the corresponsing level. On the corse grid the system is solver either by a direct solver (dense or sparse) or by a few smoothing steps.
d | the dimension of (vector valued) optimization variables |
Prolongation | defines the data type of the prolongation matrices defining the hierarchy |
Smoother | the smoother for the local QPs on all grids but the coarse one |
CoarseSmoother | the smoother for the coarse grid |
Public Types | |
using | VectorX = Dune::BlockVector< Dune::FieldVector< Real, d > > |
using | VectorB = Dune::BlockVector< Dune::FieldVector< Real, 1 > > |
using | MatrixA = NumaBCRSMatrix< EntryA > |
using | MatrixB = NumaBCRSMatrix< Dune::FieldMatrix< Real, 1, d > > |
using | Self = QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real > |
using | SmootherType = std::conditional_t< smoothersEqual, Smoother, std::variant< Smoother, CoarseSmoother > > |
Public Member Functions | |
QPMultiGridBase (MatrixA A, MatrixB const &B, std::vector< Prolongation > &&prolongations, Real smootherRegularization=0, bool blocks=true, bool directOnCoarse=true) | |
Constructs the solver, providing the matrices \( A \) and \( B \). More... | |
virtual | ~QPMultiGridBase () |
Destructor. More... | |
VectorX | step (VectorX const &x, VectorX c, VectorB b) const |
Computes a single V-cycle correction for approximately solving the QP for the given right hand side vectors. More... | |
std::tuple< VectorX, int > | solve (VectorX x, VectorX const &c, VectorB const &b, double tol, int vcycles) const |
Approximately solves the QP for the given right hand side vectors up to a given tolerance. More... | |
Self & | setCoarseLevel (int coarselevel) |
Defines the grid level to be used as coarse grid. More... | |
Self & | setSmoothings (int pre, int post) |
Self & | setCoarseCorrections (int n) |
Sets the number of coarse corrections to compute and apply. More... | |
Self & | setBulkMode (ParallelMode m) |
Defines the smoothing strategy to use for unconstrained degrees of freedom. More... | |
Self & | setLogger (MGSolverStatistics< d, Real > *logger) |
Sets the logging facility for reporting solver statistics. More... | |
Static Public Attributes | |
static constexpr bool | smoothersEqual = std::is_same<Smoother,CoarseSmoother>::value |
Protected Member Functions | |
virtual double | computeEnergy (MatrixA const &A, MatrixB const &B, VectorX const &c, VectorB const &b, VectorX const &x) const =0 |
Compute the problem energy of current iterate for logging purposes. More... | |
virtual std::tuple< VectorX, VectorB, int > | gradient (MatrixA const &A, MatrixB const &B, VectorX const &c, VectorB const &b, VectorX const &x) const =0 |
Compute the problem gradient of current iterate. More... | |
virtual std::tuple< VectorX, VectorX, VectorX, VectorB, int > | gradient_extended (MatrixA const &A, MatrixB const &B, VectorX const &c, VectorB const &b, VectorX const &x) const =0 |
Computes the gradient of current iterate, giving more information on the components of the gradient. (It might be split in several terms grad0 and grad1). More... | |
virtual double | qpLinesearch (MatrixA const &A, MatrixB const &B, VectorX const &c, VectorB const &b, VectorX const &dx) const =0 |
Compute the optimal step size in direction of dx. More... | |
MatrixB const & | B () const |
Provides a reference to the fine grid constraints. More... | |
Protected Attributes | |
std::vector< SmootherType > | smoothers |
using Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::MatrixA = NumaBCRSMatrix<EntryA> |
using Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::MatrixB = NumaBCRSMatrix<Dune::FieldMatrix<Real,1,d> > |
using Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::Self = QPMultiGridBase<d,Prolongation,Smoother,CoarseSmoother,Real> |
using Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::SmootherType = std::conditional_t<smoothersEqual, Smoother, std::variant<Smoother,CoarseSmoother> > |
using Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::VectorB = Dune::BlockVector<Dune::FieldVector<Real,1> > |
using Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::VectorX = Dune::BlockVector<Dune::FieldVector<Real,d> > |
Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::QPMultiGridBase | ( | MatrixA | A, |
MatrixB const & | B, | ||
std::vector< Prolongation > && | prolongations, | ||
Real | smootherRegularization = 0 , |
||
bool | blocks = true , |
||
bool | directOnCoarse = true |
||
) |
Constructs the solver, providing the matrices \( A \) and \( B \).
A | the Hessian. The matrix is referenced internally and has to exist during the lifetime of the smoother. |
B | the constraints. The matrix is referenced internally and has to exist during the lifetime of the smoother. |
blocks |
|
directOnCoarse |
|
|
virtual |
Destructor.
|
protected |
Provides a reference to the fine grid constraints.
|
protectedpure virtual |
Compute the problem energy of current iterate for logging purposes.
|
protectedpure virtual |
Compute the problem gradient of current iterate.
|
protectedpure virtual |
Computes the gradient of current iterate, giving more information on the components of the gradient. (It might be split in several terms grad0 and grad1).
|
protectedpure virtual |
Compute the optimal step size in direction of dx.
Self & Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::setBulkMode | ( | ParallelMode | m | ) |
Defines the smoothing strategy to use for unconstrained degrees of freedom.
In contact problems, those are usually the interior nodes, and the ones without potential contact partner. The smoother can operate in parallel (Jacobi) or sequentially (Gauss-Seidel). In the parallel way, an exact linesearch is performed after the smoothing in order to guarantee descent.
Self & Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::setCoarseCorrections | ( | int | n | ) |
Sets the number of coarse corrections to compute and apply.
n | number of coarse corrections |
Self & Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::setCoarseLevel | ( | int | coarselevel | ) |
Defines the grid level to be used as coarse grid.
coarselevel | the coarse grid level, between 0 and the max grid level |
If coarselevel is outside its allowed range, it is projected onto that range.
Self & Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::setLogger | ( | MGSolverStatistics< d, Real > * | logger | ) |
Sets the logging facility for reporting solver statistics.
logger | a pointer to the logger object. This can be a nullpointer, in which case logging is stopped. |
Ownership of the logger object remains with the caller. The logger object needs to exist for the lifetime of the QPMultiGrid class or until it is replaced by a different logger (or none).
Self & Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::setSmoothings | ( | int | pre, |
int | post | ||
) |
std::tuple< VectorX, int > Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::solve | ( | VectorX | x, |
VectorX const & | c, | ||
VectorB const & | b, | ||
double | tol, | ||
int | vcycles | ||
) | const |
Approximately solves the QP for the given right hand side vectors up to a given tolerance.
x | initial guess for the solution |
c | linear objective term |
b | constant constraint term |
tol | tolerance |
VectorX Kaskade::QPMultiGridBase< d, Prolongation, Smoother, CoarseSmoother, Real >::step | ( | VectorX const & | x, |
VectorX | c, | ||
VectorB | b | ||
) | const |
Computes a single V-cycle correction for approximately solving the QP for the given right hand side vectors.
x | initial guess for the solution |
c | linear objective term |
b | constant constraint term |
|
protected |
|
staticconstexpr |