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.18184410798147943

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.042s

──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
    Iteration      Objective           RMSE       r_primal         r_dual              ρ           Time
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
            0      5.505e+01      5.246e-01            Inf            Inf      1.000e+00         0.000
            1      5.505e+01      5.246e-01      0.000e+00      9.994e+00      1.000e+00         0.288
           20      3.755e+01      3.297e-01      5.352e-02      2.085e-01      1.000e+00         0.293
           40      3.735e+01      3.162e-01      1.053e-02      6.937e-02      1.000e+00         0.297
           60      3.730e+01      3.104e-01      7.514e-03      4.175e-02      1.000e+00         0.302
           80      3.728e+01      3.075e-01      3.500e-03      2.238e-02      1.000e+00         0.308
          100      3.728e+01      3.058e-01      2.184e-03      1.384e-02      1.000e+00         0.314
          120      3.728e+01      3.047e-01      1.377e-03      9.316e-03      1.000e+00         0.320
          140      3.728e+01      3.041e-01      8.090e-04      6.417e-03      1.000e+00         0.327
          160      3.728e+01      3.036e-01      6.316e-04      4.580e-03      1.000e+00         0.333
          180      3.728e+01      3.033e-01      4.291e-04      3.225e-03      1.000e+00         0.339
          199      3.728e+01      3.031e-01      4.722e-04      2.222e-03      1.000e+00         0.345

SOLVED in  0.345s, 199 iterations
Total time:  0.387s
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────

Final RMSE: 0.62386089

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.041s

──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
    Iteration      Objective           RMSE       Dual Gap       r_primal         r_dual              ρ           Time
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────
            0      5.505e+01      5.246e-01      6.266e+00            Inf            Inf      1.000e+00         0.000
            1      5.505e+01      5.246e-01      6.266e+00      0.000e+00      9.994e+00      1.000e+00         0.000
           20      3.755e+01      3.297e-01      1.529e-01      5.352e-02      2.085e-01      1.000e+00         0.011
           40      3.735e+01      3.162e-01      6.772e-02      1.053e-02      6.937e-02      1.000e+00         0.019
           60      3.730e+01      3.104e-01      3.585e-02      7.514e-03      4.175e-02      1.000e+00         0.025
           80      3.728e+01      3.075e-01      1.983e-02      3.500e-03      2.238e-02      1.000e+00         0.031
          100      3.728e+01      3.058e-01      1.283e-02      2.184e-03      1.384e-02      1.000e+00         0.037
          120      3.728e+01      3.047e-01      9.128e-03      1.377e-03      9.316e-03      1.000e+00         0.044
          140      3.728e+01      3.041e-01      6.302e-03      8.090e-04      6.417e-03      1.000e+00         0.050
          160      3.728e+01      3.036e-01      4.532e-03      6.316e-04      4.580e-03      1.000e+00         0.057
          180      3.728e+01      3.033e-01      3.218e-03      4.291e-04      3.225e-03      1.000e+00         0.064
          200      3.728e+01      3.031e-01      2.172e-03      2.871e-04      2.185e-03      1.000e+00         0.071
          220      3.728e+01      3.028e-01      1.272e-03      2.688e-04      1.277e-03      1.000e+00         0.078
          240      3.728e+01      3.027e-01      7.490e-04      1.567e-04      7.492e-04      1.000e+00         0.085
          260      3.728e+01      3.026e-01      4.418e-04      9.170e-05      4.402e-04      1.000e+00         0.091
          280      3.728e+01      3.026e-01      2.607e-04      5.375e-05      2.590e-04      1.000e+00         0.098
          300      3.728e+01      3.026e-01      1.538e-04      3.154e-05      1.525e-04      1.000e+00         0.104
          317      3.728e+01      3.026e-01      9.821e-05      2.007e-05      9.730e-05      1.000e+00         0.108

SOLVED in  0.108s, 317 iterations
Total time:  0.150s
──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────────

Final RMSE: 0.62163896

This page was generated using Literate.jl.