alpaqa 1.0.0a11
Nonconvex constrained optimization
Loading...
Searching...
No Matches
nuclear-norm.hpp
Go to the documentation of this file.
1#pragma once
2
4#include <Eigen/SVD>
5
6namespace alpaqa::functions {
7
8template <Config Conf>
9using DefaultSVD = Eigen::BDCSVD<typename Conf::mat,
10 Eigen::ComputeThinU | Eigen::ComputeThinV>;
11
12/// Nuclear norm (ℓ₁-norm of singular values).
13/// @ingroup grp_Functions
14template <Config Conf, class SVD = DefaultSVD<Conf>>
17
18 /// Construct without pre-allocation.
20 if (λ < 0 || !std::isfinite(λ))
21 throw std::invalid_argument("NuclearNorm::λ must be nonnegative");
22 }
23 /// Construct with pre-allocation.
25 : λ{λ}, rows{rows}, cols{cols}, svd{rows, cols},
26 singular_values{std::min(rows, cols)} {
27 if (λ < 0 || !std::isfinite(λ))
28 throw std::invalid_argument("NuclearNorm::λ must be nonnegative");
29 }
30
32 length_t rows = 0, cols = 0;
33 SVD svd;
35
36 real_t prox(crmat in, rmat out, real_t γ = 1) {
37 if (λ == 0) {
38 out = in;
39 return 0;
40 }
41 if (rows == 0 || cols == 0) { // dynamic size
42 assert(in.rows() == out.rows());
43 assert(in.cols() == out.cols());
44 svd.compute(in);
45 } else { // fixed size
46 assert(in.size() == rows * cols);
47 assert(out.size() == rows * cols);
48 svd.compute(in.reshaped(rows, cols));
49 }
50 const length_t n = svd.singularValues().size();
51 auto step = vec::Constant(n, λ * γ);
52 singular_values = vec::Zero(n).cwiseMax(svd.singularValues() - step);
53 real_t value = λ * singular_values.template lpNorm<1>();
54 auto it0 = std::find(singular_values.begin(), singular_values.end(), 0);
55 index_t rank = it0 - singular_values.begin();
56 using Eigen::placeholders::all, Eigen::seqN;
57 auto sel = seqN(0, rank);
58 auto &&U = svd.matrixU(), &&V = svd.matrixV();
59 auto &&U1 = U(all, sel);
60 auto &&Σ1 = singular_values(sel).asDiagonal();
61 auto &&V1T = V.transpose()(sel, all);
62 out.reshaped().noalias() = (U1 * Σ1 * V1T).reshaped();
63 return value;
64 }
65
67 crmat in, rmat out, real_t γ) {
68 return self.prox(std::move(in), std::move(out), γ);
69 }
70};
71
72} // namespace alpaqa::functions
#define USING_ALPAQA_CONFIG(Conf)
Definition: config.hpp:54
Eigen::BDCSVD< typename Conf::mat, Eigen::ComputeThinU|Eigen::ComputeThinV > DefaultSVD
std::decay_t< decltype(Tag)> tag_t
Definition: tag-invoke.hpp:74
typename Conf::crmat crmat
Definition: config.hpp:73
typename Conf::rmat rmat
Definition: config.hpp:72
typename Conf::real_t real_t
Definition: config.hpp:63
typename Conf::index_t index_t
Definition: config.hpp:75
typename Conf::length_t length_t
Definition: config.hpp:74
typename Conf::vec vec
Definition: config.hpp:64
Nuclear norm (ℓ₁-norm of singular values).
real_t prox(crmat in, rmat out, real_t γ=1)
NuclearNorm(real_t λ=1)
Construct without pre-allocation.
NuclearNorm(real_t λ, length_t rows, length_t cols)
Construct with pre-allocation.
friend real_t alpaqa_tag_invoke(tag_t< alpaqa::prox >, NuclearNorm &self, crmat in, rmat out, real_t γ)