scirs2-optimize 0.4.2

Optimization module for SciRS2 (scirs2-optimize)
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
//! Core implicit differentiation engine for optimization layers.
//!
//! Given an optimal solution z* satisfying KKT conditions F(z*, θ) = 0,
//! the implicit function theorem yields:
//!
//!   dz*/dθ = -(∂F/∂z)⁻¹ · (∂F/∂θ)
//!
//! This module builds the KKT Jacobian ∂F/∂z, solves the resulting linear
//! system, and supports an active-set variant that restricts differentiation
//! to active inequality constraints.

use crate::error::{OptimizeError, OptimizeResult};

// ─────────────────────────────────────────────────────────────────────────────
// KKT Jacobian construction
// ─────────────────────────────────────────────────────────────────────────────

/// Build the KKT Jacobian matrix ∂F/∂z for the QP:
///
///   min  ½ x'Qx + c'x
///   s.t. Gx ≤ h   (m inequalities)
///        Ax = b    (p equalities)
///
/// The KKT conditions (with slacks absorbed into complementarity) are:
///
///   F₁ = Qx + c + G'diag(λ) + A'ν  = 0   (stationarity, n eqs)
///   F₂ = diag(λ)(Gx - h)           = 0   (complementarity, m eqs)
///   F₃ = Ax - b                      = 0   (primal equality, p eqs)
///
/// The Jacobian w.r.t. z = (x, λ, ν) is the (n+m+p) × (n+m+p) matrix:
///
///   ┌ Q          G'diag(λ)?   A' ┐
///   │ diag(λ)G   diag(Gx-h)  0  │
///   └ A          0            0  ┘
///
/// For the complementarity row the correct linearisation is:
///   ∂F₂/∂x = diag(λ) G
///   ∂F₂/∂λ = diag(Gx - h)   (= diag(s) where s = h - Gx ≥ 0 at optimum)
///
/// # Arguments
/// * `q` – n*n cost matrix (row-major `Vec<Vec<f64>>`).
/// * `g` – m×n inequality constraint matrix.
/// * `a` – p×n equality constraint matrix.
/// * `x` – optimal primal (length n).
/// * `lam` – optimal inequality duals (length m).
/// * `nu` – optimal equality duals (length p).
///
/// # Returns
/// The full KKT Jacobian as a flat (n+m+p) × (n+m+p) row-major matrix.
pub fn compute_kkt_jacobian(
    q: &[Vec<f64>],
    g: &[Vec<f64>],
    a: &[Vec<f64>],
    x: &[f64],
    lam: &[f64],
    nu: &[f64],
) -> Vec<Vec<f64>> {
    let n = x.len();
    let m = lam.len();
    let p = nu.len();
    let dim = n + m + p;

    let mut jac = vec![vec![0.0; dim]; dim];

    // ── Block (0,0): Q  (n×n) ──────────────────────────────────────────
    for i in 0..n {
        for j in 0..n {
            jac[i][j] = if i < q.len() && j < q[i].len() {
                q[i][j]
            } else {
                0.0
            };
        }
    }

    // ── Block (0,1): G' diag(λ)  →  actually for stationarity the
    //     derivative ∂F₁/∂λ_j = G_j (the j-th row of G, transposed).
    //     But we keep the simpler form: column j of block = G[j][:] (transpose).
    //     Stationarity is Qx + c + G'λ + A'ν = 0, so ∂/∂λ_j = G[j][:] transposed.
    for j in 0..m {
        for i in 0..n {
            let g_val = if j < g.len() && i < g[j].len() {
                g[j][i]
            } else {
                0.0
            };
            jac[i][n + j] = g_val;
        }
    }

    // ── Block (0,2): A'  (n×p) ─────────────────────────────────────────
    for j in 0..p {
        for i in 0..n {
            let a_val = if j < a.len() && i < a[j].len() {
                a[j][i]
            } else {
                0.0
            };
            jac[i][n + m + j] = a_val;
        }
    }

    // ── Block (1,0): diag(λ) G  (m×n) ─────────────────────────────────
    for i in 0..m {
        let li = lam[i];
        for j in 0..n {
            let g_val = if i < g.len() && j < g[i].len() {
                g[i][j]
            } else {
                0.0
            };
            jac[n + i][j] = li * g_val;
        }
    }

    // ── Block (1,1): diag(Gx - h)  (m×m) ──────────────────────────────
    // s_i = (Gx - h)_i
    for i in 0..m {
        let mut gx_i = 0.0;
        if i < g.len() {
            for j in 0..n.min(g[i].len()) {
                gx_i += g[i][j] * x[j];
            }
        }
        // Note: h is not passed here; it cancels with the complementarity form.
        // At optimality λ_i (Gx-h)_i = 0, but the Jacobian entry is (Gx-h)_i.
        // We store the slack as-is; caller must account for h offset.
        jac[n + i][n + i] = gx_i; // will be adjusted by caller with -h_i
    }

    // ── Block (2,0): A  (p×n) ──────────────────────────────────────────
    for i in 0..p {
        for j in 0..n {
            let a_val = if i < a.len() && j < a[i].len() {
                a[i][j]
            } else {
                0.0
            };
            jac[n + m + i][j] = a_val;
        }
    }

    // Blocks (1,2), (2,1), (2,2) are zero.
    jac
}

