SOPT
Sparse OPTimisation
Classes | Public Types | Public Member Functions | Static Public Member Functions | List of all members
sopt::algorithm::PrimalDual< SCALAR > Class Template Reference

Primal Dual Algorithm. More...

#include <primal_dual.h>

Classes

struct  Diagnostic
 Values indicating how the algorithm ran. More...
 
struct  DiagnosticAndResult
 Holds result vector as well. More...
 

Public Types

using value_type = SCALAR
 Scalar type. More...
 
using Scalar = value_type
 Scalar type. More...
 
using Real = typename real_type< Scalar >::type
 Real type. More...
 
using t_Vector = Vector< Scalar >
 Type of then underlying vectors. More...
 
using t_LinearTransform = LinearTransform< t_Vector >
 Type of the Ψ and Ψ^H operations, as well as Φ and Φ^H. More...
 
using t_IsConverged = std::function< bool(const t_Vector &, const t_Vector &)>
 Type of the convergence function. More...
 
using t_Constraint = std::function< void(t_Vector &, const t_Vector &)>
 Type of the constraint function. More...
 
using t_Random_Updater = std::function< bool()>
 Type of random update function. More...
 
using t_Proximal = ProximalFunction< Scalar >
 Type of the convergence function. More...
 

Public Member Functions

template<typename DERIVED >
 PrimalDual (t_Proximal const &f_proximal, t_Proximal const &g_proximal, Eigen::MatrixBase< DERIVED > const &target)
 
virtual ~PrimalDual ()
 
 SOPT_MACRO (itermax, t_uint)
 Maximum number of iterations. More...
 
 SOPT_MACRO (update_scale, Real)
 Update parameter. More...
 
 SOPT_MACRO (regulariser_strength, Real)
 γ parameter. More...
 
 SOPT_MACRO (sigma, Real)
 sigma parameter More...
 
 SOPT_MACRO (xi, Real)
 xi parameter More...
 
 SOPT_MACRO (rho, Real)
 rho parameter More...
 
 SOPT_MACRO (tau, Real)
 tau parameter More...
 
 SOPT_MACRO (is_converged, t_IsConverged)
 A function verifying convergence. More...
 
 SOPT_MACRO (constraint, t_Constraint)
 A function applying a simple constraint. More...
 
 SOPT_MACRO (Phi, t_LinearTransform)
 Measurement operator. More...
 
 SOPT_MACRO (Psi, t_LinearTransform)
 Wavelet operator. More...
 
 SOPT_MACRO (f_proximal, t_Proximal)
 First proximal. More...
 
 SOPT_MACRO (g_proximal, t_Proximal)
 Second proximal. More...
 
 SOPT_MACRO (random_measurement_updater, t_Random_Updater)
 lambda that determines if to update measurements More...
 
 SOPT_MACRO (random_wavelet_updater, t_Random_Updater)
 lambda that determines if to update wavelets More...
 
void f_proximal (t_Vector &out, Real regulariser_strength, t_Vector const &x) const
 Simplifies calling the proximal of f. More...
 
void g_proximal (t_Vector &out, Real regulariser_strength, t_Vector const &x) const
 Simplifies calling the proximal of f. More...
 
PrimalDual< Scalar > & is_converged (std::function< bool(t_Vector const &x)> const &func)
 Convergence function that takes only the output as argument. More...
 
t_Vector const & target () const
 Vector of target measurements. More...
 
template<typename DERIVED >
PrimalDual< Scalar > & target (Eigen::MatrixBase< DERIVED > const &target)
 Sets the vector of target measurements. More...
 
bool is_converged (t_Vector const &x, t_Vector const &residual) const
 Facilitates call to user-provided convergence function. More...
 
Diagnostic operator() (t_Vector &out) const
 Calls Primal Dual. More...
 
Diagnostic operator() (t_Vector &out, std::tuple< t_Vector, t_Vector > const &guess) const
 Calls Primal Dual. More...
 
Diagnostic operator() (t_Vector &out, std::tuple< t_Vector const &, t_Vector const & > const &guess) const
 Calls Primal Dual. More...
 
DiagnosticAndResult operator() (std::tuple< t_Vector, t_Vector > const &guess) const
 Calls Primal Dual. More...
 
