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
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
//! Trust-region algorithm for constrained optimization
//!
//! This module implements the trust-region constrained optimization algorithm,
//! supporting both numerical and user-supplied gradients and Hessians.
//!
//! # Example with custom gradient
//!
//! ```
//! use scirs2_optimize::constrained::trust_constr::{
//!     minimize_trust_constr_with_derivatives, HessianUpdate,
//! };
//! use scirs2_optimize::constrained::{Constraint, Options};
//! use scirs2_core::ndarray::{array, Array1};
//!
//! // Objective: f(x) = (x[0] - 1)² + (x[1] - 2.5)²
//! fn objective(x: &[f64]) -> f64 {
//!     (x[0] - 1.0).powi(2) + (x[1] - 2.5).powi(2)
//! }
//!
//! // Analytical gradient: [2(x[0] - 1), 2(x[1] - 2.5)]
//! fn gradient(x: &[f64]) -> Array1<f64> {
//!     array![2.0 * (x[0] - 1.0), 2.0 * (x[1] - 2.5)]
//! }
//!
//! // Constraint: x[0] + x[1] <= 3  =>  g(x) = 3 - x[0] - x[1] >= 0
//! fn constraint(x: &[f64]) -> f64 {
//!     3.0 - x[0] - x[1]
//! }
//!
//! let x0 = array![0.0, 0.0];
//! let constraints = vec![Constraint::new(constraint, Constraint::INEQUALITY)];
//! let options = Options::default();
//!
//! let result = minimize_trust_constr_with_derivatives(
//!     objective,
//!     &x0,
//!     &constraints,
//!     &options,
//!     Some(gradient),           // Custom gradient
//!     None::<fn(&[f64]) -> _>,  // No custom Hessian (use BFGS)
//!     HessianUpdate::BFGS,
//! ).expect("valid input");
//!
//! assert!(result.success);
//! ```

use crate::constrained::{Constraint, ConstraintFn, ConstraintKind, Options};
use crate::error::OptimizeResult;
use crate::result::OptimizeResults;
use scirs2_core::ndarray::{Array1, Array2, ArrayBase, Axis, Data, Ix1};

/// Hessian update strategy for quasi-Newton methods
#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
pub enum HessianUpdate {
    /// BFGS update (default) - good for general problems
    #[default]
    BFGS,
    /// SR1 update - can capture negative curvature
    SR1,
    /// Use exact Hessian provided by user
    Exact,
}

/// Type alias for gradient function
pub type GradientFn = fn(&[f64]) -> Array1<f64>;

/// Type alias for Hessian function
pub type HessianFn = fn(&[f64]) -> Array2<f64>;

