The source files for all examples can be found in /examples.

Huber Fitting

This example sets up a $\ell_1$-regularized huber fitting problem using the MLSolver interface provided by GeNIOS. Huber fitting is a form of 'robust regression' that is less sensitive to outliers.

The huber loss function is defined as

\[f^\mathrm{hub}(x) = \begin{cases} \frac{1}{2}x^2 & \lvert x\rvert \leq 1 \\ |x| - \frac{1}{2} & \lvert x \rvert > 1 \end{cases}\]

We want to solve the problem

\[\begin{array}{ll} \text{minimize} & \sum_{i=1}^N f^\mathrm{hub}(a_i^T x - b_i) + \gamma \|x\|_1 \end{array}\]

using GeNIOS
using Random, LinearAlgebra, SparseArrays

Generating the problem data

Random.seed!(1)
N, n = 200, 400
A = randn(N, n)
A .-= sum(A, dims=1) ./ N
normalize!.(eachcol(A))
xstar = sprandn(n, 0.1)
b = A*xstar + 1e-3*randn(N)

# Add outliers
b += 10*collect(sprand(N, 0.05))
λ = 0.05*norm(A'*b, Inf)
0.15848896526224718

MLSolver interface

We just need to specify $f$ and the regularization parameters.

# Huber problem: min ∑ fʰᵘᵇ(aᵢᵀx - bᵢ) + λ||x||₁
f(x) = abs(x) <= 1 ? 0.5*x^2 : abs(x) - 0.5
λ1 = λ
λ2 = 0.0
solver = GeNIOS.MLSolver(f, λ1, λ2, A, b)
res = solve!(solver; options=GeNIOS.SolverOptions(use_dual_gap=false))
rmse = sqrt(1/N*norm(A*solver.zk - b, 2)^2)
println("Final RMSE: $(round(rmse, digits=8))")
Starting setup...
Setup in  0.036s

──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
    Iteration      Objective           RMSE       r_primal         r_dual              ρ           Time
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
            0      4.864e+01      4.932e-01            Inf            Inf      1.000e+00         0.000
            1      4.864e+01      4.932e-01      0.000e+00      1.083e+01      1.000e+00         0.476
           20      2.892e+01      2.713e-01      4.607e-02      2.036e-01      1.000e+00         0.481
           40      2.869e+01      2.516e-01      1.362e-02      7.071e-02      1.000e+00         0.487
           60      2.864e+01      2.465e-01      5.991e-03      4.049e-02      1.000e+00         0.493
           80      2.863e+01      2.439e-01      3.404e-03      2.073e-02      1.000e+00         0.500
          100      2.862e+01      2.426e-01      2.089e-03      1.358e-02      1.000e+00         0.507
          120      2.862e+01      2.419e-01      1.405e-03      8.822e-03      1.000e+00         0.514
          140      2.862e+01      2.416e-01      8.922e-04      5.881e-03      1.000e+00         0.521
          160      2.862e+01      2.413e-01      6.359e-04      3.967e-03      1.000e+00         0.528
          180      2.862e+01      2.412e-01      4.102e-04      2.655e-03      1.000e+00         0.536
          188      2.862e+01      2.411e-01      5.451e-04      2.200e-03      1.000e+00         0.539

SOLVED in  0.539s, 188 iterations
Total time:  0.576s
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────

Final RMSE: 0.4456551

Duality gap

However, we can supply the conjugate function if we want to use the duality gap.

fconj(y) = y > 1 ? Inf : y^2/2.0
solver = GeNIOS.MLSolver(f, λ1, λ2, A, b; fconj=fconj)
res = solve!(solver; options=GeNIOS.SolverOptions(use_dual_gap=true))
rmse = sqrt(1/N*norm(A*solver.zk - b, 2)^2)
println("Final RMSE: $(round(rmse, digits=8))")
Starting setup...
Setup in  0.035s

──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
    Iteration      Objective           RMSE       Dual Gap       r_primal         r_dual              ρ           Time
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
            0      4.864e+01      4.932e-01      1.058e+01            Inf            Inf      1.000e+00         0.000
            1      4.864e+01      4.932e-01      1.058e+01      0.000e+00      1.083e+01      1.000e+00         0.000
           20      2.892e+01      2.713e-01      1.695e-01      4.607e-02      2.036e-01      1.000e+00         0.005
           40      2.869e+01      2.516e-01      4.753e-02      1.362e-02      7.071e-02      1.000e+00         0.010
           60      2.864e+01      2.465e-01      3.121e-02      5.991e-03      4.049e-02      1.000e+00         0.017
           80      2.863e+01      2.439e-01      2.058e-02      3.404e-03      2.073e-02      1.000e+00         0.023
          100      2.862e+01      2.426e-01      1.389e-02      2.089e-03      1.358e-02      1.000e+00         0.030
          120      2.862e+01      2.419e-01      9.869e-03      1.405e-03      8.822e-03      1.000e+00         0.037
          140      2.862e+01      2.416e-01      6.505e-03      8.922e-04      5.881e-03      1.000e+00         0.045
          160      2.862e+01      2.413e-01      4.483e-03      6.359e-04      3.967e-03      1.000e+00         0.053
          180      2.862e+01      2.412e-01      2.848e-03      4.102e-04      2.655e-03      1.000e+00         0.061
          200      2.862e+01      2.411e-01      1.685e-03      3.832e-04      1.568e-03      1.000e+00         0.068
          220      2.862e+01      2.410e-01      9.638e-04      2.093e-04      8.717e-04      1.000e+00         0.076
          240      2.862e+01      2.409e-01      5.632e-04      1.162e-04      4.910e-04      1.000e+00         0.083
          260      2.862e+01      2.409e-01      3.328e-04      6.536e-05      2.798e-04      1.000e+00         0.090
          280      2.862e+01      2.409e-01      1.981e-04      3.716e-05      1.612e-04      1.000e+00         0.097
          300      2.862e+01      2.409e-01      1.186e-04      2.135e-05      9.379e-05      1.000e+00         0.103
          307      2.862e+01      2.409e-01      9.918e-05      1.763e-05      7.779e-05      1.000e+00         0.105

SOLVED in  0.105s, 307 iterations
Total time:  0.141s
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────

Final RMSE: 0.44488898

This page was generated using Literate.jl.