numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
|
Class to solve linear equations with 2nd order total generalized variation (TGV) regularization [3] using the alternating direction method of multipliers (ADMM) [2]. More...
#include <num_collect/regularization/tgv2_admm.h>
Public Types | |
using | base_type = iterative_regularized_solver_base<this_type, Data> |
Type of the base class. | |
using | coeff_type = Coeff |
Type of coefficient matrices. | |
using | derivative_matrix_type = DerivativeMatrix |
Type of matrices to compute derivatives. | |
using | this_type = tgv2_admm<Coeff, DerivativeMatrix, Data> |
This type. | |
Public Types inherited from num_collect::regularization::regularized_solver_base< tgv2_admm< Coeff, DerivativeMatrix, Data >, Data > | |
using | data_type |
Type of data. | |
using | scalar_type |
Type of scalars. |
Public Member Functions | |
tgv2_admm () | |
Constructor. | |
void | calculate_data_for (const data_type &solution, data_type &data) const |
Calculate data for a solution. | |
void | change_data (const data_type &data) |
Change data. | |
void | compute (const Coeff &coeff, const DerivativeMatrix &first_derivative_matrix, const DerivativeMatrix &second_derivative_matrix, const Data &data) |
Compute internal parameters. | |
void | configure_iteration_logger (logging::iterations::iteration_logger< this_type > &iteration_logger) const |
Configure an iteration logger. | |
auto | data_size () const -> index_type |
Get the size of data. | |
void | init (const scalar_type ¶m, data_type &solution) |
Initialize. | |
auto | is_stop_criteria_satisfied (const data_type &solution) const -> bool |
Determine if stopping criteria of the algorithm are satisfied. | |
void | iterate (const scalar_type ¶m, data_type &solution) |
Iterate the algorithm once. | |
auto | iterations () const noexcept -> index_type |
Get the number of iterations. | |
auto | max_iterations () const -> index_type |
Get the maximum number of iterations. | |
auto | max_iterations (index_type value) -> this_type & |
Set the maximum number of iterations. | |
auto | param_search_region () const -> std::pair< scalar_type, scalar_type > |
Get the default region to search for the optimal regularization parameter. | |
auto | regularization_term (const data_type &solution) const -> scalar_type |
Calculate the regularization term. | |
auto | residual_norm (const data_type &solution) const -> scalar_type |
Calculate the squared norm of the residual. | |
auto | residual_norm_rate () const -> scalar_type |
Get the rate of the last residual norm. | |
auto | tol_update_rate () const -> scalar_type |
Get the tolerance of update rate of the solution. | |
auto | tol_update_rate (scalar_type value) -> this_type & |
Set the tolerance of update rate of the solution. | |
auto | update_rate () const noexcept -> scalar_type |
Get the rate of the norm of the update of the solution in the last iteration. | |
Public Member Functions inherited from num_collect::regularization::iterative_regularized_solver_base< tgv2_admm< Coeff, DerivativeMatrix, Data >, Data > | |
void | init (const scalar_type ¶m, data_type &solution) |
Initialize. | |
auto | is_stop_criteria_satisfied (const data_type &solution) const -> bool |
Determine if stopping criteria of the algorithm are satisfied. | |
void | iterate (const scalar_type ¶m, data_type &solution) |
Iterate the algorithm once. | |
void | solve (const scalar_type ¶m, data_type &solution) |
Solve for a regularization parameter. | |
Public Member Functions inherited from num_collect::regularization::implicit_regularized_solver_base< tgv2_admm< Coeff, DerivativeMatrix, Data >, Data > | |
void | calculate_data_for (const data_type &solution, data_type &data) const |
Calculate data for a solution. | |
void | change_data (const data_type &data) |
Change data. | |
auto | regularization_term (const data_type &solution) const -> scalar_type |
Calculate the regularization term. | |
auto | residual_norm (const data_type &solution) const -> scalar_type |
Calculate the squared norm of the residual. | |
Public Member Functions inherited from num_collect::regularization::regularized_solver_base< tgv2_admm< Coeff, DerivativeMatrix, Data >, Data > | |
auto | data_size () const -> index_type |
Get the size of data. | |
auto | param_search_region () const -> std::pair< scalar_type, scalar_type > |
Get the default region to search for the optimal regularization parameter. | |
void | solve (const scalar_type ¶m, data_type &solution) |
Solve for a regularization parameter. | |
Public Member Functions inherited from num_collect::logging::logging_mixin | |
logging_mixin (log_tag_view tag) | |
Constructor. | |
template<typename Child> | |
void | configure_child_algorithm_logger_if_exists (Child &child) |
Configure a logger of a child algorithm if exists. | |
auto | logger () const noexcept -> const num_collect::logging::logger & |
Access to the logger. | |
auto | logger () noexcept -> num_collect::logging::logger & |
Access to the logger. | |
Public Member Functions inherited from num_collect::logging::iterations::iteration_logger_mixin< tgv2_admm< Coeff, DerivativeMatrix, Data > > | |
iteration_logger_mixin ()=default | |
Constructor. | |
void | configure_iteration_logger (num_collect::logging::iterations::iteration_logger< tgv2_admm< Coeff, DerivativeMatrix, Data > > &iteration_logger) const |
Configure an iteration logger. | |
auto | initialize_iteration_logger () -> num_collect::logging::iterations::iteration_logger< tgv2_admm< Coeff, DerivativeMatrix, Data > > & |
Get the iteration logger. |
Private Member Functions | |
void | update_p (const scalar_type ¶m, const data_type &solution) |
Update p_. | |
void | update_s (const scalar_type ¶m, const data_type &solution) |
Update s_. | |
void | update_solution (const scalar_type ¶m, data_type &solution) |
Update the solution. | |
void | update_t (const scalar_type ¶m, const data_type &solution) |
Update t_. | |
void | update_u (const scalar_type ¶m, const data_type &solution) |
Update u_. | |
void | update_z (const scalar_type ¶m, const data_type &solution) |
Update z_. |
Private Attributes | |
const coeff_type * | coeff_ {nullptr} |
Coefficient matrix to compute data vector. | |
linear::impl::operator_conjugate_gradient< data_type > | conjugate_gradient_solution_ {} |
Conjugate gradient solver for update of the solution. | |
linear::impl::operator_conjugate_gradient< data_type > | conjugate_gradient_z_ {} |
Conjugate gradient solver for update of the 1st order derivative. | |
scalar_type | constraint_coeff_ {default_constraint_coeff} |
Coefficient of the constraint for the derivative. | |
const data_type * | data_ {nullptr} |
Data vector. | |
const derivative_matrix_type * | first_derivative_matrix_ {nullptr} |
Matrix to compute the first order derivative. | |
index_type | iterations_ {} |
Number of iterations. | |
index_type | max_iterations_ {default_max_iterations} |
Maximum number of iterations. | |
data_type | p_ {} |
Lagrange multiplier for the 1st order derivative constraint. | |
data_type | p_update_ {} |
Update of p_. | |
data_type | previous_s_ {} |
Previous value of s_. | |
data_type | previous_solution_ {} |
Previous solution. | |
data_type | previous_t_ {} |
Previous value of t_. | |
data_type | previous_z_ {} |
Previous value of z_. | |
scalar_type | rate_of_cg_tol_rate_to_tol_update_rate_ |
Rate of the ratio of the rate of tolerance in CG method to the tolerance of update rate of the solution. | |
data_type | residual_ {} |
Residual vector. | |
data_type | s_ {} |
Component of the 1st order derivative for regularization. | |
const derivative_matrix_type * | second_derivative_matrix_ {nullptr} |
Matrix to compute the second order derivative. | |
scalar_type | second_derivative_ratio_ {default_second_derivative_ratio} |
Ratio of regularization term for the 2nd order derivative. | |
data_type | t_ {} |
2nd order derivative of the solution. | |
data_type | temp_data_ {} |
Temporary vector with the size of data. | |
data_type | temp_solution_ {} |
Temporary vector for the update of the solution. | |
data_type | temp_t_ {} |
Temporary vector for the update of the 2nd order derivative. | |
data_type | temp_z_ {} |
Temporary vector for the update of the 1st order derivative. | |
scalar_type | tol_update_rate_ {default_tol_update_rate} |
Tolerance of update rate of the solution. | |
data_type | u_ {} |
Lagrange multiplier for the 2nd order derivative constraint. | |
data_type | u_update_ {} |
Update of u_. | |
scalar_type | update_rate_ {} |
Rate of norm of the update of the variables in the last iteration. | |
data_type | z_ {} |
Component of the 1st order derivative. |
Static Private Attributes | |
static constexpr auto | default_constraint_coeff |
Default coefficient of the constraint for the derivative. | |
static constexpr index_type | default_max_iterations = 100000 |
Default maximum number of iterations. | |
static constexpr auto | default_rate_of_cg_tol_rate_to_tol_update_rate |
Default value of the ratio of the rate of tolerance in CG method to the tolerance of update rate of the solution. | |
static constexpr auto | default_second_derivative_ratio |
Default ratio of regularization term for the 2nd order derivative. | |
static constexpr auto | default_tol_update_rate |
Default tolerance of update rate of the solution. |
Additional Inherited Members | |
Protected Member Functions inherited from num_collect::regularization::iterative_regularized_solver_base< tgv2_admm< Coeff, DerivativeMatrix, Data >, Data > | |
iterative_regularized_solver_base (logging::log_tag_view tag) | |
Constructor. | |
Protected Member Functions inherited from num_collect::regularization::implicit_regularized_solver_base< tgv2_admm< Coeff, DerivativeMatrix, Data >, Data > | |
implicit_regularized_solver_base (logging::log_tag_view tag) | |
Constructor. | |
auto | derived () noexcept -> tgv2_admm< Coeff, DerivativeMatrix, Data > & |
Access derived object. | |
Protected Member Functions inherited from num_collect::regularization::regularized_solver_base< tgv2_admm< Coeff, DerivativeMatrix, Data >, Data > | |
regularized_solver_base (logging::log_tag_view tag) | |
Constructor. | |
auto | derived () noexcept -> tgv2_admm< Coeff, DerivativeMatrix, Data > & |
Access derived object. |
Class to solve linear equations with 2nd order total generalized variation (TGV) regularization [3] using the alternating direction method of multipliers (ADMM) [2].
Coeff | Type of coefficient matrices. |
DerivativeMatrix | Type of matrices to compute derivatives. |
Data | Type of data vectors. |
This class minimizes the following evaluation function [20] :
\[E(\boldsymbol{x}, \boldsymbol{z}) = \| A \boldsymbol{x} - \boldsymbol{y} \|_2^2 + \lambda (\| D \boldsymbol{x} - \boldsymbol{z} \|_1 + \alpha \| E \boldsymbol{z} \|_1) \]
where variables are defined as follows:
For solving the above problem using ADMM, this class solves the following optimization problem:
\[\begin{aligned} \text{minimize} \hspace{1em} & \|A \boldsymbol{x} - \boldsymbol{y}\|_2^2 + \lambda \left( \|\boldsymbol{s}\|_1 + \alpha \|\boldsymbol{t}\|_1 \right) \\ \text{s.t.} \hspace{1em} & D \boldsymbol{x} - \boldsymbol{z} = \boldsymbol{s} \\ & E \boldsymbol{z} = \boldsymbol{t} \end{aligned} \]
where additional variables are defined as follows:
This class uses the following augmented Lagrangian function:
\[ L_{\rho}(\boldsymbol{x}, \boldsymbol{z}, \boldsymbol{s}, \boldsymbol{t}, \boldsymbol{p}, \boldsymbol{u}) \equiv \|A \boldsymbol{x} - \boldsymbol{y}\|_2^2 + \lambda \left( \|\boldsymbol{s}\|_1 + \alpha \|\boldsymbol{t}\|_1 \right) + \boldsymbol{p}^\top (D \boldsymbol{x} - \boldsymbol{z} - \boldsymbol{s}) + \boldsymbol{u}^\top (E \boldsymbol{z} - \boldsymbol{t}) + \frac{\rho}{2} \|D \boldsymbol{x} - \boldsymbol{z} - \boldsymbol{s}\|_2^2 + \frac{\rho}{2} \|E \boldsymbol{z} - \boldsymbol{t}\|_2^2 \]
where additional variables are defined as follows:
This class uses the following update rules:
\[\begin{aligned} \boldsymbol{x}_{k+1} & = (2 A^\top A + \rho D^\top D)^{-1} \left( 2 A^\top \boldsymbol{y} - D^\top \boldsymbol{p}_k + \rho D^\top \boldsymbol{z}_k + \rho D^\top \boldsymbol{s}_k \right) \\ \boldsymbol{z}_{k+1} & = (\rho I + \rho E^\top E)^{-1} \left( \boldsymbol{p}_k - E^\top \boldsymbol{u}_k + \rho D \boldsymbol{x}_{k+1} - \rho \boldsymbol{s}_k + \rho E^\top \boldsymbol{t}_k \right) \\ \boldsymbol{s}_{k+1} & = \mathcal{T}_{\lambda/\rho} \left( D \boldsymbol{x}_{k+1} - \boldsymbol{z}_{k+1} + \frac{\boldsymbol{p}_k}{\rho} \right) \\ \boldsymbol{t}_{k+1} & = \mathcal{T}_{\lambda/\rho \alpha} \left( E \boldsymbol{z}_{k+1} + \frac{\boldsymbol{u}_k}{\rho} \right) \\ \boldsymbol{p}_{k+1} & = \boldsymbol{p}_k + \rho (D \boldsymbol{x}_{k+1} - \boldsymbol{z}_{k+1} - \boldsymbol{s}_{k+1}) \\ \boldsymbol{u}_{k+1} & = \boldsymbol{u}_k + \rho (E \boldsymbol{z}_{k+1} - \boldsymbol{t}_{k+1}) \end{aligned} \]
where \(\mathcal{T}_\lambda\) is the thresholding operator implemented in num_collect::regularization::impl::apply_shrinkage_operator function.
This class uses the following variable names:
In Formula | In C++ |
---|---|
\(\boldsymbol{x}\) | solution |
\(\boldsymbol{y}\) | data |
\(\boldsymbol{z}\) | z |
\(\boldsymbol{s}\) | s |
\(\boldsymbol{t}\) | t |
\(\boldsymbol{p}\) | p |
\(\boldsymbol{u}\) | u |
\(A\) | coeff |
\(D\) | first_derivative_matrix |
\(E\) | second_derivative_matrix |
\(\lambda\) | param |
\(\alpha\) | second_derivative_ratio |
\(\rho\) | constraint_coeff |
Definition at line 179 of file tgv2_admm.h.
using num_collect::regularization::tgv2_admm< Coeff, DerivativeMatrix, Data >::base_type = iterative_regularized_solver_base<this_type, Data> |
Type of the base class.
Definition at line 186 of file tgv2_admm.h.
using num_collect::regularization::tgv2_admm< Coeff, DerivativeMatrix, Data >::coeff_type = Coeff |
Type of coefficient matrices.
Definition at line 192 of file tgv2_admm.h.
using num_collect::regularization::tgv2_admm< Coeff, DerivativeMatrix, Data >::derivative_matrix_type = DerivativeMatrix |
Type of matrices to compute derivatives.
Definition at line 195 of file tgv2_admm.h.
using num_collect::regularization::tgv2_admm< Coeff, DerivativeMatrix, Data >::this_type = tgv2_admm<Coeff, DerivativeMatrix, Data> |
This type.
Definition at line 183 of file tgv2_admm.h.
|
inline |
Constructor.
Definition at line 200 of file tgv2_admm.h.
|
inline |
Calculate data for a solution.
[in] | solution | Solution. |
[out] | data | Data. |
Definition at line 345 of file tgv2_admm.h.
|
inline |
|
inline |
Compute internal parameters.
[in] | coeff | Coefficient matrix. |
[in] | first_derivative_matrix | Matrix to compute the first derivative. |
[in] | second_derivative_matrix | Matrix to compute the second derivative from the first order derivative. |
[in] | data | Data vector. |
Definition at line 220 of file tgv2_admm.h.
|
inline |
Configure an iteration logger.
[in] | iteration_logger | Iteration logger. |
Definition at line 318 of file tgv2_admm.h.
|
inlinenodiscard |
|
inline |
Initialize.
[in] | param | Regularization parameter. |
[in,out] | solution | Solution. (Assumed to be the initial solution.) |
Definition at line 231 of file tgv2_admm.h.
|
inlinenodiscard |
Determine if stopping criteria of the algorithm are satisfied.
[in] | solution | Solution. |
Definition at line 310 of file tgv2_admm.h.
|
inline |
Iterate the algorithm once.
[in] | param | Regularization parameter. |
[in,out] | solution | Solution. (Assumed to be the last solution.) |
Definition at line 296 of file tgv2_admm.h.
|
inlinenodiscardnoexcept |
|
inlinenodiscard |
|
inline |
Set the maximum number of iterations.
[in] | value | Value. |
Definition at line 414 of file tgv2_admm.h.
|
inlinenodiscard |
Get the default region to search for the optimal regularization parameter.
Definition at line 353 of file tgv2_admm.h.
|
inlinenodiscard |
Calculate the regularization term.
[in] | solution | Solution. |
Definition at line 336 of file tgv2_admm.h.
|
inlinenodiscard |
Calculate the squared norm of the residual.
[in] | solution | Solution. |
Definition at line 330 of file tgv2_admm.h.
|
inlinenodiscard |
Get the rate of the last residual norm.
Definition at line 395 of file tgv2_admm.h.
|
inlinenodiscard |
Get the tolerance of update rate of the solution.
Definition at line 426 of file tgv2_admm.h.
|
inline |
Set the tolerance of update rate of the solution.
[in] | value | Value. |
Definition at line 436 of file tgv2_admm.h.
|
inlineprivate |
Update p_.
[in] | param | Regularization parameter. |
[in] | solution | Solution vector. |
Definition at line 544 of file tgv2_admm.h.
|
inlinenodiscardnoexcept |
Get the rate of the norm of the update of the solution in the last iteration.
Definition at line 386 of file tgv2_admm.h.
|
inlineprivate |
Update s_.
[in] | param | Regularization parameter. |
[in] | solution | Solution vector. |
Definition at line 511 of file tgv2_admm.h.
|
inlineprivate |
Update the solution.
[in] | param | Regularization parameter. |
[in,out] | solution | Solution vector. |
Definition at line 451 of file tgv2_admm.h.
|
inlineprivate |
Update t_.
[in] | param | Regularization parameter. |
[in] | solution | Solution vector. |
Definition at line 526 of file tgv2_admm.h.
|
inlineprivate |
Update u_.
[in] | param | Regularization parameter. |
[in] | solution | Solution vector. |
Definition at line 561 of file tgv2_admm.h.
|
inlineprivate |
Update z_.
[in] | param | Regularization parameter. |
[in] | solution | Solution vector. |
Definition at line 482 of file tgv2_admm.h.
|
private |
Coefficient matrix to compute data vector.
Definition at line 574 of file tgv2_admm.h.
|
private |
Conjugate gradient solver for update of the solution.
Definition at line 641 of file tgv2_admm.h.
|
private |
Conjugate gradient solver for update of the 1st order derivative.
Definition at line 645 of file tgv2_admm.h.
|
private |
Coefficient of the constraint for the derivative.
Definition at line 659 of file tgv2_admm.h.
|
private |
Data vector.
Definition at line 583 of file tgv2_admm.h.
|
staticconstexprprivate |
Default coefficient of the constraint for the derivative.
Definition at line 655 of file tgv2_admm.h.
|
staticconstexprprivate |
Default maximum number of iterations.
Definition at line 662 of file tgv2_admm.h.
|
staticconstexprprivate |
Default value of the ratio of the rate of tolerance in CG method to the tolerance of update rate of the solution.
Definition at line 675 of file tgv2_admm.h.
|
staticconstexprprivate |
Default ratio of regularization term for the 2nd order derivative.
Definition at line 648 of file tgv2_admm.h.
|
staticconstexprprivate |
Default tolerance of update rate of the solution.
Definition at line 668 of file tgv2_admm.h.
|
private |
Matrix to compute the first order derivative.
Definition at line 577 of file tgv2_admm.h.
|
private |
Number of iterations.
Definition at line 586 of file tgv2_admm.h.
|
private |
Maximum number of iterations.
Definition at line 665 of file tgv2_admm.h.
|
private |
Lagrange multiplier for the 1st order derivative constraint.
Definition at line 598 of file tgv2_admm.h.
|
private |
Update of p_.
Definition at line 628 of file tgv2_admm.h.
|
private |
Previous value of s_.
Definition at line 622 of file tgv2_admm.h.
|
private |
Previous solution.
Definition at line 616 of file tgv2_admm.h.
|
private |
Previous value of t_.
Definition at line 625 of file tgv2_admm.h.
|
private |
Previous value of z_.
Definition at line 619 of file tgv2_admm.h.
|
private |
Rate of the ratio of the rate of tolerance in CG method to the tolerance of update rate of the solution.
Definition at line 679 of file tgv2_admm.h.
|
private |
Residual vector.
Definition at line 634 of file tgv2_admm.h.
|
private |
Component of the 1st order derivative for regularization.
Definition at line 592 of file tgv2_admm.h.
|
private |
Matrix to compute the second order derivative.
Definition at line 580 of file tgv2_admm.h.
|
private |
Ratio of regularization term for the 2nd order derivative.
Definition at line 652 of file tgv2_admm.h.
|
private |
2nd order derivative of the solution.
Definition at line 595 of file tgv2_admm.h.
|
private |
Temporary vector with the size of data.
Definition at line 607 of file tgv2_admm.h.
|
private |
Temporary vector for the update of the solution.
Definition at line 604 of file tgv2_admm.h.
|
private |
Temporary vector for the update of the 2nd order derivative.
Definition at line 613 of file tgv2_admm.h.
|
private |
Temporary vector for the update of the 1st order derivative.
Definition at line 610 of file tgv2_admm.h.
|
private |
Tolerance of update rate of the solution.
Definition at line 672 of file tgv2_admm.h.
|
private |
Lagrange multiplier for the 2nd order derivative constraint.
Definition at line 601 of file tgv2_admm.h.
|
private |
Update of u_.
Definition at line 631 of file tgv2_admm.h.
|
private |
Rate of norm of the update of the variables in the last iteration.
Definition at line 637 of file tgv2_admm.h.
|
private |
Component of the 1st order derivative.
Definition at line 589 of file tgv2_admm.h.