KASKADE 7 development version
|
A multigrid solver for convex QPs. More...
#include <qpmg.hh>
A multigrid solver for convex QPs.
This solves QPs of the form
\[ \min \frac{1}{2} x^T A x + c^T x \quad \text{s.t.} \quad Bx\le b. \]
On the finest level, this works with the augmented Lagrangian formulation
\[ \min \frac{1}{2} x^T A x + c^T x + \frac{\gamma}{2} \|(Bx-(b-\lambda/\gamma))_+\|^2 \]
for some multiplier \( \lambda \ge 0 \), and an appropriate multiplier update.
The multilevel minimization works with level-dependent penalty factors. On each level, the following computations ( \(\text{MGM}(A,B,x,c,b,\gamma)\) are performed:
\[ (\delta x, \delta \lambda) = \text{MGM}(P^TAP,BP,0,P^Tc_r,b_r,\gamma/4) \]
Rationale: The level-dependence of the penalty factor is based on the following idea. Enforcing the constraint by the penalty should allow sufficient slack such as not to impede convergence of Gauß-Seidel. On coarser levels, GS is supposed to take larger steps (the whole reason for multigrid) approximately by a factor of two for each grid level, and so the penalty should increase the permitted steps by a factor of two as well. This means decreasing the penalty factor by four.
The restriction of multipliers is not necessary for asymptotically optimal complexity, as the number of constraints is \( \mathcal{O}(N^{1-1/d}) \) and therefore the computational work of a single V-cycle is \( \mathcal{O}(N + N^{1-1/d}\log N) = \mathcal{O}(N) \). Nevertheless, restricting the constraints in an appropriate way (open how to do this) might be beneficial for performance.
d | the dimension of (vector valued) optimization variables |
Prolongation | defines the data type of the prolongation matrices defining the hierarchy |
Real | floating point format |
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 = QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother > |
using | SmootherType = std::conditional_t< smoothersEqual, QPPenaltySmoother< d, double >, std::variant< QPPenaltySmoother< d, double >, QPPenaltySmoother< d, double > > > |
Public Member Functions | |
QPMultiGrid (MatrixA const &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... | |
~QPMultiGrid () | |
Destructor. More... | |
std::tuple< VectorX, VectorB, int > | solve (VectorX x, VectorX const &c, VectorB const &b, VectorB const &lambda, double tol, int vcycles) const |
Approximately solves the QP for the given right hand side vectors up to a given tolerance. More... | |
VectorB | firstOrderUpdate (VectorX const &x, VectorB const &b, VectorB const &lambda) const |
A first order update for the Lagrange multiplier. More... | |
Self & | setPenalty (Real gamma, bool levelDependent=false) |
Sets the penalty factor \( \gamma \). More... | |
Self & | setConstraintsMode (ParallelMode m) |
Self & | setGlobalMode (ParallelMode m) |
double | getPenalty () |
Returns the current penalty factor \( \gamma \). 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, double > *logger) |
Sets the logging facility for reporting solver statistics. More... | |
Static Public Attributes | |
static constexpr bool | smoothersEqual |
Protected Member Functions | |
virtual double | computeEnergy (MatrixA const &A, MatrixB const &B, VectorX const &c, VectorB const &b, VectorX const &x) const |
Computes the energy, given by the penalized energy function. More... | |
virtual std::tuple< VectorX, VectorB, int > | gradient (MatrixA const &A, MatrixB const &B, VectorX const &c, VectorB const &b, VectorX const &x) const |
Computes the gradient of the penalized energy function. 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 |
Computes the gradient of the penalized energy function. More... | |
virtual double | qpLinesearch (MatrixA const &A, MatrixB const &B, VectorX const &c, VectorB const &b, VectorX const &dx) const |
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::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::MatrixA = NumaBCRSMatrix<EntryA> |
using Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::MatrixB = NumaBCRSMatrix<Dune::FieldMatrix<Real,1,d> > |
using Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::Self = QPMultiGrid<d,Prolongation,Real,Smoother,CoarseSmoother> |
|
inherited |
using Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::VectorB = Dune::BlockVector<Dune::FieldVector<Real,1> > |
using Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::VectorX = Dune::BlockVector<Dune::FieldVector<Real,d> > |
Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::QPMultiGrid | ( | MatrixA const & | 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 copied internally and need not exist after solver construction. |
B | the constraints. The matrix is referenced internally and has to exist during the lifetime of the solver. |
blocks |
|
directOnCoarse |
|
Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::~QPMultiGrid | ( | ) |
Destructor.
|
protectedinherited |
Provides a reference to the fine grid constraints.
|
protectedvirtual |
Computes the energy, given by the penalized energy function.
\[ \frac{1}{2} x^T A x + c^T x + \frac{\gamma}{2} \|(Bx-b)_+\|^2 \]
|
protectedpure virtualinherited |
Compute the problem energy of current iterate for logging purposes.
VectorB Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::firstOrderUpdate | ( | VectorX const & | x, |
VectorB const & | b, | ||
VectorB const & | lambda | ||
) | const |
A first order update for the Lagrange multiplier.
Compute first order multiplier update \lambda <- (\lambda + \gamma*(B*x-b))_+, or equivalently \lambda <- \lambda + max(-\lambda, \gamma*(B*x-b)). This makes only sense, however, if we have actually minimized - otherwise there is no new information in the constraint violation. Thus the accuracy of minimization should be sufficient.
x | primal variable |
b | constant constraint term |
lambda | approximate multiplier |
double Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::getPenalty | ( | ) |
Returns the current penalty factor \( \gamma \).
|
protectedvirtual |
Computes the gradient of the penalized energy function.
\[ \frac{1}{2} x^T A x + c^T x + \frac{\gamma}{2} \|(Bx-b)_+\|^2 \]
\[ Ax + c + \gamma*B^T(Bx-b)_+ \]
\[ (Bx-b)_+ \]
\[ | \{ i \colon (Bx-b)_i >= 0 \} | \]
|
protectedpure virtualinherited |
Compute the problem gradient of current iterate.
|
protectedvirtual |
Computes the gradient of the penalized energy function.
\[ \frac{1}{2} x^T A x + c^T x + \frac{\gamma}{2} \|(Bx-b)_+\|^2 \]
\[ Ax + c + \gamma*B^T(Bx-b)_+ \]
\[ Ax + c \]
\[ \gamma*B^T(Bx-b)_+ \]
\[ (Bx-b)_+ \]
\[ | \{ i \colon (Bx-b)_i >= 0 \} | \]
|
protectedpure virtualinherited |
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).
|
protectedvirtual |
|
protectedpure virtualinherited |
Compute the optimal step size in direction of dx.
|
inherited |
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.
|
inherited |
Sets the number of coarse corrections to compute and apply.
n | number of coarse corrections |
|
inherited |
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::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::setConstraintsMode | ( | ParallelMode | m | ) |
Self & Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::setGlobalMode | ( | ParallelMode | m | ) |
|
inherited |
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::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::setPenalty | ( | Real | gamma, |
bool | levelDependent = false |
||
) |
Sets the penalty factor \( \gamma \).
The default on construction is \( \gamma = 10^3 \).
gamma | the penalty parameter |
levelDependent | if set true, we use level dependent penaltys, i.e., gamma on finest level and gamma/4 recursively on coarser levels |
|
inherited |
|
inherited |
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 |
std::tuple< VectorX, VectorB, int > Kaskade::QPMultiGrid< d, Prolongation, Real, Smoother, CoarseSmoother >::solve | ( | VectorX | x, |
VectorX const & | c, | ||
VectorB const & | b, | ||
VectorB const & | lambda, | ||
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 |
lambda | approximate multiplier |
tol | tolerance |
|
inherited |
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 |
|
protectedinherited |
|
staticconstexprinherited |