So-Bogus
A c++ sparse block matrix library aimed at Second Order cone problems
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Pages
ProjectedGradient.hpp
1 /*
2  * This file is part of bogus, a C++ sparse block matrix library.
3  *
4  * Copyright 2013 Gilles Daviet <gdaviet@gmail.com>
5  *
6  * This Source Code Form is subject to the terms of the Mozilla Public
7  * License, v. 2.0. If a copy of the MPL was not distributed with this
8  * file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 */
10 
11 
12 #ifndef BOGUS_BLOCK_PROJECTED_GRADIENT_HPP
13 #define BOGUS_BLOCK_PROJECTED_GRADIENT_HPP
14 
15 #include "ConstrainedSolverBase.hpp"
16 
17 #include <vector>
18 
19 namespace bogus
20 {
21 
23 namespace projected_gradient {
25  enum Variant {
36  } ;
37 }
38 
40 template < typename BlockMatrixType >
41 class ProjectedGradient : public ConstrainedSolverBase< ProjectedGradient<BlockMatrixType >, BlockMatrixType >
42 {
43 public:
45  typedef typename Base::Scalar Scalar ;
46 
48  ProjectedGradient( ) : Base() { init() ; }
51  { init() ; Base::setMatrix( matrix ) ; }
52 
54  template < typename NSLaw, typename RhsT, typename ResT >
55  Scalar solve( const NSLaw &law, const RhsT &b, ResT &x ) const ;
56 
58 
67  template < projected_gradient::Variant variant, typename NSLaw, typename RhsT, typename ResT >
68  Scalar solve( const NSLaw &law, const RhsT &b, ResT &x ) const ;
69 
71 
73  {
74  m_matrix = &matrix ;
75  Base::updateScalings() ;
76  return *this ;
77  }
78 
80  void setLineSearchIterations( const unsigned lsIterations )
81  { m_lsIters = lsIterations ; }
82 
85  void setLineSearchOptimisticFactor( const Scalar lsOptimisticFactor )
86  { m_lsOptimisticFactor = lsOptimisticFactor ; }
87 
90  void setLineSearchPessimisticFactor( const Scalar lsPessimisticFactor )
91  { m_lsPessimisticFactor = lsPessimisticFactor ; }
92 
95  void setLineSearchArmijoCoefficient( const Scalar lsArmijoCoefficient )
96  { m_lsArmijoCoefficient = lsArmijoCoefficient ; }
97 
100  { m_defaultVariant = variant ; }
101 
102  unsigned lineSearchIterations() const { return m_lsIters ; }
103  Scalar lineSearchOptimisticFactor() const { return m_lsOptimisticFactor ; }
104  Scalar lineSearchPessimisticFactor() const { return m_lsPessimisticFactor ; }
105  Scalar lineSearchArmijoCoefficient() const { return m_lsArmijoCoefficient ; }
106 
107 protected:
108 
109  typedef typename Base::Index Index ;
110 
112  void init()
113  {
114  m_tol = 1.e-6 ;
115  m_maxIters = 300 ;
116  m_lsIters = 8 ;
117  m_lsOptimisticFactor = 1.25 ;
118  m_lsPessimisticFactor = .5 ;
119  m_lsArmijoCoefficient = .5 ;
120  m_defaultVariant = projected_gradient::APGD ;
121  }
122 
123  using Base::m_matrix ;
124  using Base::m_maxIters ;
125  using Base::m_tol ;
126 
127  unsigned m_lsIters ;
128  Scalar m_lsOptimisticFactor ;
129  Scalar m_lsPessimisticFactor ;
130  Scalar m_lsArmijoCoefficient ;
131 
132  projected_gradient::Variant m_defaultVariant ;
133 
134 } ;
135 
136 } //namespace bogus
137 
138 
139 #endif
140 
void setLineSearchArmijoCoefficient(const Scalar lsArmijoCoefficient)
Definition: ProjectedGradient.hpp:95
void setLineSearchPessimisticFactor(const Scalar lsPessimisticFactor)
Definition: ProjectedGradient.hpp:90
void setDefaultVariant(projected_gradient::Variant variant)
Sets the variant that will be used when calling solve() without template arguments.
Definition: ProjectedGradient.hpp:99
Scalar m_tol
See setTol()
Definition: BlockSolverBase.hpp:70
Accelerated Projected Gradient Descent based on and developed in .
Definition: ProjectedGradient.hpp:33
Scalar solve(const NSLaw &law, const RhsT &b, ResT &x) const
Finds an approximate minimum for a constrained quadratic problem.
Projected Gradient iterative solver.
Definition: ProjectedGradient.hpp:41
void setLineSearchOptimisticFactor(const Scalar lsOptimisticFactor)
Definition: ProjectedGradient.hpp:85
Spectral Projected Gradient, loosely adapted from .
Definition: ProjectedGradient.hpp:35
Projected gradient descent.
Definition: ProjectedGradient.hpp:29
Derived & setMatrix(const BlockObjectBase< BlockMatrixType > &matrix)
Sets the system matrix and initializes internal structures.
Definition: ConstrainedSolverBase.hpp:62
void setLineSearchIterations(const unsigned lsIterations)
Sets the maximum number of line-search iterations.
Definition: ProjectedGradient.hpp:80
const BlockObjectBase< BlockMatrixType > * m_matrix
Pointer to the matrix of the system.
Definition: BlockSolverBase.hpp:65
Projected gradient with conjugation of search direction.
Definition: ProjectedGradient.hpp:31
void init()
Sets up the default values for all parameters.
Definition: ProjectedGradient.hpp:112
unsigned m_maxIters
See setMaxIters()
Definition: BlockSolverBase.hpp:68
ProjectedGradient()
Default constructor – you will have to call setMatrix() before using the solve() function.
Definition: ProjectedGradient.hpp:48
Variant
Variants of Projected Gradient algorithm.
Definition: ProjectedGradient.hpp:25
Definition: ConstrainedSolverBase.hpp:20
Standard projected gradient.
Definition: ProjectedGradient.hpp:27
ProjectedGradient(const BlockObjectBase< BlockMatrixType > &matrix)
Constructor with the system matrix.
Definition: ProjectedGradient.hpp:50
ProjectedGradient & setMatrix(const BlockObjectBase< BlockMatrixType > &matrix)
Sets the matrix M defining the quadratic objective function.
Definition: ProjectedGradient.hpp:72