1use std::sync::Arc;
22
23use scirs2_symbolic::eml::eval::{eval_real, EvalCtx};
24use scirs2_symbolic::eml::{grad, hessian, LoweredOp};
25
26#[derive(Debug, Clone)]
28pub struct SymbolicNewtonResult {
29 pub x: Vec<f64>,
31 pub f_final: f64,
33 pub iters: usize,
35 pub converged: bool,
37}
38
39#[derive(Debug, thiserror::Error)]
41pub enum SymbolicNewtonError {
42 #[error("evaluation error: {0}")]
45 EvalError(String),
46 #[error("singular Hessian — cannot invert")]
49 SingularHessian,
50 #[error("dimension mismatch: x0 has {x0_len}, cost expects {n_vars}")]
52 DimMismatch {
53 x0_len: usize,
55 n_vars: usize,
57 },
58}
59
60pub fn newton(
75 cost: &LoweredOp,
76 x0: &[f64],
77 max_iter: usize,
78 tol: f64,
79) -> Result<SymbolicNewtonResult, SymbolicNewtonError> {
80 let n_vars = cost.count_vars();
81 if x0.len() < n_vars {
82 return Err(SymbolicNewtonError::DimMismatch {
83 x0_len: x0.len(),
84 n_vars,
85 });
86 }
87
88 let grad_ops: Vec<LoweredOp> = (0..n_vars).map(|i| grad(cost, i)).collect();
90 let hess_ops: Vec<Vec<LoweredOp>> = hessian(cost, n_vars);
91
92 let mut x = x0.to_vec();
93 let mut converged = false;
94 let mut iters = 0;
95
96 for k in 0..max_iter {
97 iters = k + 1;
98 let ctx = EvalCtx::new(&x);
99
100 let mut grad_vec: Vec<f64> = Vec::with_capacity(n_vars);
102 for g_op in &grad_ops {
103 let v =
104 eval_real(g_op, &ctx).map_err(|e| SymbolicNewtonError::EvalError(e.to_string()))?;
105 grad_vec.push(v);
106 }
107
108 let grad_norm: f64 = grad_vec.iter().map(|g| g * g).sum::<f64>().sqrt();
110 if grad_norm < tol {
111 converged = true;
112 break;
113 }
114
115 let mut hess_mat: Vec<Vec<f64>> = Vec::with_capacity(n_vars);
117 for row in &hess_ops {
118 let mut hess_row: Vec<f64> = Vec::with_capacity(n_vars);
119 for h_op in row {
120 let v = eval_real(h_op, &ctx)
121 .map_err(|e| SymbolicNewtonError::EvalError(e.to_string()))?;
122 hess_row.push(v);
123 }
124 hess_mat.push(hess_row);
125 }
126
127 let dx = solve_linear(&hess_mat, &grad_vec)?;
129
130 for (xi, dxi) in x.iter_mut().zip(dx.iter()) {
132 *xi -= *dxi;
133 }
134 }
135
136 let final_ctx = EvalCtx::new(&x);
137 let f_final =
138 eval_real(cost, &final_ctx).map_err(|e| SymbolicNewtonError::EvalError(e.to_string()))?;
139
140 Ok(SymbolicNewtonResult {
141 x,
142 f_final,
143 iters,
144 converged,
145 })
146}
147
148#[derive(Debug, Clone)]
154pub struct SymbolicOptResult {
155 pub x: scirs2_core::ndarray::Array1<f64>,
157 pub f_val: f64,
159 pub grad_norm: f64,
161 pub iters: usize,
163 pub converged: bool,
165}
166
167#[derive(Debug, thiserror::Error)]
169pub enum SymbolicOptError {
170 #[error("eval error: {0}")]
172 EvalError(String),
173 #[error("gradient eval error: {0}")]
175 GradEvalError(String),
176 #[error("Hessian eval error: {0}")]
178 HessEvalError(String),
179 #[error("dimension mismatch: expected {expected} variables, got {got}")]
181 DimMismatch {
182 expected: usize,
184 got: usize,
186 },
187 #[error("not converged after {iters} iters (grad_norm = {grad_norm:.6e})")]
189 NotConverged {
190 iters: usize,
192 grad_norm: f64,
194 },
195 #[error("line search failed")]
197 LineSearchFailed,
198}
199
200fn eval_gradient(grad_ops: &[LoweredOp], x: &[f64]) -> Result<Vec<f64>, SymbolicOptError> {
206 let ctx = EvalCtx::new(x);
207 let mut g = Vec::with_capacity(grad_ops.len());
208 for op in grad_ops {
209 let v = eval_real(op, &ctx).map_err(|e| SymbolicOptError::GradEvalError(e.to_string()))?;
210 g.push(v);
211 }
212 Ok(g)
213}
214
215fn eval_objective(obj: &LoweredOp, x: &[f64]) -> Result<f64, SymbolicOptError> {
217 let ctx = EvalCtx::new(x);
218 eval_real(obj, &ctx).map_err(|e| SymbolicOptError::EvalError(e.to_string()))
219}
220
221fn l2_norm(v: &[f64]) -> f64 {
223 v.iter().map(|x| x * x).sum::<f64>().sqrt()
224}
225
226fn dot(a: &[f64], b: &[f64]) -> f64 {
228 a.iter().zip(b.iter()).map(|(ai, bi)| ai * bi).sum()
229}
230
231pub fn lbfgs_symbolic(
272 objective: &Arc<LoweredOp>,
273 x0: scirs2_core::ndarray::ArrayView1<f64>,
274 max_iter: usize,
275 tol: f64,
276 memory: usize,
277) -> Result<SymbolicOptResult, SymbolicOptError> {
278 use scirs2_core::ndarray::Array1;
279
280 let n_vars = objective.count_vars();
281 if x0.len() != n_vars {
282 return Err(SymbolicOptError::DimMismatch {
283 expected: n_vars,
284 got: x0.len(),
285 });
286 }
287
288 let grad_ops: Vec<LoweredOp> = (0..n_vars).map(|i| grad(objective.as_ref(), i)).collect();
290
291 let mut x: Vec<f64> = x0.iter().copied().collect();
292
293 let mut g = eval_gradient(&grad_ops, &x)?;
295 let mut grad_norm = l2_norm(&g);
296
297 if max_iter == 0 {
299 return Err(SymbolicOptError::NotConverged {
300 iters: 0,
301 grad_norm,
302 });
303 }
304
305 if grad_norm < tol {
306 let f_val = eval_objective(objective.as_ref(), &x)?;
307 return Ok(SymbolicOptResult {
308 x: Array1::from_vec(x),
309 f_val,
310 grad_norm,
311 iters: 0,
312 converged: true,
313 });
314 }
315
316 let mem = memory.max(1);
319 let mut history: std::collections::VecDeque<(Vec<f64>, Vec<f64>, f64)> =
320 std::collections::VecDeque::with_capacity(mem);
321
322 let mut converged = false;
323 let mut iters = 0;
324
325 for k in 0..max_iter {
326 iters = k + 1;
327
328 let direction = lbfgs_direction(&g, &history);
330
331 let f_k = eval_objective(objective.as_ref(), &x)?;
333 let dg = dot(&g, &direction); const C1: f64 = 1e-4; const C2: f64 = 0.9; const MAX_LS: usize = 20;
338
339 let mut alpha = 1.0_f64;
340 let mut x_new: Vec<f64> = Vec::with_capacity(n_vars);
341 let mut g_new: Vec<f64>;
342 let mut armijo_ok = false;
343
344 let mut ls_iter = 0;
346 loop {
347 x_new = x
348 .iter()
349 .zip(&direction)
350 .map(|(xi, di)| xi + alpha * di)
351 .collect();
352 let f_new = eval_objective(objective.as_ref(), &x_new)?;
353 let suff_dec = f_new <= f_k + C1 * alpha * dg;
354
355 if suff_dec {
356 g_new = eval_gradient(&grad_ops, &x_new)?;
357 let curvature = dot(&g_new, &direction).abs() <= C2 * dg.abs();
358 if curvature {
359 armijo_ok = true;
360 break;
361 }
362 armijo_ok = true;
364 } else {
366 g_new = Vec::new(); }
368
369 ls_iter += 1;
370 if ls_iter >= MAX_LS {
371 break;
372 }
373 alpha *= 0.5;
374 }
375
376 if !armijo_ok {
378 x_new = x
380 .iter()
381 .zip(&direction)
382 .map(|(xi, di)| xi + alpha * di)
383 .collect();
384 let f_new = eval_objective(objective.as_ref(), &x_new)?;
385 if f_new > f_k + C1 * alpha * dg {
386 return Err(SymbolicOptError::LineSearchFailed);
387 }
388 g_new = eval_gradient(&grad_ops, &x_new)?;
389 } else if g_new.is_empty() {
390 g_new = eval_gradient(&grad_ops, &x_new)?;
392 }
393
394 let s_k: Vec<f64> = x_new.iter().zip(&x).map(|(xn, xo)| xn - xo).collect();
396 let y_k: Vec<f64> = g_new.iter().zip(&g).map(|(gn, go)| gn - go).collect();
397 let sy = dot(&s_k, &y_k);
398 if sy.abs() > 1e-20 {
399 let rho = 1.0 / sy;
400 if history.len() == mem {
401 history.pop_front();
402 }
403 history.push_back((s_k, y_k, rho));
404 }
405
406 x = x_new;
408 g = g_new;
409 grad_norm = l2_norm(&g);
410
411 if grad_norm < tol {
412 converged = true;
413 break;
414 }
415 }
416
417 let f_val = eval_objective(objective.as_ref(), &x)?;
418
419 if converged {
420 Ok(SymbolicOptResult {
421 x: Array1::from_vec(x),
422 f_val,
423 grad_norm,
424 iters,
425 converged: true,
426 })
427 } else {
428 Err(SymbolicOptError::NotConverged { iters, grad_norm })
429 }
430}
431
432fn lbfgs_direction(
436 g: &[f64],
437 history: &std::collections::VecDeque<(Vec<f64>, Vec<f64>, f64)>,
438) -> Vec<f64> {
439 let n = g.len();
440 let m = history.len();
441 let mut q: Vec<f64> = g.to_vec();
442 let mut alphas: Vec<f64> = vec![0.0; m];
443
444 for (idx, (s, y, rho)) in history.iter().rev().enumerate() {
446 let a = rho * dot(s, &q);
447 alphas[m - 1 - idx] = a;
448 for j in 0..n {
449 q[j] -= a * y[j];
450 }
451 }
452
453 let mut r: Vec<f64> = if let Some((s, y, _)) = history.back() {
455 let sy = dot(s, y);
456 let yy = dot(y, y);
457 let gamma = if yy > 1e-20 { sy / yy } else { 1.0 };
458 q.iter().map(|qi| gamma * qi).collect()
459 } else {
460 q.clone() };
462
463 for (idx, (s, y, rho)) in history.iter().enumerate() {
465 let beta = rho * dot(y, &r);
466 let a = alphas[idx];
467 for j in 0..n {
468 r[j] += s[j] * (a - beta);
469 }
470 }
471
472 r.iter().map(|ri| -ri).collect()
474}
475
476pub fn trust_region_symbolic(
518 objective: &Arc<LoweredOp>,
519 x0: scirs2_core::ndarray::ArrayView1<f64>,
520 max_iter: usize,
521 tol: f64,
522 initial_radius: f64,
523) -> Result<SymbolicOptResult, SymbolicOptError> {
524 use scirs2_core::ndarray::Array1;
525
526 let n_vars = objective.count_vars();
527 if x0.len() != n_vars {
528 return Err(SymbolicOptError::DimMismatch {
529 expected: n_vars,
530 got: x0.len(),
531 });
532 }
533
534 let grad_ops: Vec<LoweredOp> = (0..n_vars).map(|i| grad(objective.as_ref(), i)).collect();
536 let hess_ops: Vec<Vec<LoweredOp>> = hessian(objective.as_ref(), n_vars);
537
538 let mut x: Vec<f64> = x0.iter().copied().collect();
539
540 let mut g = eval_gradient(&grad_ops, &x)?;
542 let mut grad_norm = l2_norm(&g);
543
544 if max_iter == 0 {
545 return Err(SymbolicOptError::NotConverged {
546 iters: 0,
547 grad_norm,
548 });
549 }
550
551 if grad_norm < tol {
552 let f_val = eval_objective(objective.as_ref(), &x)?;
553 return Ok(SymbolicOptResult {
554 x: Array1::from_vec(x),
555 f_val,
556 grad_norm,
557 iters: 0,
558 converged: true,
559 });
560 }
561
562 let mut delta = initial_radius.max(1e-8);
563 let mut converged = false;
564 let mut iters = 0;
565
566 for k in 0..max_iter {
567 iters = k + 1;
568
569 let f_k = eval_objective(objective.as_ref(), &x)?;
570
571 let ctx = EvalCtx::new(&x);
573 let mut h_mat: Vec<Vec<f64>> = Vec::with_capacity(n_vars);
574 for row in &hess_ops {
575 let mut h_row: Vec<f64> = Vec::with_capacity(n_vars);
576 for h_op in row {
577 let v = eval_real(h_op, &ctx)
578 .map_err(|e| SymbolicOptError::HessEvalError(e.to_string()))?;
579 h_row.push(v);
580 }
581 h_mat.push(h_row);
582 }
583
584 let p = dogleg_step_sym(&g, &h_mat, delta);
586
587 let x_new: Vec<f64> = x.iter().zip(&p).map(|(xi, pi)| xi + pi).collect();
589 let f_new = eval_objective(objective.as_ref(), &x_new)?;
590 let actual_red = f_k - f_new;
591
592 let gtp = dot(&g, &p);
594 let htp: Vec<f64> = matvec(&h_mat, &p);
595 let phtp = dot(&p, &htp);
596 let pred_red = -(gtp + 0.5 * phtp);
597
598 let rho = if pred_red.abs() < 1e-20 {
600 1.0
601 } else {
602 actual_red / pred_red
603 };
604
605 if rho >= 0.1 {
606 x = x_new;
608 g = eval_gradient(&grad_ops, &x)?;
609 grad_norm = l2_norm(&g);
610
611 if grad_norm < tol {
612 converged = true;
613 if rho >= 0.75 {
615 let p_norm = l2_norm(&p);
616 if p_norm >= 0.8 * delta {
617 delta = (2.0 * delta).min(100.0);
618 }
619 }
620 break;
621 }
622 }
623 if rho >= 0.75 {
627 let p_norm = l2_norm(&p);
628 if p_norm >= 0.8 * delta {
629 delta = (2.0 * delta).min(100.0);
630 }
631 } else if rho < 0.25 {
632 delta *= 0.25;
633 }
634
635 if delta < 1e-14 {
636 break;
638 }
639 }
640
641 let f_val = eval_objective(objective.as_ref(), &x)?;
642
643 if converged {
644 Ok(SymbolicOptResult {
645 x: Array1::from_vec(x),
646 f_val,
647 grad_norm,
648 iters,
649 converged: true,
650 })
651 } else {
652 Err(SymbolicOptError::NotConverged { iters, grad_norm })
653 }
654}
655
656fn dogleg_step_sym(g: &[f64], h: &[Vec<f64>], delta: f64) -> Vec<f64> {
662 let n = g.len();
663
664 let neg_g: Vec<f64> = g.iter().map(|gi| -gi).collect();
666 let p_n_opt = solve_linear_opt(h, &neg_g);
667
668 let gtg = dot(g, g);
670 let hg = matvec(h, g);
671 let gthg = dot(g, &hg); let p_sd: Vec<f64> = if gthg > 1e-20 {
674 let scale = gtg / gthg;
675 g.iter().map(|gi| -scale * gi).collect()
676 } else {
677 let gn = gtg.sqrt().max(1e-20);
679 g.iter().map(|gi| -gi / gn).collect()
680 };
681
682 let p_sd_norm = l2_norm(&p_sd);
683
684 if let Some(p_n) = p_n_opt {
686 let p_n_norm = l2_norm(&p_n);
687
688 if p_n_norm <= delta {
689 return p_n;
691 }
692
693 if p_sd_norm >= delta {
694 return p_sd.iter().map(|pi| delta / p_sd_norm * pi).collect();
696 }
697
698 let d: Vec<f64> = p_n.iter().zip(&p_sd).map(|(pni, psi)| pni - psi).collect();
702 let dd = dot(&d, &d);
703 let psd_d = dot(&p_sd, &d);
704 let psd2 = p_sd_norm * p_sd_norm;
705 let discriminant = psd_d * psd_d - dd * (psd2 - delta * delta);
706 let tau = if dd < 1e-20 || discriminant < 0.0 {
707 1.0
708 } else {
709 (-psd_d + discriminant.sqrt()) / dd
710 };
711 let tau = tau.clamp(0.0, 1.0);
712 p_sd.iter()
713 .zip(&d)
714 .map(|(psi, di)| psi + tau * di)
715 .collect()
716 } else {
717 if p_sd_norm <= delta {
719 p_sd
720 } else {
721 p_sd.iter().map(|pi| delta / p_sd_norm * pi).collect()
722 }
723 }
724}
725
726fn matvec(h: &[Vec<f64>], v: &[f64]) -> Vec<f64> {
728 h.iter().map(|row| dot(row, v)).collect()
729}
730
731fn solve_linear_opt(a: &[Vec<f64>], b: &[f64]) -> Option<Vec<f64>> {
733 let n = b.len();
734 if n == 0 {
735 return Some(Vec::new());
736 }
737
738 let mut mat: Vec<Vec<f64>> = a
739 .iter()
740 .zip(b.iter())
741 .map(|(row, &bi)| {
742 let mut r = row.clone();
743 r.push(bi);
744 r
745 })
746 .collect();
747
748 for k in 0..n {
749 let mut max_idx = k;
751 let mut max_val = mat[k][k].abs();
752 for i in (k + 1)..n {
753 let v = mat[i][k].abs();
754 if v > max_val {
755 max_val = v;
756 max_idx = i;
757 }
758 }
759 if max_val < 1e-12 {
760 return None; }
762 mat.swap(k, max_idx);
763 for i in (k + 1)..n {
764 let factor = mat[i][k] / mat[k][k];
765 for j in k..(n + 1) {
766 mat[i][j] -= factor * mat[k][j];
767 }
768 }
769 }
770
771 let mut x = vec![0.0; n];
772 for i in (0..n).rev() {
773 let mut sum = mat[i][n];
774 for j in (i + 1)..n {
775 sum -= mat[i][j] * x[j];
776 }
777 x[i] = sum / mat[i][i];
778 }
779 Some(x)
780}
781
782fn solve_linear(a: &[Vec<f64>], b: &[f64]) -> Result<Vec<f64>, SymbolicNewtonError> {
784 let n = b.len();
785 if n == 0 {
786 return Ok(Vec::new());
787 }
788
789 let mut mat: Vec<Vec<f64>> = a
791 .iter()
792 .zip(b.iter())
793 .map(|(row, &bi)| {
794 let mut r = row.clone();
795 r.push(bi);
796 r
797 })
798 .collect();
799
800 for k in 0..n {
802 let mut max_idx = k;
804 let mut max_val = mat[k][k].abs();
805 for i in (k + 1)..n {
806 let v = mat[i][k].abs();
807 if v > max_val {
808 max_val = v;
809 max_idx = i;
810 }
811 }
812 if max_val < 1e-12 {
813 return Err(SymbolicNewtonError::SingularHessian);
814 }
815 mat.swap(k, max_idx);
817 for i in (k + 1)..n {
819 let factor = mat[i][k] / mat[k][k];
820 for j in k..(n + 1) {
821 mat[i][j] -= factor * mat[k][j];
822 }
823 }
824 }
825
826 let mut x = vec![0.0; n];
828 for i in (0..n).rev() {
829 let mut sum = mat[i][n];
830 for j in (i + 1)..n {
831 sum -= mat[i][j] * x[j];
832 }
833 x[i] = sum / mat[i][i];
834 }
835
836 Ok(x)
837}
838
839pub struct KktSystem {
848 pub lagrangian: LoweredOp,
850 pub stationarity: Vec<LoweredOp>,
852 pub constraint_residuals: Vec<LoweredOp>,
854 pub n_vars: usize,
856 pub n_constraints: usize,
858}
859
860#[derive(Debug)]
862pub enum LagrangianError {
863 DimMismatch {
865 expected: usize,
867 got: usize,
869 },
870 GradError(String),
872}
873
874pub fn build_kkt(
888 objective: &Arc<LoweredOp>,
889 constraints: &[Arc<LoweredOp>],
890 n_vars: usize,
891) -> Result<KktSystem, LagrangianError> {
892 let m = constraints.len();
893
894 let obj_vars = objective.count_vars();
897 if obj_vars > n_vars {
898 return Err(LagrangianError::DimMismatch {
899 expected: n_vars,
900 got: obj_vars,
901 });
902 }
903
904 let mut lagrangian: LoweredOp = objective.as_ref().clone();
907 for (i, g) in constraints.iter().enumerate() {
908 let lam_i = LoweredOp::Var(n_vars + i);
909 let term = LoweredOp::Mul(Box::new(lam_i), Box::new(g.as_ref().clone()));
910 lagrangian = LoweredOp::Add(Box::new(lagrangian), Box::new(term));
911 }
912
913 let stationarity: Vec<LoweredOp> = (0..n_vars).map(|j| grad(&lagrangian, j)).collect();
915
916 let constraint_residuals: Vec<LoweredOp> =
918 constraints.iter().map(|g| g.as_ref().clone()).collect();
919
920 Ok(KktSystem {
921 lagrangian,
922 stationarity,
923 constraint_residuals,
924 n_vars,
925 n_constraints: m,
926 })
927}
928
929pub fn solve_lagrangian_symbolic(
948 objective: &Arc<LoweredOp>,
949 constraints: &[Arc<LoweredOp>],
950 x0: scirs2_core::ndarray::ArrayView1<f64>,
951 lambda0: scirs2_core::ndarray::ArrayView1<f64>,
952 max_iter: usize,
953 tol: f64,
954) -> Result<SymbolicOptResult, SymbolicOptError> {
955 use scirs2_core::ndarray::Array1;
956
957 let m = constraints.len();
958
959 let n_vars = std::iter::once(objective.count_vars())
962 .chain(constraints.iter().map(|c| c.count_vars()))
963 .max()
964 .unwrap_or(0);
965
966 if x0.len() != n_vars {
967 return Err(SymbolicOptError::DimMismatch {
968 expected: n_vars,
969 got: x0.len(),
970 });
971 }
972 if lambda0.len() != m {
973 return Err(SymbolicOptError::DimMismatch {
974 expected: m,
975 got: lambda0.len(),
976 });
977 }
978
979 let kkt = build_kkt(objective, constraints, n_vars).map_err(|e| match e {
981 LagrangianError::DimMismatch { expected, got } => {
982 SymbolicOptError::DimMismatch { expected, got }
983 }
984 LagrangianError::GradError(s) => SymbolicOptError::GradEvalError(s),
985 })?;
986
987 let big_n = n_vars + m;
989
990 let mut z: Vec<f64> = x0.iter().copied().chain(lambda0.iter().copied()).collect();
992
993 let eqs: Vec<&LoweredOp> = kkt
995 .stationarity
996 .iter()
997 .chain(kkt.constraint_residuals.iter())
998 .collect();
999
1000 let jac_ops: Vec<Vec<LoweredOp>> = eqs
1002 .iter()
1003 .map(|eq| (0..big_n).map(|j| grad(eq, j)).collect())
1004 .collect();
1005
1006 let initial_residual = {
1008 let ctx = EvalCtx::new(&z);
1009 let mut f_vec = Vec::with_capacity(big_n);
1010 for eq in &eqs {
1011 let v = eval_real(eq, &ctx).map_err(|e| SymbolicOptError::EvalError(e.to_string()))?;
1012 f_vec.push(v);
1013 }
1014 f_vec
1015 };
1016 let mut residual_norm = l2_norm(&initial_residual);
1017
1018 if max_iter == 0 {
1019 return Err(SymbolicOptError::NotConverged {
1020 iters: 0,
1021 grad_norm: residual_norm,
1022 });
1023 }
1024
1025 if residual_norm < tol {
1027 let f_val = eval_objective(objective.as_ref(), &z[..n_vars])?;
1028 return Ok(SymbolicOptResult {
1029 x: Array1::from_vec(z[..n_vars].to_vec()),
1030 f_val,
1031 grad_norm: residual_norm,
1032 iters: 0,
1033 converged: true,
1034 });
1035 }
1036
1037 let mut converged = false;
1038 let mut iters = 0;
1039
1040 for k in 0..max_iter {
1041 iters = k + 1;
1042
1043 let ctx = EvalCtx::new(&z);
1045 let mut f_vec: Vec<f64> = Vec::with_capacity(big_n);
1046 for eq in &eqs {
1047 let v = eval_real(eq, &ctx).map_err(|e| SymbolicOptError::EvalError(e.to_string()))?;
1048 f_vec.push(v);
1049 }
1050
1051 residual_norm = l2_norm(&f_vec);
1052 if residual_norm < tol {
1053 converged = true;
1054 break;
1055 }
1056
1057 let mut jac_mat: Vec<Vec<f64>> = Vec::with_capacity(big_n);
1059 for row_ops in &jac_ops {
1060 let mut row: Vec<f64> = Vec::with_capacity(big_n);
1061 for op in row_ops {
1062 let v =
1063 eval_real(op, &ctx).map_err(|e| SymbolicOptError::EvalError(e.to_string()))?;
1064 row.push(v);
1065 }
1066 jac_mat.push(row);
1067 }
1068
1069 let neg_f: Vec<f64> = f_vec.iter().map(|v| -v).collect();
1071 let dz = match solve_linear_opt(&jac_mat, &neg_f) {
1072 Some(d) => d,
1073 None => {
1074 neg_f
1076 }
1077 };
1078
1079 for (zi, dzi) in z.iter_mut().zip(dz.iter()) {
1081 *zi += dzi;
1082 }
1083 }
1084
1085 let f_val = eval_objective(objective.as_ref(), &z[..n_vars])?;
1086
1087 if converged {
1088 Ok(SymbolicOptResult {
1089 x: Array1::from_vec(z[..n_vars].to_vec()),
1090 f_val,
1091 grad_norm: residual_norm,
1092 iters,
1093 converged: true,
1094 })
1095 } else {
1096 Err(SymbolicOptError::NotConverged {
1097 iters,
1098 grad_norm: residual_norm,
1099 })
1100 }
1101}
1102
1103pub mod line_search;
1108pub use line_search::{OptLineSearchError, SymbolicLineSearch};
1109
1110#[cfg(test)]
1111mod tests {
1112 use super::*;
1113
1114 fn var(i: usize) -> LoweredOp {
1115 LoweredOp::Var(i)
1116 }
1117 fn c(v: f64) -> LoweredOp {
1118 LoweredOp::Const(v)
1119 }
1120
1121 #[test]
1122 fn newton_x_squared_one_step() {
1123 let cost = LoweredOp::Mul(Box::new(var(0)), Box::new(var(0)));
1125 let result = newton(&cost, &[5.0], 10, 1e-8).expect("converge");
1126 assert!(result.converged);
1127 assert!(result.x[0].abs() < 1e-6, "x = {}", result.x[0]);
1128 assert!(result.iters <= 2);
1129 }
1130
1131 #[test]
1132 fn newton_quadratic_2d() {
1133 let cost = LoweredOp::Add(
1135 Box::new(LoweredOp::Mul(Box::new(var(0)), Box::new(var(0)))),
1136 Box::new(LoweredOp::Mul(Box::new(var(1)), Box::new(var(1)))),
1137 );
1138 let result = newton(&cost, &[3.0, -4.0], 10, 1e-8).expect("converge");
1139 assert!(result.converged);
1140 assert!(result.x[0].abs() < 1e-6);
1141 assert!(result.x[1].abs() < 1e-6);
1142 }
1143
1144 #[test]
1145 fn newton_shifted_quadratic() {
1146 let inner = LoweredOp::Sub(Box::new(var(0)), Box::new(c(3.0)));
1148 let cost = LoweredOp::Mul(Box::new(inner.clone()), Box::new(inner));
1149 let result = newton(&cost, &[0.0], 10, 1e-8).expect("converge");
1150 assert!(result.converged);
1151 assert!((result.x[0] - 3.0).abs() < 1e-6, "x = {}", result.x[0]);
1152 }
1153
1154 #[test]
1155 fn newton_dim_mismatch_returns_err() {
1156 let cost = LoweredOp::Mul(Box::new(var(0)), Box::new(var(0)));
1157 let result = newton(&cost, &[], 10, 1e-8);
1158 assert!(matches!(
1159 result,
1160 Err(SymbolicNewtonError::DimMismatch { .. })
1161 ));
1162 }
1163
1164 #[test]
1165 fn solve_linear_2x2() {
1166 let a = vec![vec![2.0, 1.0], vec![1.0, 3.0]];
1168 let b = vec![5.0, 10.0];
1169 let x = solve_linear(&a, &b).expect("solve");
1170 assert!((x[0] - 1.0).abs() < 1e-10);
1171 assert!((x[1] - 3.0).abs() < 1e-10);
1172 }
1173
1174 #[test]
1175 fn solve_linear_singular_returns_err() {
1176 let a = vec![vec![1.0, 2.0], vec![2.0, 4.0]]; let b = vec![3.0, 6.0];
1178 assert!(matches!(
1179 solve_linear(&a, &b),
1180 Err(SymbolicNewtonError::SingularHessian)
1181 ));
1182 }
1183}