DiagnosticAndResult operator() (std::tuple< t_Vector const &, t_Vector const & > const &guess) const
 Calls Primal Dual. More...
 
DiagnosticAndResult operator() () const
 Calls Primal Dual. More...
 
DiagnosticAndResult operator() (DiagnosticAndResult const &warmstart) const
 Makes it simple to chain different calls to PD. More...
 
PrimalDual &::type Phi (ARGS &&... args)
 
PrimalDual &::type Psi (ARGS &&... args)
 
std::tuple< t_Vector, t_Vectorinitial_guess () const
 Computes initial guess for x and the residual using the targets. More...
 

Static Public Member Functions

static std::tuple< t_Vector, t_Vectorinitial_guess (t_Vector const &target, t_LinearTransform const &phi)
 Computes initial guess for x and the residual using the targets. More...
 

Detailed Description

template<typename SCALAR>
class sopt::algorithm::PrimalDual< SCALAR >

Primal Dual Algorithm.

\(\min_{x, z} f(x) + h(z)\) subject to \(Φx + z = y\). \(y\) is a target vector.

Definition at line 24 of file primal_dual.h.

Member Typedef Documentation

◆ Real

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::Real = typename real_type<Scalar>::type

Real type.

Definition at line 31 of file primal_dual.h.

◆ Scalar

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::Scalar = value_type

Scalar type.

Definition at line 29 of file primal_dual.h.

◆ t_Constraint

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::t_Constraint = std::function<void (t_Vector &, const t_Vector &)>

Type of the constraint function.

Definition at line 39 of file primal_dual.h.

◆ t_IsConverged

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::t_IsConverged = std::function<bool (const t_Vector &, const t_Vector &)>

Type of the convergence function.

Definition at line 37 of file primal_dual.h.

◆ t_LinearTransform

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::t_LinearTransform = LinearTransform<t_Vector>

Type of the Ψ and Ψ^H operations, as well as Φ and Φ^H.

Definition at line 35 of file primal_dual.h.

◆ t_Proximal

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::t_Proximal = ProximalFunction<Scalar>

Type of the convergence function.

Definition at line 43 of file primal_dual.h.

◆ t_Random_Updater

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::t_Random_Updater = std::function<bool ()>

Type of random update function.

Definition at line 41 of file primal_dual.h.

◆ t_Vector

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::t_Vector = Vector<Scalar>

Type of then underlying vectors.

Definition at line 33 of file primal_dual.h.

◆ value_type

template<typename SCALAR >
using sopt::algorithm::PrimalDual< SCALAR >::value_type = SCALAR

Scalar type.

Definition at line 27 of file primal_dual.h.

Constructor & Destructor Documentation

◆ PrimalDual()

template<typename SCALAR >
template<typename DERIVED >
sopt::algorithm::PrimalDual< SCALAR >::PrimalDual ( t_Proximal const &  f_proximal,
t_Proximal const &  g_proximal,
Eigen::MatrixBase< DERIVED > const &  target 
)
inline

Setups PrimalDual

Parameters
[in]f_proximalproximal operator of the \(f\) function.
[in]g_proximalproximal operator of the \(g\) function

Definition at line 69 of file primal_dual.h.

71  : itermax_(std::numeric_limits<t_uint>::max()),
72  sigma_(1),
73  tau_(0.5),
74  regulariser_strength_(0.5),
75  update_scale_(1),
76  xi_(1),
77  rho_(1),
78  is_converged_(),
79  constraint_([](t_Vector &out, t_Vector const &x) { out = x; }),
80  Phi_(linear_transform_identity<Scalar>()),
81  Psi_(linear_transform_identity<Scalar>()),
82  f_proximal_(f_proximal),
83  g_proximal_(g_proximal),
84  random_measurement_updater_([]() { return true; }),
85  random_wavelet_updater_([]() { return true; }),
86 #ifdef SOPT_MPI
87  v_all_sum_all_comm_(mpi::Communicator()),
88  u_all_sum_all_comm_(mpi::Communicator()),
89 #endif
90  target_(target) {
91  }
sopt::Vector< Scalar > t_Vector
t_Vector const & target() const
Vector of target measurements.
Definition: primal_dual.h:161
void f_proximal(t_Vector &out, Real regulariser_strength, t_Vector const &x) const
Simplifies calling the proximal of f.
Definition: primal_dual.h:147
void g_proximal(t_Vector &out, Real regulariser_strength, t_Vector const &x) const
Simplifies calling the proximal of f.
Definition: primal_dual.h:151
#define SOPT_MPI
Whether or not to include mpi.
Definition: config.in.h:17