#[allow(clippy::many_single_char_names)]
#[allow(dead_code)]
pub fn minimize_trust_constr<F, S>(
    func: F,
    x0: &ArrayBase<S, Ix1>,
    constraints: &[Constraint<ConstraintFn>],
    options: &Options,
) -> OptimizeResult<OptimizeResults<f64>>
where
    F: Fn(&[f64]) -> f64,
    S: Data<Elem = f64>,
{
    // Get options or use defaults
    let ftol = options.ftol.unwrap_or(1e-8);
    let gtol = options.gtol.unwrap_or(1e-8);
    let ctol = options.ctol.unwrap_or(1e-8);
    let maxiter = options.maxiter.unwrap_or(100 * x0.len());
    let eps = options.eps.unwrap_or(1e-8);

    // Initialize variables
    let n = x0.len();
    let mut x = x0.to_owned();
    let mut f = func(x.as_slice().expect("Operation failed"));
    let mut nfev = 1;

    // Initialize the Lagrange multipliers
    let mut lambda = Array1::zeros(constraints.len());

    // Calculate initial gradient using finite differences
    let mut g = Array1::zeros(n);
    for i in 0..n {
        let mut x_h = x.clone();
        x_h[i] += eps;
        let f_h = func(x_h.as_slice().expect("Operation failed"));
        g[i] = (f_h - f) / eps;
        nfev += 1;
    }

    // Evaluate initial constraints
    let mut c = Array1::zeros(constraints.len());
    for (i, constraint) in constraints.iter().enumerate() {
        if !constraint.is_bounds() {
            let val = (constraint.fun)(x.as_slice().expect("Operation failed"));

            match constraint.kind {
                ConstraintKind::Inequality => {
                    c[i] = val; // g(x) >= 0 constraint
                }
                ConstraintKind::Equality => {
                    c[i] = val; // h(x) = 0 constraint
                }
            }
        }
    }

    // Calculate constraint Jacobian
    let mut a = Array2::zeros((constraints.len(), n));
    for (i, constraint) in constraints.iter().enumerate() {
        if !constraint.is_bounds() {
            for j in 0..n {
                let mut x_h = x.clone();
                x_h[j] += eps;
                let c_h = (constraint.fun)(x_h.as_slice().expect("Operation failed"));
                a[[i, j]] = (c_h - c[i]) / eps;
                nfev += 1;
            }
        }
    }

    // Initialize trust region radius
    let mut delta = 1.0;

    // Initialize approximation of the Hessian of the Lagrangian
    let mut b = Array2::eye(n);

    // Main optimization loop
    let mut iter = 0;

    while iter < maxiter {
        // Check constraint violation
        let mut max_constraint_violation = 0.0;
        for (i, &ci) in c.iter().enumerate() {
            match constraints[i].kind {
                ConstraintKind::Inequality => {
                    if ci < -ctol {
                        max_constraint_violation = f64::max(max_constraint_violation, -ci);
                    }
                }
                ConstraintKind::Equality => {
                    // For equality constraints, violation is |h(x)|
                    let violation = ci.abs();
                    if violation > ctol {
                        max_constraint_violation = f64::max(max_constraint_violation, violation);
                    }
                }
            }
        }

        // Check convergence on gradient and constraints
        if g.iter().all(|&gi| gi.abs() < gtol) && max_constraint_violation < ctol {
            break;
        }

        // Compute the Lagrangian gradient
        let mut lag_grad = g.clone();
        for (i, &li) in lambda.iter().enumerate() {
            let include_constraint = match constraints[i].kind {
                ConstraintKind::Inequality => {
                    // For inequality constraints: include if active (lambda > 0) or violated (c < 0)
                    li > 0.0 || c[i] < -ctol
                }
                ConstraintKind::Equality => {
                    // For equality constraints: always include (always active)
                    true
                }
            };

            if include_constraint {
                for j in 0..n {
                    // L = f - sum(lambda_i * c_i) for constraints
                    // So gradient of L is grad(f) - sum(lambda_i * grad(c_i))
                    lag_grad[j] -= li * a[[i, j]];
                }
            }
        }

        // Compute the step using a trust-region approach
        // Solve the constrained quadratic subproblem:
        // min 0.5 * p^T B p + g^T p  s.t. ||p|| <= delta and linearized constraints

        let (p, predicted_reduction) =
            compute_trust_region_step_constrained(&lag_grad, &b, &a, &c, delta, constraints, ctol);

        // Try the step
        let x_new = &x + &p;

        // Evaluate function and constraints at new point
        let f_new = func(x_new.as_slice().expect("Operation failed"));
        nfev += 1;

        let mut c_new = Array1::zeros(constraints.len());
        for (i, constraint) in constraints.iter().enumerate() {
            if !constraint.is_bounds() {
                c_new[i] = (constraint.fun)(x_new.as_slice().expect("Operation failed"));
                nfev += 1;
            }
        }

        // Compute change in merit function (includes constraint violation)
        let mut merit = f;
        let mut merit_new = f_new;

        // Add constraint violation penalty
        let penalty = 10.0; // Simple fixed penalty parameter
        for (i, &ci) in c.iter().enumerate() {
            match constraints[i].kind {
                ConstraintKind::Inequality => {
                    // For inequality constraints: penalize only violations (g(x) < 0)
                    merit += penalty * f64::max(0.0, -ci);
                    merit_new += penalty * f64::max(0.0, -c_new[i]);
                }
                ConstraintKind::Equality => {
                    // For equality constraints: penalize any deviation from zero
                    merit += penalty * ci.abs();
                    merit_new += penalty * c_new[i].abs();
                }
            }
        }

        // Compute actual reduction in merit function
        let actual_reduction = merit - merit_new;

        // Compute ratio of actual to predicted reduction
        let rho = if predicted_reduction > 0.0 {
            actual_reduction / predicted_reduction
        } else {
            0.0
        };

        // Update trust region radius based on the quality of the step
        if rho < 0.25 {
            delta *= 0.5;
        } else if rho > 0.75 && p.iter().map(|&pi| pi * pi).sum::<f64>().sqrt() >= 0.9 * delta {
            delta *= 2.0;
        }

        // Accept or reject the step
        if rho > 0.1 {
            // Accept the step
            x = x_new;
            f = f_new;
            c = c_new;

            // Check convergence on function value
            if (merit - merit_new).abs() < ftol * (1.0 + merit.abs()) {
                break;
            }

            // Compute new gradient
            let mut g_new = Array1::zeros(n);
            for i in 0..n {
                let mut x_h = x.clone();
                x_h[i] += eps;
                let f_h = func(x_h.as_slice().expect("Operation failed"));
                g_new[i] = (f_h - f) / eps;
                nfev += 1;
            }

            // Compute new constraint Jacobian
            let mut a_new = Array2::zeros((constraints.len(), n));
            for (i, constraint) in constraints.iter().enumerate() {
                if !constraint.is_bounds() {
                    for j in 0..n {
                        let mut x_h = x.clone();
                        x_h[j] += eps;
                        let c_h = (constraint.fun)(x_h.as_slice().expect("Operation failed"));
                        a_new[[i, j]] = (c_h - c[i]) / eps;
                        nfev += 1;
                    }
                }
            }

            // Update Lagrange multipliers using projected gradient method
            for (i, constraint) in constraints.iter().enumerate() {
                match constraint.kind {
                    ConstraintKind::Inequality => {
                        if c[i] < -ctol {
                            // Active or violated constraint
                            // Increase multiplier if constraint is violated
                            lambda[i] = f64::max(0.0, lambda[i] - c[i] * penalty);
                        } else {
                            // Decrease multiplier towards zero
                            lambda[i] = f64::max(0.0, lambda[i] - 0.1 * lambda[i]);
                        }
                    }
                    ConstraintKind::Equality => {
                        // For equality constraints, multipliers can be positive or negative
                        // Update based on constraint violation to drive h(x) -> 0
                        let step_size = 0.1;
                        lambda[i] -= step_size * c[i] * penalty;

                        // Optional: add some damping to prevent oscillations
                        lambda[i] *= 0.9;
                    }
                }
            }

            // Update Hessian approximation using BFGS or SR1
            let s = &p;
            let y = &g_new - &g;

            // Simple BFGS update for the Hessian approximation
            let s_dot_y = s.dot(&y);
            if s_dot_y > 1e-10 {
                let s_col = s.clone().insert_axis(Axis(1));
                let s_row = s.clone().insert_axis(Axis(0));

                let bs = b.dot(s);
                let bs_col = bs.clone().insert_axis(Axis(1));
                let bs_row = bs.clone().insert_axis(Axis(0));

                let term1 = s_dot_y + s.dot(&bs);
                let term2 = &s_col.dot(&s_row) * (term1 / (s_dot_y * s_dot_y));

                let term3 = &bs_col.dot(&s_row) + &s_col.dot(&bs_row);

                b = &b + &term2 - &(&term3 / s_dot_y);
            }

            // Update variables for next iteration
            g = g_new;
            a = a_new;
        }

        iter += 1;
    }

    // Prepare constraint values for the result
    let mut c_result = Array1::zeros(constraints.len());
    for (i, constraint) in constraints.iter().enumerate() {
        if !constraint.is_bounds() {
            c_result[i] = c[i];
        }
    }

    // Create and return result
    let mut result = OptimizeResults::default();
    result.x = x;
    result.fun = f;
    result.jac = Some(g.into_raw_vec_and_offset().0);
    result.constr = Some(c_result);
    result.nfev = nfev;
    result.nit = iter;
    result.success = iter < maxiter;

    if result.success {
        result.message = "Optimization terminated successfully.".to_string();
    } else {
        result.message = "Maximum iterations reached.".to_string();
    }

    Ok(result)
}