/// Adjust the complementarity diagonal block with the rhs h.
///
/// Call after `compute_kkt_jacobian` to set diag(Gx - h) correctly.
pub fn adjust_complementarity_diagonal(jac: &mut [Vec<f64>], h: &[f64], n: usize) {
    for (i, &h_i) in h.iter().enumerate() {
        jac[n + i][n + i] -= h_i;
    }
}

// ─────────────────────────────────────────────────────────────────────────────
// Linear system solver (Gaussian elimination with partial pivoting)
// ─────────────────────────────────────────────────────────────────────────────

/// Solve the linear system `A x = rhs` via Gaussian elimination with partial
/// pivoting.
///
/// Both `mat` (square, row-major) and `rhs` are consumed / mutated.
///
/// # Errors
/// Returns `OptimizeError::ComputationError` if the matrix is singular.
pub fn solve_implicit_system(mat: &[Vec<f64>], rhs: &[f64]) -> OptimizeResult<Vec<f64>> {
    let n = rhs.len();
    if mat.len() != n {
        return Err(OptimizeError::InvalidInput(format!(
            "KKT matrix rows ({}) != rhs length ({})",
            mat.len(),
            n
        )));
    }

    // Build augmented matrix
    let mut aug: Vec<Vec<f64>> = mat
        .iter()
        .enumerate()
        .map(|(i, row)| {
            let mut r = row.clone();
            r.push(rhs[i]);
            r
        })
        .collect();

    // Forward elimination with partial pivoting
    for col in 0..n {
        // Find pivot
        let mut max_val = aug[col][col].abs();
        let mut max_row = col;
        for row in (col + 1)..n {
            let v = aug[row][col].abs();
            if v > max_val {
                max_val = v;
                max_row = row;
            }
        }

        if max_val < 1e-30 {
            return Err(OptimizeError::ComputationError(
                "Singular KKT matrix in implicit differentiation".to_string(),
            ));
        }

        if max_row != col {
            aug.swap(col, max_row);
        }

        let pivot = aug[col][col];
        for row in (col + 1)..n {
            let factor = aug[row][col] / pivot;
            for j in col..=n {
                let val = aug[col][j];
                aug[row][j] -= factor * val;
            }
        }
    }

    // Back substitution
    let mut solution = vec![0.0; n];
    for i in (0..n).rev() {
        let mut sum = aug[i][n];
        for j in (i + 1)..n {
            sum -= aug[i][j] * solution[j];
        }
        let diag = aug[i][i];
        if diag.abs() < 1e-30 {
            return Err(OptimizeError::ComputationError(
                "Zero diagonal in back substitution".to_string(),
            ));
        }
        solution[i] = sum / diag;
    }

    Ok(solution)
}

/// Solve the system `mat * X = rhs_matrix` where rhs_matrix has `k` columns.
/// Returns the solution matrix (n × k) in row-major form.
pub fn solve_implicit_system_multi(
    mat: &[Vec<f64>],
    rhs_cols: &[Vec<f64>],
) -> OptimizeResult<Vec<Vec<f64>>> {
    rhs_cols
        .iter()
        .map(|rhs| solve_implicit_system(mat, rhs))
        .collect()
}