◆ ~PrimalDual()

template<typename SCALAR >
virtual sopt::algorithm::PrimalDual< SCALAR >::~PrimalDual ( )
inlinevirtual

Definition at line 92 of file primal_dual.h.

92 {}

Member Function Documentation

◆ f_proximal()

template<typename SCALAR >
void sopt::algorithm::PrimalDual< SCALAR >::f_proximal ( t_Vector out,
Real  regulariser_strength,
t_Vector const &  x 
) const
inline

Simplifies calling the proximal of f.

Definition at line 147 of file primal_dual.h.

147  {
148  f_proximal()(out, regulariser_strength, x);
149  }

◆ g_proximal()

template<typename SCALAR >
void sopt::algorithm::PrimalDual< SCALAR >::g_proximal ( t_Vector out,
Real  regulariser_strength,
t_Vector const &  x 
) const
inline

Simplifies calling the proximal of f.

Definition at line 151 of file primal_dual.h.

151  {
152  g_proximal()(out, regulariser_strength, x);
153  }

◆ initial_guess() [1/2]

template<typename SCALAR >
std::tuple<t_Vector, t_Vector> sopt::algorithm::PrimalDual< SCALAR >::initial_guess ( ) const
inline

Computes initial guess for x and the residual using the targets.

with y the vector of measurements

  • x = Φ^T y / nu = Φ^T y / (Φ_norm^2)
  • residuals = Φ x - y

Definition at line 233 of file primal_dual.h.

233  {
235  }
PrimalDual &::type Phi(ARGS &&... args)
Definition: primal_dual.h:218
std::tuple< t_Vector, t_Vector > initial_guess() const
Computes initial guess for x and the residual using the targets.
Definition: primal_dual.h:233

References sopt::algorithm::PrimalDual< SCALAR >::Phi(), and sopt::algorithm::PrimalDual< SCALAR >::target().

Referenced by sopt::algorithm::ImagingPrimalDual< SCALAR >::operator()(), sopt::algorithm::PrimalDual< SCALAR >::operator()(), and sopt::algorithm::TVPrimalDual< SCALAR >::operator()().

◆ initial_guess() [2/2]

template<typename SCALAR >
static std::tuple<t_Vector, t_Vector> sopt::algorithm::PrimalDual< SCALAR >::initial_guess ( t_Vector const &  target,
t_LinearTransform const &  phi 
)
inlinestatic

Computes initial guess for x and the residual using the targets.

with y the vector of measurements

  • x = Φ^T y / nu = Φ^T y / (Φ_norm^2)
  • residuals = Φ x - y

This function simplifies creating overloads for operator() in PD wrappers.

Definition at line 243 of file primal_dual.h.

244  {
245  std::tuple<t_Vector, t_Vector> guess;
246  std::get<0>(guess) = static_cast<t_Vector>(phi.adjoint() * target) / phi.sq_norm();
247  std::get<1>(guess) = target;
248  return guess;
249  }

References sopt::LinearTransform< VECTOR >::adjoint(), sopt::LinearTransform< VECTOR >::sq_norm(), and sopt::algorithm::PrimalDual< SCALAR >::target().

◆ is_converged() [1/2]

template<typename SCALAR >
PrimalDual<Scalar>& sopt::algorithm::PrimalDual< SCALAR >::is_converged ( std::function< bool(t_Vector const &x)> const &  func)
inline

Convergence function that takes only the output as argument.

Definition at line 156 of file primal_dual.h.

156  {
157  return is_converged([func](t_Vector const &x, t_Vector const &) { return func(x); });
158  }
PrimalDual< Scalar > & is_converged(std::function< bool(t_Vector const &x)> const &func)
Convergence function that takes only the output as argument.
Definition: primal_dual.h:156

Referenced by sopt::algorithm::PrimalDual< SCALAR >::is_converged(), and TEST_CASE().

◆ is_converged() [2/2]