/// Minimizes a function with constraints using trust-region method with custom derivatives.
///
/// This function extends [`minimize_trust_constr`] by allowing user-supplied gradient
/// and Hessian functions, which can significantly improve performance and accuracy
/// compared to finite-difference approximations.
///
/// # Arguments
///
/// * `func` - The objective function to minimize
/// * `x0` - Initial guess for the optimization
/// * `constraints` - Vector of constraints (equality and/or inequality)
/// * `options` - Optimization options (tolerances, max iterations, etc.)
/// * `jac` - Optional gradient function. If `None`, finite differences are used.
/// * `hess` - Optional Hessian function. Only used when `hess_update` is [`HessianUpdate::Exact`].
/// * `hess_update` - Strategy for Hessian updates (BFGS, SR1, or Exact)
///
/// # Returns
///
/// * `OptimizeResults` containing the optimization solution
///
/// # Example
///
/// ```
/// use scirs2_optimize::constrained::trust_constr::{
///     minimize_trust_constr_with_derivatives, HessianUpdate,
/// };
/// use scirs2_optimize::constrained::{Constraint, Options};
/// use scirs2_core::ndarray::{array, Array1, Array2};
///
/// // Rosenbrock function
/// fn rosenbrock(x: &[f64]) -> f64 {
///     (1.0 - x[0]).powi(2) + 100.0 * (x[1] - x[0].powi(2)).powi(2)
/// }
///
/// // Analytical gradient
/// fn rosenbrock_grad(x: &[f64]) -> Array1<f64> {
///     array![
///         -2.0 * (1.0 - x[0]) - 400.0 * x[0] * (x[1] - x[0].powi(2)),
///         200.0 * (x[1] - x[0].powi(2))
///     ]
/// }
///
/// // Constraint: x[0]^2 + x[1]^2 <= 2
/// fn circle_constraint(x: &[f64]) -> f64 {
///     2.0 - x[0].powi(2) - x[1].powi(2)
/// }
///
/// let x0 = array![0.0, 0.0];
/// let constraints = vec![Constraint::new(circle_constraint, Constraint::INEQUALITY)];
///
/// let result = minimize_trust_constr_with_derivatives(
///     rosenbrock,
///     &x0,
///     &constraints,
///     &Options::default(),
///     Some(rosenbrock_grad),
///     None::<fn(&[f64]) -> Array2<f64>>,
///     HessianUpdate::BFGS,
/// ).expect("valid input");
/// ```
#[allow(clippy::many_single_char_names)]
#[allow(clippy::too_many_arguments)]
#[allow(dead_code)]
pub fn minimize_trust_constr_with_derivatives<F, S, G, H>(
    func: F,
    x0: &ArrayBase<S, Ix1>,
    constraints: &[Constraint<ConstraintFn>],
    options: &Options,
    jac: Option<G>,
    hess: Option<H>,
    hess_update: HessianUpdate,
) -> OptimizeResult<OptimizeResults<f64>>
where
    F: Fn(&[f64]) -> f64,
    S: Data<Elem = f64>,
    G: Fn(&[f64]) -> Array1<f64>,
    H: Fn(&[f64]) -> Array2<f64>,
{
    // Get options or use defaults
    let ftol = options.ftol.unwrap_or(1e-8);
    let gtol = options.gtol.unwrap_or(1e-8);
    let ctol = options.ctol.unwrap_or(1e-8);
    let maxiter = options.maxiter.unwrap_or(100 * x0.len());
    let eps = options.eps.unwrap_or(1e-8);

    // Initialize variables
    let n = x0.len();
    let mut x = x0.to_owned();
    let mut f = func(x.as_slice().expect("Operation failed"));
    let mut nfev = 1;
    let mut njev = 0;
    let mut nhev = 0;

    // Initialize the Lagrange multipliers
    let mut lambda = Array1::zeros(constraints.len());

    // Calculate initial gradient (using custom jac if provided)
    let mut g = if let Some(ref grad_fn) = jac {
        njev += 1;
        grad_fn(x.as_slice().expect("Operation failed"))
    } else {
        // Finite difference gradient
        let mut grad = Array1::zeros(n);
        for i in 0..n {
            let mut x_h = x.clone();
            x_h[i] += eps;
            let f_h = func(x_h.as_slice().expect("Operation failed"));
            grad[i] = (f_h - f) / eps;
            nfev += 1;
        }
        grad
    };

    // Evaluate initial constraints
    let mut c = Array1::zeros(constraints.len());
    for (i, constraint) in constraints.iter().enumerate() {
        if !constraint.is_bounds() {
            let val = (constraint.fun)(x.as_slice().expect("Operation failed"));
            match constraint.kind {
                ConstraintKind::Inequality => c[i] = val,
                ConstraintKind::Equality => c[i] = val,
            }
        }
    }

    // Calculate constraint Jacobian
    let mut a = Array2::zeros((constraints.len(), n));
    for (i, constraint) in constraints.iter().enumerate() {
        if !constraint.is_bounds() {
            for j in 0..n {
                let mut x_h = x.clone();
                x_h[j] += eps;
                let c_h = (constraint.fun)(x_h.as_slice().expect("Operation failed"));
                a[[i, j]] = (c_h - c[i]) / eps;
                nfev += 1;
            }
        }
    }

    // Initialize trust region radius
    let mut delta = 1.0;

    // Initialize Hessian or its approximation
    let mut b = if hess_update == HessianUpdate::Exact {
        if let Some(ref hess_fn) = hess {
            nhev += 1;
            hess_fn(x.as_slice().expect("Operation failed"))
        } else {
            // Fallback to identity if exact requested but no Hessian provided
            Array2::eye(n)
        }
    } else {
        Array2::eye(n)
    };

    // Main optimization loop
    let mut iter = 0;

    while iter < maxiter {
        // Check constraint violation
        let mut max_constraint_violation = 0.0;
        for (i, &ci) in c.iter().enumerate() {
            match constraints[i].kind {
                ConstraintKind::Inequality => {
                    if ci < -ctol {
                        max_constraint_violation = f64::max(max_constraint_violation, -ci);
                    }
                }
                ConstraintKind::Equality => {
                    let violation = ci.abs();
                    if violation > ctol {
                        max_constraint_violation = f64::max(max_constraint_violation, violation);
                    }
                }
            }
        }

        // Check convergence
        if g.iter().all(|&gi| gi.abs() < gtol) && max_constraint_violation < ctol {
            break;
        }

        // Compute the Lagrangian gradient
        let mut lag_grad = g.clone();
        for (i, &li) in lambda.iter().enumerate() {
            let include_constraint = match constraints[i].kind {
                ConstraintKind::Inequality => li > 0.0 || c[i] < -ctol,
                ConstraintKind::Equality => true,
            };

            if include_constraint {
                for j in 0..n {
                    lag_grad[j] -= li * a[[i, j]];
                }
            }
        }

        // Compute the step
        let (p, predicted_reduction) =
            compute_trust_region_step_constrained(&lag_grad, &b, &a, &c, delta, constraints, ctol);

        // Try the step
        let x_new = &x + &p;
        let f_new = func(x_new.as_slice().expect("Operation failed"));
        nfev += 1;

        let mut c_new = Array1::zeros(constraints.len());
        for (i, constraint) in constraints.iter().enumerate() {
            if !constraint.is_bounds() {
                c_new[i] = (constraint.fun)(x_new.as_slice().expect("Operation failed"));
                nfev += 1;
            }
        }

        // Compute merit function
        let mut merit = f;
        let mut merit_new = f_new;
        let penalty = 10.0;

        for (i, &ci) in c.iter().enumerate() {
            match constraints[i].kind {
                ConstraintKind::Inequality => {
                    merit += penalty * f64::max(0.0, -ci);
                    merit_new += penalty * f64::max(0.0, -c_new[i]);
                }
                ConstraintKind::Equality => {
                    merit += penalty * ci.abs();
                    merit_new += penalty * c_new[i].abs();
                }
            }
        }

        let actual_reduction = merit - merit_new;
        let rho = if predicted_reduction > 0.0 {
            actual_reduction / predicted_reduction
        } else {
            0.0
        };

        // Update trust region radius
        if rho < 0.25 {
            delta *= 0.5;
        } else if rho > 0.75 && p.iter().map(|&pi| pi * pi).sum::<f64>().sqrt() >= 0.9 * delta {
            delta *= 2.0;
        }

        // Accept or reject the step
        if rho > 0.1 {
            x = x_new;
            f = f_new;
            c = c_new;

            if (merit - merit_new).abs() < ftol * (1.0 + merit.abs()) {
                break;
            }

            // Compute new gradient
            let g_new = if let Some(ref grad_fn) = jac {
                njev += 1;
                grad_fn(x.as_slice().expect("Operation failed"))
            } else {
                let mut grad = Array1::zeros(n);
                for i in 0..n {
                    let mut x_h = x.clone();
                    x_h[i] += eps;
                    let f_h = func(x_h.as_slice().expect("Operation failed"));
                    grad[i] = (f_h - f) / eps;
                    nfev += 1;
                }
                grad
            };

            // Compute new constraint Jacobian
            let mut a_new = Array2::zeros((constraints.len(), n));
            for (i, constraint) in constraints.iter().enumerate() {
                if !constraint.is_bounds() {
                    for j in 0..n {
                        let mut x_h = x.clone();
                        x_h[j] += eps;
                        let c_h = (constraint.fun)(x_h.as_slice().expect("Operation failed"));
                        a_new[[i, j]] = (c_h - c[i]) / eps;
                        nfev += 1;
                    }
                }
            }

            // Update Lagrange multipliers
            for (i, constraint) in constraints.iter().enumerate() {
                match constraint.kind {
                    ConstraintKind::Inequality => {
                        if c[i] < -ctol {
                            lambda[i] = f64::max(0.0, lambda[i] - c[i] * penalty);
                        } else {
                            lambda[i] = f64::max(0.0, lambda[i] - 0.1 * lambda[i]);
                        }
                    }
                    ConstraintKind::Equality => {
                        let step_size = 0.1;
                        lambda[i] -= step_size * c[i] * penalty;
                        lambda[i] *= 0.9;
                    }
                }
            }

            // Update Hessian based on strategy
            match hess_update {
                HessianUpdate::Exact => {
                    if let Some(ref hess_fn) = hess {
                        nhev += 1;
                        b = hess_fn(x.as_slice().expect("Operation failed"));
                    }
                }
                HessianUpdate::BFGS => {
                    let s = &p;
                    let y = &g_new - &g;
                    let s_dot_y = s.dot(&y);
                    if s_dot_y > 1e-10 {
                        let s_col = s.clone().insert_axis(Axis(1));
                        let s_row = s.clone().insert_axis(Axis(0));
                        let bs = b.dot(s);
                        let bs_col = bs.clone().insert_axis(Axis(1));
                        let bs_row = bs.clone().insert_axis(Axis(0));
                        let term1 = s_dot_y + s.dot(&bs);
                        let term2 = &s_col.dot(&s_row) * (term1 / (s_dot_y * s_dot_y));
                        let term3 = &bs_col.dot(&s_row) + &s_col.dot(&bs_row);
                        b = &b + &term2 - &(&term3 / s_dot_y);
                    }
                }
                HessianUpdate::SR1 => {
                    let s = &p;
                    let y = &g_new - &g;
                    let bs = b.dot(s);
                    let diff = &y - &bs;
                    let s_dot_diff = s.dot(&diff);
                    // SR1 update with skipping condition
                    if s_dot_diff.abs() > 1e-8 * s.dot(s).sqrt() * diff.dot(&diff).sqrt() {
                        let diff_col = diff.clone().insert_axis(Axis(1));
                        let diff_row = diff.clone().insert_axis(Axis(0));
                        let update = &diff_col.dot(&diff_row) / s_dot_diff;
                        b = &b + &update;
                    }
                }
            }

            g = g_new;
            a = a_new;
        }

        iter += 1;
    }

    // Prepare result
    let mut c_result = Array1::zeros(constraints.len());
    for (i, constraint) in constraints.iter().enumerate() {
        if !constraint.is_bounds() {
            c_result[i] = c[i];
        }
    }

    let mut result = OptimizeResults::default();
    result.x = x;
    result.fun = f;
    result.jac = Some(g.into_raw_vec_and_offset().0);
    result.constr = Some(c_result);
    result.nfev = nfev;
    result.njev = njev;
    result.nhev = nhev;
    result.nit = iter;
    result.success = iter < maxiter;

    if result.success {
        result.message = "Optimization terminated successfully.".to_string();
    } else {
        result.message = "Maximum iterations reached.".to_string();
    }

    Ok(result)
}