// ─────────────────────────────────────────────────────────────────────────────
// Active constraint identification
// ─────────────────────────────────────────────────────────────────────────────

/// Identify active inequality constraints: those where h_i - G_i x ≤ tol
/// (i.e., the slack is near zero).
///
/// Returns indices of the active constraints.
pub fn identify_active_constraints(g: &[Vec<f64>], h: &[f64], x: &[f64], tol: f64) -> Vec<usize> {
    let m = h.len();
    let n = x.len();
    let mut active = Vec::new();

    for i in 0..m {
        let mut gx_i = 0.0;
        if i < g.len() {
            for j in 0..n.min(g[i].len()) {
                gx_i += g[i][j] * x[j];
            }
        }
        let slack = h[i] - gx_i; // slack ≥ 0 at feasibility
        if slack.abs() <= tol {
            active.push(i);
        }
    }

    active
}

/// Extract rows from G and h corresponding to the given active indices.
pub fn extract_active_constraints(
    g: &[Vec<f64>],
    h: &[f64],
    active: &[usize],
) -> (Vec<Vec<f64>>, Vec<f64>) {
    let g_active: Vec<Vec<f64>> = active.iter().filter_map(|&i| g.get(i).cloned()).collect();
    let h_active: Vec<f64> = active.iter().filter_map(|&i| h.get(i).copied()).collect();
    (g_active, h_active)
}

// ─────────────────────────────────────────────────────────────────────────────
// Full implicit backward pass
// ─────────────────────────────────────────────────────────────────────────────

/// Compute the full implicit gradient dz*/dθ for a QP.
///
/// Given dl/dx (the upstream gradient from the loss), compute dl/dθ for
/// θ = (Q, c, G, h, A, b) by solving:
///
///   (∂F/∂z)' dz = -[dl/dx, 0, 0]'
///
/// and then computing dl/dθ = (∂F/∂θ)' dz.
///
/// # Arguments
/// * `q` – n×n cost matrix.
/// * `g` – m×n inequality constraint matrix.
/// * `h` – m inequality rhs.
/// * `a` – p×n equality constraint matrix.
/// * `x` – optimal primal solution.
/// * `lam` – optimal inequality duals.
/// * `nu` – optimal equality duals.
/// * `dl_dx` – upstream gradient dl/dx (length n).
///
/// # Returns
/// The implicit gradients for all parameters.
pub fn compute_full_implicit_gradient(
    q: &[Vec<f64>],
    g: &[Vec<f64>],
    h: &[f64],
    a: &[Vec<f64>],
    x: &[f64],
    lam: &[f64],
    nu: &[f64],
    dl_dx: &[f64],
) -> OptimizeResult<super::types::ImplicitGradient> {
    let n = x.len();
    let m = lam.len();
    let p = nu.len();
    let dim = n + m + p;

    // Build the KKT Jacobian
    let mut kkt = compute_kkt_jacobian(q, g, a, x, lam, nu);
    adjust_complementarity_diagonal(&mut kkt, h, n);

    // Transpose the KKT Jacobian (we solve the adjoint system)
    let mut kkt_t = vec![vec![0.0; dim]; dim];
    for i in 0..dim {
        for j in 0..dim {
            kkt_t[i][j] = kkt[j][i];
        }
    }

    // RHS = -[dl/dx, 0_m, 0_p]
    let mut rhs = vec![0.0; dim];
    for i in 0..n {
        rhs[i] = -dl_dx[i];
    }

    // Solve: kkt_t * dz = rhs
    let dz = solve_implicit_system(&kkt_t, &rhs)?;

    // Extract components: dz = (dx, dlam, dnu)
    let dx = &dz[..n];
    let dlam = &dz[n..n + m];
    let dnu = &dz[n + m..];

    // ── Compute dl/dθ from dz ──────────────────────────────────────────
    // dl/dc = dx  (from ∂F₁/∂c = I)
    let dl_dc = dx.to_vec();

    // dl/dh = -dlam  (from ∂F₂/∂h = -diag(λ), contracted with dlam gives -λ·dlam,
    //  but more directly: ∂F₂/∂h_i = -λ_i, so dl/dh_i = -λ_i * dlam_i / λ_i = -dlam_i
    //  when λ_i ≠ 0, and 0 otherwise; we use dlam directly.)
    let dl_dh: Vec<f64> = dlam.iter().map(|&v| -v).collect();

    // dl/db = -dnu  (from ∂F₃/∂b = -I)
    let dl_db: Vec<f64> = dnu.iter().map(|&v| -v).collect();

    // dl/dQ = dx * x' (outer product, symmetric part)
    let mut dl_dq = vec![vec![0.0; n]; n];
    for i in 0..n {
        for j in 0..n {
            // ∂F₁/∂Q_{ij} · x_j, contracted with dx_i
            // = 0.5 * (dx_i * x_j + dx_j * x_i)  for symmetric Q
            dl_dq[i][j] = 0.5 * (dx[i] * x[j] + dx[j] * x[i]);
        }
    }

    // dl/dG: from stationarity (G'λ term) and complementarity (diag(λ)G term)
    let mut dl_dg = vec![vec![0.0; n]; m];
    for i in 0..m {
        for j in 0..n {
            // From stationarity: ∂F₁_j / ∂G_{ij} = λ_i, contracted with dx_j
            // From complementarity: ∂F₂_i / ∂G_{ij} = λ_i * x_j, contracted with dlam_i
            dl_dg[i][j] = dx[j] * lam[i] + dlam[i] * lam[i] * x[j];
        }
    }

    // dl/dA: from stationarity (A'ν) and primal equality (Ax-b)
    let mut dl_da = vec![vec![0.0; n]; p];
    for i in 0..p {
        for j in 0..n {
            // From stationarity: ∂F₁_j / ∂A_{ij} = ν_i, contracted with dx_j
            // From primal eq: ∂F₃_i / ∂A_{ij} = x_j, contracted with dnu_i
            dl_da[i][j] = dx[j] * nu[i] + dnu[i] * x[j];
        }
    }

    Ok(super::types::ImplicitGradient {
        dl_dq: Some(dl_dq),
        dl_dc,
        dl_dg: Some(dl_dg),
        dl_dh,
        dl_da: if p > 0 { Some(dl_da) } else { None },
        dl_db,
    })
}