template<typename SCALAR >
bool sopt::algorithm::PrimalDual< SCALAR >::is_converged ( t_Vector const &  x,
t_Vector const &  residual 
) const
inline

Facilitates call to user-provided convergence function.

Definition at line 170 of file primal_dual.h.

170  {
171  return static_cast<bool>(is_converged()) and is_converged()(x, residual);
172  }

References sopt::algorithm::PrimalDual< SCALAR >::is_converged().

◆ operator()() [1/7]

template<typename SCALAR >
DiagnosticAndResult sopt::algorithm::PrimalDual< SCALAR >::operator() ( ) const
inline

Calls Primal Dual.

Parameters
[in]guessinitial guess

Definition at line 205 of file primal_dual.h.

205  {
206  DiagnosticAndResult result;
207  static_cast<Diagnostic &>(result) = operator()(result.x, initial_guess());
208  return result;
209  }

References sopt::algorithm::PrimalDual< SCALAR >::initial_guess(), and sopt::algorithm::PrimalDual< SCALAR >::DiagnosticAndResult::x.

Referenced by sopt::algorithm::PrimalDual< SCALAR >::operator()().

◆ operator()() [2/7]

template<typename SCALAR >
DiagnosticAndResult sopt::algorithm::PrimalDual< SCALAR >::operator() ( DiagnosticAndResult const &  warmstart) const
inline

Makes it simple to chain different calls to PD.

Definition at line 211 of file primal_dual.h.

211  {
212  DiagnosticAndResult result = warmstart;
213  static_cast<Diagnostic &>(result) = operator()(result.x, warmstart.x, warmstart.residual);
214  return result;
215  }

References sopt::algorithm::PrimalDual< SCALAR >::Diagnostic::residual, and sopt::algorithm::PrimalDual< SCALAR >::DiagnosticAndResult::x.

◆ operator()() [3/7]

template<typename SCALAR >
DiagnosticAndResult sopt::algorithm::PrimalDual< SCALAR >::operator() ( std::tuple< t_Vector const &, t_Vector const & > const &  guess) const
inline

Calls Primal Dual.

Parameters
[in]guessinitial guess

Definition at line 197 of file primal_dual.h.

198  {
199  DiagnosticAndResult result;
200  static_cast<Diagnostic &>(result) = operator()(result.x, guess);
201  return result;
202  }

References sopt::algorithm::PrimalDual< SCALAR >::DiagnosticAndResult::x.

◆ operator()() [4/7]

template<typename SCALAR >
DiagnosticAndResult sopt::algorithm::PrimalDual< SCALAR >::operator() ( std::tuple< t_Vector, t_Vector > const &  guess) const
inline

Calls Primal Dual.

Parameters
[in]guessinitial guess

Definition at line 192 of file primal_dual.h.

192  {
193  return operator()(std::tie(std::get<0>(guess), std::get<1>(guess)));
194  }
DiagnosticAndResult operator()() const
Calls Primal Dual.
Definition: primal_dual.h:205

References sopt::algorithm::PrimalDual< SCALAR >::operator()().

◆ operator()() [5/7]

template<typename SCALAR >
Diagnostic sopt::algorithm::PrimalDual< SCALAR >::operator() ( t_Vector out) const
inline

Calls Primal Dual.

Parameters
[out]outOutput vector x

Definition at line 176 of file primal_dual.h.

176 { return operator()(out, initial_guess()); }

References sopt::algorithm::PrimalDual< SCALAR >::initial_guess(), and sopt::algorithm::PrimalDual< SCALAR >::operator()().

Referenced by sopt::algorithm::PrimalDual< SCALAR >::operator()().

◆ operator()() [6/7]

template<typename SCALAR >
Diagnostic sopt::algorithm::PrimalDual< SCALAR >::operator() ( t_Vector out,
std::tuple< t_Vector const &, t_Vector const & > const &  guess 
) const
inline

Calls Primal Dual.

Parameters
[out]outOutput vector x
[in]guessinitial guess

Definition at line 186 of file primal_dual.h.

187  {
188  return operator()(out, std::get<0>(guess), std::get<1>(guess));
189  }

References sopt::algorithm::PrimalDual< SCALAR >::operator()().

◆ operator()() [7/7]