/// Compute a trust-region step for constrained optimization
#[allow(clippy::many_single_char_names)]
#[allow(dead_code)]
fn compute_trust_region_step_constrained(
    g: &Array1<f64>,
    b: &Array2<f64>,
    a: &Array2<f64>,
    c: &Array1<f64>,
    delta: f64,
    constraints: &[Constraint<ConstraintFn>],
    ctol: f64,
) -> (Array1<f64>, f64) {
    let n = g.len();
    let n_constr = constraints.len();

    // Compute the unconstrained Cauchy point (steepest descent direction)
    let p_unconstrained = compute_unconstrained_cauchy_point(g, b, delta);

    // Check if unconstrained Cauchy point satisfies linearized constraints
    let mut constraint_violated = false;
    for i in 0..n_constr {
        let grad_c_dot_p = (0..n).map(|j| a[[i, j]] * p_unconstrained[j]).sum::<f64>();

        match constraints[i].kind {
            ConstraintKind::Inequality => {
                // For inequality constraints: check if g(x) + grad_g^T p >= -ctol
                if c[i] + grad_c_dot_p < -ctol {
                    constraint_violated = true;
                    break;
                }
            }
            ConstraintKind::Equality => {
                // For equality constraints: check if |h(x) + grad_h^T p| <= ctol
                if (c[i] + grad_c_dot_p).abs() > ctol {
                    constraint_violated = true;
                    break;
                }
            }
        }
    }

    // If unconstrained point is feasible, use it
    if !constraint_violated {
        // Compute predicted reduction
        let g_dot_p = g.dot(&p_unconstrained);
        let bp = b.dot(&p_unconstrained);
        let p_dot_bp = p_unconstrained.dot(&bp);
        let predicted_reduction = -g_dot_p - 0.5 * p_dot_bp;

        return (p_unconstrained, predicted_reduction);
    }

    // Otherwise, project onto the linearized feasible region
    // This is a simplified approach - in practice, you would solve a CQP

    // Start with the steepest descent direction
    let mut p = Array1::zeros(n);
    for i in 0..n {
        p[i] = -g[i];
    }

    // Normalize to trust region radius
    let p_norm = p.iter().map(|&pi| pi * pi).sum::<f64>().sqrt();
    if p_norm > 1e-10 {
        p = &p * (delta / p_norm);
    }

    // Project onto each constraint
    for _iter in 0..5 {
        // Limited iterations for projection
        let mut max_viol = 0.0;
        let mut most_violated = 0;

        // Find most violated constraint
        for i in 0..n_constr {
            let grad_c_dot_p = (0..n).map(|j| a[[i, j]] * p[j]).sum::<f64>();

            let viol = match constraints[i].kind {
                ConstraintKind::Inequality => {
                    // For inequality constraints: violation is max(0, -(g + grad_g^T p))
                    f64::max(0.0, -(c[i] + grad_c_dot_p))
                }
                ConstraintKind::Equality => {
                    // For equality constraints: violation is |h + grad_h^T p|
                    (c[i] + grad_c_dot_p).abs()
                }
            };

            if viol > max_viol {
                max_viol = viol;
                most_violated = i;
            }
        }

        if max_viol < ctol {
            break;
        }

        // Project p onto the constraint
        let mut a_norm_sq = 0.0;
        for j in 0..n {
            a_norm_sq += a[[most_violated, j]] * a[[most_violated, j]];
        }

        if a_norm_sq > 1e-10 {
            let grad_c_dot_p = (0..n).map(|j| a[[most_violated, j]] * p[j]).sum::<f64>();

            let proj_dist = match constraints[most_violated].kind {
                ConstraintKind::Inequality => {
                    // For inequality constraints: project to boundary when violated
                    if c[most_violated] + grad_c_dot_p < 0.0 {
                        -(c[most_violated] + grad_c_dot_p) / a_norm_sq
                    } else {
                        0.0
                    }
                }
                ConstraintKind::Equality => {
                    // For equality constraints: always project to satisfy h + grad_h^T p = 0
                    -(c[most_violated] + grad_c_dot_p) / a_norm_sq
                }
            };

            // Project p
            for j in 0..n {
                p[j] += a[[most_violated, j]] * proj_dist;
            }

            // Rescale to trust region
            let p_norm = p.iter().map(|&pi| pi * pi).sum::<f64>().sqrt();
            if p_norm > delta {
                p = &p * (delta / p_norm);
            }
        }
    }

    // Compute predicted reduction
    let g_dot_p = g.dot(&p);
    let bp = b.dot(&p);
    let p_dot_bp = p.dot(&bp);
    let predicted_reduction = -g_dot_p - 0.5 * p_dot_bp;

    (p, predicted_reduction)
}

/// Compute the unconstrained Cauchy point (steepest descent to trust region boundary)
#[allow(dead_code)]
fn compute_unconstrained_cauchy_point(g: &Array1<f64>, b: &Array2<f64>, delta: f64) -> Array1<f64> {
    let n = g.len();

    // Compute the steepest descent direction: -g
    let mut p = Array1::zeros(n);
    for i in 0..n {
        p[i] = -g[i];
    }

    // Compute g^T B g and ||g||^2
    let bg = b.dot(g);
    let g_dot_bg = g.dot(&bg);
    let g_norm_sq = g.dot(g);

    // Check if the gradient is practically zero
    if g_norm_sq < 1e-10 {
        // If gradient is practically zero, don't move
        return Array1::zeros(n);
    }

    // Compute tau (step to the boundary)
    let tau = if g_dot_bg <= 0.0 {
        // Negative curvature or zero curvature case
        delta / g_norm_sq.sqrt()
    } else {
        // Positive curvature case
        f64::min(delta / g_norm_sq.sqrt(), g_norm_sq / g_dot_bg)
    };

    // Scale the direction by tau
    for i in 0..n {
        p[i] *= tau;
    }

    p
}