/// Compute implicit gradient using only active constraints.
///
/// This is faster than full differentiation when many inequality constraints
/// are inactive, since those constraints have zero dual variables and do not
/// contribute to the gradient.
pub fn compute_active_set_implicit_gradient(
    q: &[Vec<f64>],
    g: &[Vec<f64>],
    h: &[f64],
    a: &[Vec<f64>],
    x: &[f64],
    lam: &[f64],
    nu: &[f64],
    dl_dx: &[f64],
    active_tol: f64,
) -> OptimizeResult<super::types::ImplicitGradient> {
    let m = lam.len();

    // Identify active constraints
    let active = identify_active_constraints(g, h, x, active_tol);
    let (g_active, h_active) = extract_active_constraints(g, h, &active);
    let lam_active: Vec<f64> = active
        .iter()
        .filter_map(|&i| if i < m { Some(lam[i]) } else { None })
        .collect();

    // Solve reduced system
    let grad =
        compute_full_implicit_gradient(q, &g_active, &h_active, a, x, &lam_active, nu, dl_dx)?;

    // Expand gradient back to full dimension
    let m_full = lam.len();
    let n = x.len();

    let mut dl_dh_full = vec![0.0; m_full];
    for (idx, &ai) in active.iter().enumerate() {
        if ai < m_full && idx < grad.dl_dh.len() {
            dl_dh_full[ai] = grad.dl_dh[idx];
        }
    }

    let dl_dg_full = if let Some(ref dg) = grad.dl_dg {
        let mut full = vec![vec![0.0; n]; m_full];
        for (idx, &ai) in active.iter().enumerate() {
            if ai < m_full && idx < dg.len() {
                full[ai] = dg[idx].clone();
            }
        }
        Some(full)
    } else {
        None
    };

    Ok(super::types::ImplicitGradient {
        dl_dq: grad.dl_dq,
        dl_dc: grad.dl_dc,
        dl_dg: dl_dg_full,
        dl_dh: dl_dh_full,
        dl_da: grad.dl_da,
        dl_db: grad.dl_db,
    })
}

#[cfg(test)]
mod tests {
    use super::*;

    #[test]
    fn test_kkt_jacobian_structure_2x2() {
        // Simple 2-var QP with 1 inequality, 0 equalities
        // Q = [[2, 0], [0, 2]], G = [[1, 1]], h = [1]
        let q = vec![vec![2.0, 0.0], vec![0.0, 2.0]];
        let g = vec![vec![1.0, 1.0]];
        let a: Vec<Vec<f64>> = vec![];
        let x = vec![0.25, 0.25];
        let lam = vec![0.5];
        let nu: Vec<f64> = vec![];

        let jac = compute_kkt_jacobian(&q, &g, &a, &x, &lam, &nu);

        // Dimension: 2 + 1 + 0 = 3
        assert_eq!(jac.len(), 3);
        assert_eq!(jac[0].len(), 3);

        // Block (0,0): Q
        assert!((jac[0][0] - 2.0).abs() < 1e-12);
        assert!((jac[1][1] - 2.0).abs() < 1e-12);

        // Block (0,1): G' (transposed)
        assert!((jac[0][2] - 1.0).abs() < 1e-12);
        assert!((jac[1][2] - 1.0).abs() < 1e-12);

        // Block (1,0): diag(λ) G
        assert!((jac[2][0] - 0.5).abs() < 1e-12);
        assert!((jac[2][1] - 0.5).abs() < 1e-12);
    }

    #[test]
    fn test_active_constraint_identification() {
        let g = vec![vec![1.0, 0.0], vec![0.0, 1.0], vec![1.0, 1.0]];
        let h = vec![1.0, 2.0, 1.5];
        let x = vec![1.0, 0.5]; // Gx = [1.0, 0.5, 1.5], slacks = [0.0, 1.5, 0.0]

        let active = identify_active_constraints(&g, &h, &x, 1e-6);
        assert_eq!(active, vec![0, 2]);
    }

    #[test]
    fn test_solve_implicit_system_simple() {
        // Solve [[2, 1], [1, 3]] x = [5, 7]  →  x = [8/5, 9/5]
        let mat = vec![vec![2.0, 1.0], vec![1.0, 3.0]];
        let rhs = vec![5.0, 7.0];
        let sol = solve_implicit_system(&mat, &rhs).expect("solve failed");
        assert!((sol[0] - 1.6).abs() < 1e-10);
        assert!((sol[1] - 1.8).abs() < 1e-10);
    }

    #[test]
    fn test_solve_singular_matrix() {
        let mat = vec![vec![1.0, 2.0], vec![2.0, 4.0]];
        let rhs = vec![3.0, 6.0];
        let result = solve_implicit_system(&mat, &rhs);
        assert!(result.is_err());
    }

    #[test]
    fn test_extract_active_constraints() {
        let g = vec![vec![1.0], vec![2.0], vec![3.0]];
        let h = vec![10.0, 20.0, 30.0];
        let active = vec![0, 2];

        let (ga, ha) = extract_active_constraints(&g, &h, &active);
        assert_eq!(ga.len(), 2);
        assert!((ga[0][0] - 1.0).abs() < 1e-12);
        assert!((ga[1][0] - 3.0).abs() < 1e-12);
        assert!((ha[0] - 10.0).abs() < 1e-12);
        assert!((ha[1] - 30.0).abs() < 1e-12);
    }

    #[test]
    fn test_empty_constraints() {
        let q = vec![vec![2.0, 0.0], vec![0.0, 2.0]];
        let g: Vec<Vec<f64>> = vec![];
        let a: Vec<Vec<f64>> = vec![];
        let x = vec![1.0, 2.0];
        let lam: Vec<f64> = vec![];
        let nu: Vec<f64> = vec![];

        let jac = compute_kkt_jacobian(&q, &g, &a, &x, &lam, &nu);
        assert_eq!(jac.len(), 2);
        assert!((jac[0][0] - 2.0).abs() < 1e-12);
    }
}