template<typename SCALAR >
Diagnostic sopt::algorithm::PrimalDual< SCALAR >::operator() ( t_Vector out,
std::tuple< t_Vector, t_Vector > const &  guess 
) const
inline

Calls Primal Dual.

Parameters
[out]outOutput vector x
[in]guessinitial guess

Definition at line 180 of file primal_dual.h.

180  {
181  return operator()(out, std::get<0>(guess), std::get<1>(guess));
182  }

References sopt::algorithm::PrimalDual< SCALAR >::operator()().

◆ Phi()

template<typename SCALAR >
PrimalDual& ::type sopt::algorithm::PrimalDual< SCALAR >::Phi ( ARGS &&...  args)
inline

Definition at line 218 of file primal_dual.h.

218  {
219  Phi_ = linear_transform(std::forward<ARGS>(args)...);
220  return *this;
221  }
LinearTransform< VECTOR > linear_transform(OperatorFunction< VECTOR > const &direct, OperatorFunction< VECTOR > const &indirect, std::array< t_int, 3 > const &sizes={{1, 1, 0}})

References sopt::linear_transform().

Referenced by sopt::algorithm::PrimalDual< SCALAR >::initial_guess().

◆ Psi()

template<typename SCALAR >
PrimalDual& ::type sopt::algorithm::PrimalDual< SCALAR >::Psi ( ARGS &&...  args)
inline

Definition at line 224 of file primal_dual.h.

224  {
225  Psi_ = linear_transform(std::forward<ARGS>(args)...);
226  return *this;
227  }

References sopt::linear_transform().

◆ SOPT_MACRO() [1/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( constraint  ,
t_Constraint   
)

A function applying a simple constraint.

◆ SOPT_MACRO() [2/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( f_proximal  ,
t_Proximal   
)

First proximal.

◆ SOPT_MACRO() [3/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( g_proximal  ,
t_Proximal   
)

Second proximal.

◆ SOPT_MACRO() [4/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( is_converged  ,
t_IsConverged   
)

A function verifying convergence.

It takes as input two arguments: the current solution x and the current residual.

◆ SOPT_MACRO() [5/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( itermax  ,
t_uint   
)

Maximum number of iterations.

◆ SOPT_MACRO() [6/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( Phi  ,
t_LinearTransform   
)

Measurement operator.

◆ SOPT_MACRO() [7/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( Psi  ,
t_LinearTransform   
)

Wavelet operator.

◆ SOPT_MACRO() [8/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( random_measurement_updater  ,
t_Random_Updater   
)

lambda that determines if to update measurements

◆ SOPT_MACRO() [9/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( random_wavelet_updater  ,
t_Random_Updater   
)

lambda that determines if to update wavelets

◆ SOPT_MACRO() [10/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( regulariser_strength  ,
Real   
)

γ parameter.

◆ SOPT_MACRO() [11/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( rho  ,
Real   
)

rho parameter

◆ SOPT_MACRO() [12/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( sigma  ,
Real   
)

sigma parameter

◆ SOPT_MACRO() [13/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( tau  ,
Real   
)

tau parameter

◆ SOPT_MACRO() [14/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( update_scale  ,
Real   
)

Update parameter.

◆ SOPT_MACRO() [15/15]

template<typename SCALAR >
sopt::algorithm::PrimalDual< SCALAR >::SOPT_MACRO ( xi  ,
Real   
)

xi parameter

◆ target() [1/2]

template<typename SCALAR >
t_Vector const& sopt::algorithm::PrimalDual< SCALAR >::target ( ) const
inline

Vector of target measurements.

Definition at line 161 of file primal_dual.h.

161 { return target_; }

Referenced by sopt::algorithm::PrimalDual< SCALAR >::initial_guess(), and sopt::algorithm::PrimalDual< SCALAR >::target().

◆ target() [2/2]

template<typename SCALAR >
template<typename DERIVED >
PrimalDual<Scalar>& sopt::algorithm::PrimalDual< SCALAR >::target ( Eigen::MatrixBase< DERIVED > const &  target)
inline

Sets the vector of target measurements.

Definition at line 164 of file primal_dual.h.

164  {
165  target_ = target;
166  return *this;
167  }

References sopt::algorithm::PrimalDual< SCALAR >::target().


The documentation for this class was generated from the following file: