Skip to main content

pounce_cli/minima/
mod.rs

1//! Multistart / find-minima driver for the `pounce` CLI (`--minima`).
2//!
3//! A pure-Rust port of `pounce.find_minima` (`python/pounce/_minima.py`):
4//! drive the same local IPM solver in a loop, escaping already-found minima
5//! by one of six strategies, and collect the distinct local minima into a
6//! deduplicated archive. The strategies and their references:
7//!
8//! * `multistart` — random / Sobol' box sampling.
9//! * `mlsl` — Multi-Level Single Linkage clustering (Rinnooy Kan & Timmer 1987).
10//! * `basinhopping` — Metropolis chain over minima (Wales & Doye 1997).
11//! * `flooding` — repulsive Gaussian bumps (filled-function; Ge 1990).
12//! * `deflation` — softened `1/‖x−x*‖^p` poles (Farrell et al. 2015).
13//! * `tunneling` — equal-height tunnel between descents (Levy & Montalvo 1985).
14//!
15//! The local solver is reused across starts on a single `IpoptApplication`
16//! (no rebuild per start): each start wraps the base TNLP in a
17//! [`SeededTnlp`] (and, for the repulsion strategies, a penalty wrapper).
18//! Acceptance mirrors `_minima.py`: solve succeeded ∧ finite ∧ in-bounds ∧
19//! (Hessian PSD within `psd_tol`) ∧ not already in the archive.
20
21pub mod archive;
22pub mod penalty_tnlp;
23pub mod sampling;
24
25use crate::cli::{Args, MinimaArgs, MinimaMethod, ProblemSource};
26use crate::seeded_tnlp::SeededTnlp;
27use crate::solve_report::{status_to_solve_result_num, InputDescriptor, ReportBuilder};
28use archive::{scaled_distance, Archive};
29use penalty_tnlp::{Kernel, PenaltyTnlp, TunnelTnlp};
30use pounce_algorithm::application::IpoptApplication;
31use pounce_common::types::{Index, Number};
32use pounce_nlp::return_codes::ApplicationReturnStatus;
33use pounce_nlp::tnlp::{BoundsInfo, IndexStyle, SparsityRequest, StartingPoint, TNLP};
34use sampling::{clip, Sampler};
35use std::cell::RefCell;
36use std::path::{Path, PathBuf};
37use std::process::ExitCode;
38use std::rc::Rc;
39
40/// AMPL bound sentinel: a bound beyond this magnitude counts as ±∞.
41const BOUND_INF: Number = 1e19;
42/// Above this dimension a dense symmetric eigendecomposition (cyclic
43/// Jacobi, O(n³)) is too slow for the per-acceptance saddle-rejection
44/// check, so we skip it and accept (matching `find_minima`'s `hess=None`).
45const PSD_MAX_N: usize = 256;
46
47/// Why the search stopped (mirrors `MinimaResult.status`).
48#[derive(Clone, Copy, Debug)]
49enum Stop {
50    TargetReached,
51    Converged,
52    BudgetExhausted,
53}
54
55impl Stop {
56    fn as_str(self) -> &'static str {
57        match self {
58            Stop::TargetReached => "target_reached",
59            Stop::Converged => "converged",
60            Stop::BudgetExhausted => "budget_exhausted",
61        }
62    }
63}
64
65/// One local solve's outcome (the captured minimizer + whether it converged).
66struct SolveOut {
67    success: bool,
68    x: Vec<Number>,
69}
70
71/// The find-minima driver. Holds the single application and base problem and
72/// runs the chosen strategy until a [`Stop`].
73struct Driver<'a> {
74    app: &'a mut IpoptApplication,
75    base: Rc<RefCell<dyn TNLP>>,
76    /// Filled by the `on_converged` hook with the converged primal (full
77    /// length); cleared before each solve, taken after.
78    capture: Rc<RefCell<Option<Vec<Number>>>>,
79    cfg: &'a MinimaArgs,
80    n: usize,
81    m: usize,
82    nnz_h: usize,
83    /// 1 when the TNLP emits Fortran (1-based) triplet indices, else 0.
84    index_offset: usize,
85    x0: Vec<Number>,
86    x_l: Vec<Number>,
87    x_u: Vec<Number>,
88    has_box: bool,
89    /// Per-dimension scale `L` (box width, 1.0 for unbounded dims).
90    l_scale: Vec<Number>,
91    sampler: Sampler,
92    archive: Archive,
93    stagnant: usize,
94    n_solves: usize,
95    max_solves: usize,
96    /// Sampled points drawn so far (only MLSL counts against this).
97    n_samples: usize,
98    /// Hard ceiling on sampled points for solve-gated strategies (MLSL),
99    /// so `--max-solves` bounds wall-clock even when the clustering filter
100    /// rejects every sample (pounce#103).
101    max_samples: usize,
102    psd_skipped_logged: bool,
103}
104
105impl<'a> Driver<'a> {
106    // ---- shared local-solver ops -------------------------------------
107
108    /// Run one solve of `solve_tnlp`, toggling the Hessian mode, and return
109    /// the converged minimizer (captured via the `on_converged` hook).
110    /// Returns `Err(Stop::BudgetExhausted)` once the solve budget is spent.
111    fn run_solve(
112        &mut self,
113        solve_tnlp: Rc<RefCell<dyn TNLP>>,
114        exact_hessian: bool,
115    ) -> Result<SolveOut, Stop> {
116        if self.n_solves >= self.max_solves {
117            return Err(Stop::BudgetExhausted);
118        }
119        self.n_solves += 1;
120        // Penalty (repulsion) solves go quasi-Newton; clean / polish solves
121        // keep the exact Hessian. The IPM rereads this option each solve.
122        let line = if exact_hessian {
123            "hessian_approximation exact\n"
124        } else {
125            "hessian_approximation limited-memory\n"
126        };
127        let _ = self.app.options_mut().read_from_str(line, true);
128        *self.capture.borrow_mut() = None;
129        let status = self.app.optimize_tnlp(solve_tnlp);
130        let success = matches!(
131            status,
132            ApplicationReturnStatus::SolveSucceeded
133                | ApplicationReturnStatus::SolvedToAcceptableLevel
134        );
135        match self.capture.borrow_mut().take() {
136            Some(x) if success => Ok(SolveOut { success: true, x }),
137            // A failed solve has no usable captured point; acceptance needs
138            // success anyway, so the empty x is never read.
139            _ => Ok(SolveOut {
140                success: false,
141                x: Vec::new(),
142            }),
143        }
144    }
145
146    fn solve_seeded(&mut self, seed_x: &[Number], exact: bool) -> Result<SolveOut, Stop> {
147        let t: Rc<RefCell<dyn TNLP>> = Rc::new(RefCell::new(SeededTnlp::new(
148            Rc::clone(&self.base),
149            seed_x.to_vec(),
150        )));
151        self.run_solve(t, exact)
152    }
153
154    fn solve_penalty(&mut self, seed_x: &[Number], kernel: Kernel) -> Result<SolveOut, Stop> {
155        let pen: Rc<RefCell<dyn TNLP>> = Rc::new(RefCell::new(PenaltyTnlp::new(
156            Rc::clone(&self.base),
157            kernel,
158        )));
159        let t: Rc<RefCell<dyn TNLP>> = Rc::new(RefCell::new(SeededTnlp::new(pen, seed_x.to_vec())));
160        self.run_solve(t, false)
161    }
162
163    fn solve_tunnel(
164        &mut self,
165        seed_x: &[Number],
166        f_ref: Number,
167        pole: Kernel,
168    ) -> Result<SolveOut, Stop> {
169        let tun: Rc<RefCell<dyn TNLP>> = Rc::new(RefCell::new(TunnelTnlp::new(
170            Rc::clone(&self.base),
171            f_ref,
172            pole,
173        )));
174        let t: Rc<RefCell<dyn TNLP>> = Rc::new(RefCell::new(SeededTnlp::new(tun, seed_x.to_vec())));
175        self.run_solve(t, false)
176    }
177
178    /// Clean objective value at `x` (the un-augmented problem).
179    fn clean_f(&mut self, x: &[Number]) -> Option<Number> {
180        self.base.borrow_mut().eval_f(x, true)
181    }
182
183    /// Is `x` inside the box, allowing for the solver's bound relaxation?
184    ///
185    /// The IPM lets a converged primal sit slightly *outside* a bound — by up
186    /// to `bound_relax_factor · max(1, |bound|)` (the Ipopt default factor is
187    /// `1e-8`). On problems whose optimum binds a large-magnitude limit (e.g.
188    /// ACOPF generator/flow limits in the hundreds), that legal slack exceeds
189    /// any fixed absolute tolerance, so a purely absolute test wrongly rejects
190    /// every minimum (pounce#101). Use a bound-magnitude-relative tolerance
191    /// comfortably above the relaxation but far below any real basin spacing.
192    fn in_bounds(&self, x: &[Number]) -> bool {
193        x.iter()
194            .zip(&self.x_l)
195            .zip(&self.x_u)
196            .all(|((&xi, &lo), &hi)| coord_in_bounds(xi, lo, hi))
197    }
198
199    /// Dense objective Hessian at `x` (row-major n×n), or `None` when no
200    /// exact Hessian is available / the problem is too large.
201    fn obj_hessian_dense(&mut self, x: &[Number]) -> Option<Vec<Number>> {
202        if self.n > PSD_MAX_N || self.nnz_h == 0 {
203            return None;
204        }
205        let nnz = self.nnz_h;
206        let mut irow = vec![0 as Index; nnz];
207        let mut jcol = vec![0 as Index; nnz];
208        {
209            let mut b = self.base.borrow_mut();
210            if !b.eval_h(
211                None,
212                false,
213                1.0,
214                None,
215                false,
216                SparsityRequest::Structure {
217                    irow: &mut irow,
218                    jcol: &mut jcol,
219                },
220            ) {
221                return None;
222            }
223        }
224        let lam = vec![0.0; self.m];
225        let mut vals = vec![0.0; nnz];
226        {
227            let mut b = self.base.borrow_mut();
228            if !b.eval_h(
229                Some(x),
230                true,
231                1.0,
232                Some(&lam),
233                true,
234                SparsityRequest::Values { values: &mut vals },
235            ) {
236                return None;
237            }
238        }
239        let n = self.n;
240        let mut dense = vec![0.0; n * n];
241        for k in 0..nnz {
242            let i = irow[k] as usize - self.index_offset;
243            let j = jcol[k] as usize - self.index_offset;
244            dense[i * n + j] += vals[k];
245            if i != j {
246                dense[j * n + i] += vals[k];
247            }
248        }
249        Some(dense)
250    }
251
252    /// Reject saddles/maxima via the clean Hessian's smallest eigenvalue
253    /// (accept when no Hessian is available — matching `find_minima`).
254    fn is_minimum(&mut self, x: &[Number]) -> bool {
255        if self.n > PSD_MAX_N {
256            if !self.psd_skipped_logged {
257                eprintln!(
258                    "pounce: --minima saddle-rejection (PSD) check skipped — n={} exceeds the \
259                     dense-eigendecomposition cap ({PSD_MAX_N}); accepting converged points as minima.",
260                    self.n
261                );
262                self.psd_skipped_logged = true;
263            }
264            return true;
265        }
266        let dense = match self.obj_hessian_dense(x) {
267            Some(d) => d,
268            None => return true,
269        };
270        let n = self.n;
271        let mut w = vec![0.0; n];
272        let mut v = vec![0.0; n * n];
273        if !pounce_sensitivity::symmetric_eigen(&dense, n, &mut w, &mut v) {
274            return true;
275        }
276        let min_eig = w.iter().cloned().fold(f64::INFINITY, f64::min);
277        min_eig >= -self.cfg.psd_tol
278    }
279
280    /// Per-dimension width vector from a knob spec, mirroring
281    /// `_resolve_lengths`: a scalar ⇒ isotropic; `None` ("auto") ⇒
282    /// `frac · L` when a box is known, else `fallback`.
283    fn resolve_lengths(
284        &self,
285        spec: Option<f64>,
286        frac_default: f64,
287        frac_override: Option<f64>,
288        fallback: f64,
289    ) -> Vec<Number> {
290        match spec {
291            Some(s) => vec![s; self.n],
292            None => {
293                let frac = frac_override.unwrap_or(frac_default);
294                if self.has_box {
295                    self.l_scale.iter().map(|&l| frac * l).collect()
296                } else {
297                    vec![fallback; self.n]
298                }
299            }
300        }
301    }
302
303    /// Curvature-based escape height for a flooding bump at `center`
304    /// (`margin · μ_min` of `diag(σ)·H·diag(σ)`); `None` when no Hessian.
305    fn auto_amplitude(&mut self, center: &[Number], sigma: &[Number], margin: f64) -> Option<f64> {
306        let h = self.obj_hessian_dense(center)?;
307        let n = self.n;
308        let mut s_mat = vec![0.0; n * n];
309        for i in 0..n {
310            for j in 0..n {
311                s_mat[i * n + j] = sigma[i] * sigma[j] * h[i * n + j];
312            }
313        }
314        let mut w = vec![0.0; n];
315        let mut v = vec![0.0; n * n];
316        if !pounce_sensitivity::symmetric_eigen(&s_mat, n, &mut w, &mut v) {
317            return None;
318        }
319        let mu_min = w.iter().cloned().fold(f64::INFINITY, f64::min);
320        Some(margin * mu_min.max(1e-12))
321    }
322
323    /// Draw a fresh start from the box (Sobol'/uniform) or jitter around x0.
324    fn sample(&mut self, jitter: f64) -> Vec<Number> {
325        let x0 = self.x0.clone();
326        let lo = self.x_l.clone();
327        let hi = self.x_u.clone();
328        self.sampler.sample(&x0, &lo, &hi, self.has_box, jitter)
329    }
330
331    // ---- acceptance --------------------------------------------------
332
333    /// Consider a candidate for the archive. With `polish`, first re-solve
334    /// the clean objective from the candidate (exact Hessian) — the
335    /// repulsion strategies escape on the augmented objective, then polish
336    /// back onto the true one. Returns whether it was accepted.
337    fn consider(
338        &mut self,
339        mut x: Vec<Number>,
340        mut success: bool,
341        polish: bool,
342    ) -> Result<bool, Stop> {
343        if success && polish {
344            let r = self.solve_seeded(&x, true)?;
345            success = r.success;
346            if success {
347                x = r.x;
348            }
349        }
350        if !success {
351            return self.reject();
352        }
353        let fval = match self.clean_f(&x) {
354            Some(f) => f,
355            None => return self.reject(),
356        };
357        let finite = x.iter().all(|v| v.is_finite()) && fval.is_finite();
358        let accepted =
359            finite && self.in_bounds(&x) && self.is_minimum(&x) && !self.archive.is_known(&x);
360        if accepted {
361            self.archive.add(x, fval);
362            self.stagnant = 0;
363            if self.archive.len() >= self.cfg.n_minima {
364                return Err(Stop::TargetReached);
365            }
366            Ok(true)
367        } else {
368            self.reject()
369        }
370    }
371
372    fn reject(&mut self) -> Result<bool, Stop> {
373        self.stagnant += 1;
374        if self.stagnant >= self.cfg.patience {
375            return Err(Stop::Converged);
376        }
377        Ok(false)
378    }
379
380    /// Count one drawn sample against the sampling budget. MLSL's expensive
381    /// work is *sampling* (an O(N²) single-linkage scan over a growing pool),
382    /// not solving, so on a problem where the clustering filter rejects
383    /// almost every sample no solve ever fires and `max_solves` cannot bound
384    /// the loop (pounce#103). The sample budget gives it a hard ceiling.
385    fn note_sample(&mut self) -> Result<(), Stop> {
386        if self.n_samples >= self.max_samples {
387            return Err(Stop::BudgetExhausted);
388        }
389        self.n_samples += 1;
390        Ok(())
391    }
392
393    // ---- strategy loops (each runs until a Stop) ---------------------
394
395    fn run(&mut self) -> Stop {
396        let res = match self.cfg.method {
397            MinimaMethod::Multistart => self.run_multistart(),
398            MinimaMethod::Mlsl => self.run_mlsl(),
399            MinimaMethod::Basinhopping => self.run_basinhopping(),
400            MinimaMethod::Flooding => self.run_flooding(),
401            MinimaMethod::Deflation => self.run_deflation(),
402            MinimaMethod::Tunneling => self.run_tunneling(),
403        };
404        // The loops only exit by propagating a `Stop`; `Ok` is unreachable.
405        res.err().unwrap_or(Stop::BudgetExhausted)
406    }
407
408    fn run_multistart(&mut self) -> Result<(), Stop> {
409        let jitter = self.cfg.restart_jitter.unwrap_or(1.0);
410        let x0 = self.x0.clone();
411        let r = self.solve_seeded(&x0, true)?;
412        self.consider(r.x, r.success, false)?;
413        loop {
414            let s = self.sample(jitter);
415            let r = self.solve_seeded(&s, true)?;
416            self.consider(r.x, r.success, false)?;
417        }
418    }
419
420    fn run_mlsl(&mut self) -> Result<(), Stop> {
421        let batch = self.cfg.samples_per_round.unwrap_or(20);
422        let gamma = self.cfg.gamma.unwrap_or(2.0);
423        let jitter = self.cfg.restart_jitter.unwrap_or(1.0);
424        let n = self.n;
425        let diag = (n as f64).sqrt();
426        let mut pool_x: Vec<Vec<Number>> = Vec::new();
427        let mut pool_f: Vec<Number> = Vec::new();
428        let x0 = self.x0.clone();
429        let r = self.solve_seeded(&x0, true)?;
430        self.consider(r.x, r.success, false)?;
431        loop {
432            // Grow the pool; each draw counts against the sample budget, so a
433            // round that solves nothing still drives the loop to terminate.
434            for _ in 0..batch {
435                self.note_sample()?;
436                let s = self.sample(jitter);
437                let f = self.clean_f(&s).unwrap_or(f64::INFINITY);
438                pool_x.push(s);
439                pool_f.push(f);
440            }
441            let bign = pool_x.len();
442            let ne = bign.max(2) as f64;
443            let radius = gamma * diag * (ne.ln() / ne).powf(1.0 / n as f64);
444            let mut order: Vec<usize> = (0..bign).collect();
445            order.sort_by(|&a, &b| {
446                pool_f[a]
447                    .partial_cmp(&pool_f[b])
448                    .unwrap_or(std::cmp::Ordering::Equal)
449            });
450            for i in order {
451                let si = pool_x[i].clone();
452                let fi = pool_f[i];
453                // Single-linkage: skip if a *better* sample is within radius.
454                let better_near = (0..bign).any(|j| {
455                    j != i
456                        && pool_f[j] < fi
457                        && scaled_distance(&si, &pool_x[j], &self.l_scale) < radius
458                });
459                if better_near || self.archive.near_any(&si, radius) {
460                    continue;
461                }
462                let r = self.solve_seeded(&si, true)?;
463                self.consider(r.x, r.success, false)?;
464            }
465        }
466    }
467
468    fn run_basinhopping(&mut self) -> Result<(), Stop> {
469        let step = self.cfg.step.unwrap_or(0.5);
470        let temperature = self.cfg.temperature.unwrap_or(1.0);
471        let x0 = self.x0.clone();
472        let r = self.solve_seeded(&x0, true)?;
473        let mut cur = if r.success { r.x.clone() } else { x0.clone() };
474        let mut cur_f = self.clean_f(&cur).unwrap_or(f64::INFINITY);
475        self.consider(r.x, r.success, false)?;
476        loop {
477            let mut trial = self.sampler.perturb(&cur, &[step]);
478            clip(&mut trial, &self.x_l, &self.x_u, self.has_box);
479            let r = self.solve_seeded(&trial, true)?;
480            if !r.success {
481                self.consider(r.x, false, false)?;
482                continue;
483            }
484            let new_x = r.x.clone();
485            self.consider(r.x, true, false)?;
486            let new_f = self.clean_f(&new_x).unwrap_or(f64::INFINITY);
487            let accept_uphill = self.sampler.uniform() < (-(new_f - cur_f) / temperature).exp();
488            if new_f < cur_f || accept_uphill {
489                cur = new_x;
490                cur_f = new_f;
491            }
492        }
493    }
494
495    fn run_flooding(&mut self) -> Result<(), Stop> {
496        let sigma = self.resolve_lengths(self.cfg.sigma, 0.1, self.cfg.sigma_frac, 0.5);
497        let inv_sigma2: Vec<f64> = sigma.iter().map(|&s| 1.0 / (s * s)).collect();
498        let amp_spec = self.cfg.amplitude;
499        let margin = self.cfg.amp_margin.unwrap_or(2.0);
500        let bump_factor = 3.0;
501        let bump_cap = 1e3;
502        let fallback_amp = 2.0;
503        let jitter = self.cfg.restart_jitter.unwrap_or(0.5);
504        let x0 = self.x0.clone();
505        let mut base_amp: Vec<f64> = Vec::new();
506        let mut mult: Vec<f64> = Vec::new();
507        let mut start = x0.clone();
508        let mut last_center: Option<usize> = None;
509        let mut fails = 0usize;
510        loop {
511            let centers = self.archive.xs.clone();
512            while base_amp.len() < centers.len() {
513                let k = base_amp.len();
514                let a = match amp_spec {
515                    Some(a) => a,
516                    None => self
517                        .auto_amplitude(&centers[k], &sigma, margin)
518                        .unwrap_or(fallback_amp),
519                };
520                base_amp.push(a);
521                mult.push(1.0);
522            }
523            let eff: Vec<f64> = (0..centers.len()).map(|k| base_amp[k] * mult[k]).collect();
524            let polish = !centers.is_empty();
525            let solve_out = if centers.is_empty() {
526                self.solve_seeded(&start, true)?
527            } else {
528                let kernel = Kernel::Gauss {
529                    centers: centers.clone(),
530                    amps: eff,
531                    inv_sigma2: inv_sigma2.clone(),
532                };
533                self.solve_penalty(&start, kernel)?
534            };
535            let accepted = self.consider(solve_out.x, solve_out.success, polish)?;
536            if accepted {
537                if let Some(last) = self.archive.xs.last() {
538                    start = last.clone();
539                }
540                last_center = Some(self.archive.xs.len() - 1);
541                fails = 0;
542            } else if let Some(lc) = last_center {
543                if mult[lc] < bump_cap && fails < 8 {
544                    // Under-flooded the basin we started from: bump and retry.
545                    mult[lc] *= bump_factor;
546                    let scale: Vec<f64> = sigma.iter().map(|&s| 0.05 * s).collect();
547                    start = self.sampler.perturb(&centers[lc], &scale);
548                    fails += 1;
549                    continue;
550                }
551                start = self.sample(jitter);
552                last_center = None;
553                fails = 0;
554            } else {
555                start = self.sample(jitter);
556                last_center = None;
557                fails = 0;
558            }
559        }
560    }
561
562    fn run_deflation(&mut self) -> Result<(), Stop> {
563        let eta = self.cfg.eta.unwrap_or(1.0);
564        let power = self.cfg.power.unwrap_or(2.0);
565        let soft = self.cfg.soft.unwrap_or(1e-3);
566        let length = self.resolve_lengths(self.cfg.length, 0.1, self.cfg.length_frac, 0.5);
567        let inv_len2: Vec<f64> = length.iter().map(|&l| 1.0 / (l * l)).collect();
568        let jitter = self.cfg.restart_jitter.unwrap_or(0.5);
569        let q = power / 2.0;
570        let mut start = self.x0.clone();
571        loop {
572            let centers = self.archive.xs.clone();
573            let polish = !centers.is_empty();
574            // Step a little off the pole so the first gradient is finite.
575            let mut s = start.clone();
576            if !centers.is_empty() && self.archive.is_known(&s) {
577                let scale: Vec<f64> = length.iter().map(|&l| 0.1 * l).collect();
578                s = self.sampler.perturb(&s, &scale);
579            }
580            let solve_out = if centers.is_empty() {
581                self.solve_seeded(&s, true)?
582            } else {
583                let kernel = Kernel::Pole {
584                    centers: centers.clone(),
585                    eta,
586                    q,
587                    soft,
588                    inv_len2: inv_len2.clone(),
589                };
590                self.solve_penalty(&s, kernel)?
591            };
592            let accepted = self.consider(solve_out.x, solve_out.success, polish)?;
593            if accepted {
594                if let Some(last) = self.archive.xs.last() {
595                    start = last.clone();
596                }
597            } else {
598                start = self.sample(jitter);
599            }
600        }
601    }
602
603    fn run_tunneling(&mut self) -> Result<(), Stop> {
604        let eta = self.cfg.eta.unwrap_or(1.0);
605        let power = self.cfg.power.unwrap_or(2.0);
606        let soft = self.cfg.soft.unwrap_or(1e-3);
607        let length = self.resolve_lengths(self.cfg.length, 0.1, self.cfg.length_frac, 0.5);
608        let inv_len2: Vec<f64> = length.iter().map(|&l| 1.0 / (l * l)).collect();
609        let jitter = self.cfg.restart_jitter.unwrap_or(0.75);
610        let q = power / 2.0;
611        let x0 = self.x0.clone();
612        // Seed: one clean descent.
613        let r = self.solve_seeded(&x0, true)?;
614        self.consider(r.x, r.success, false)?;
615        loop {
616            let centers = self.archive.xs.clone();
617            // Tunnel at the height of the most-recent minimum, away from all
618            // known minima — the classic monotone-descending tunnel.
619            let f_ref = match self.archive.fs.last() {
620                Some(&f) => f,
621                None => self.clean_f(&x0).unwrap_or(0.0),
622            };
623            let anchor = self
624                .archive
625                .xs
626                .last()
627                .cloned()
628                .unwrap_or_else(|| x0.clone());
629            let jit = vec![jitter; self.n];
630            let mut start = self.sampler.perturb(&anchor, &jit);
631            clip(&mut start, &self.x_l, &self.x_u, self.has_box);
632            let kernel = Kernel::Pole {
633                centers: centers.clone(),
634                eta,
635                q,
636                soft,
637                inv_len2: inv_len2.clone(),
638            };
639            let r = self.solve_tunnel(&start, f_ref, kernel)?;
640            self.consider(r.x, r.success, true)?;
641        }
642    }
643}
644
645/// A single found minimum (used for output / JSON).
646struct Minimum {
647    x: Vec<Number>,
648    objective: Number,
649}
650
651/// Entry point: run the `--minima` search on `base` (the raw problem TNLP —
652/// presolve / counting wrappers are intentionally bypassed so coordinates
653/// match the original problem and the clean objective is evaluated directly).
654/// Returns the process exit code.
655pub fn run(
656    app: &mut IpoptApplication,
657    base: &Rc<RefCell<dyn TNLP>>,
658    cfg: &MinimaArgs,
659    args: &Args,
660    sol_path: Option<&Path>,
661) -> ExitCode {
662    let info = match base.borrow_mut().get_nlp_info() {
663        Some(i) => i,
664        None => {
665            eprintln!("pounce: --minima could not read problem dimensions");
666            return ExitCode::from(2);
667        }
668    };
669    let n = info.n as usize;
670    let m = info.m as usize;
671    let nnz_h = info.nnz_h_lag as usize;
672    let index_offset = match info.index_style {
673        IndexStyle::Fortran => 1,
674        IndexStyle::C => 0,
675    };
676
677    // Bounds + starting point straight from the TNLP (works uniformly for
678    // built-ins and `.nl` files).
679    let mut x_l = vec![0.0; n];
680    let mut x_u = vec![0.0; n];
681    let mut g_l = vec![0.0; m];
682    let mut g_u = vec![0.0; m];
683    base.borrow_mut().get_bounds_info(BoundsInfo {
684        x_l: &mut x_l,
685        x_u: &mut x_u,
686        g_l: &mut g_l,
687        g_u: &mut g_u,
688    });
689    let mut x0 = vec![0.0; n];
690    {
691        let mut z_l = vec![0.0; n];
692        let mut z_u = vec![0.0; n];
693        let mut lambda = vec![0.0; m];
694        base.borrow_mut().get_starting_point(StartingPoint {
695            init_x: true,
696            x: &mut x0,
697            init_z: false,
698            z_l: &mut z_l,
699            z_u: &mut z_u,
700            init_lambda: false,
701            lambda: &mut lambda,
702        });
703    }
704
705    // Per-dimension scale and box availability (mirror `_scale_from_bounds`).
706    let has_box = (0..n).all(|i| x_l[i] > -BOUND_INF && x_u[i] < BOUND_INF);
707    let l_scale: Vec<Number> = (0..n)
708        .map(|i| {
709            let w = x_u[i] - x_l[i];
710            if has_box && w > 0.0 {
711                w
712            } else {
713                1.0
714            }
715        })
716        .collect();
717
718    // Capture the converged primal of each solve via the on-converged hook.
719    let capture: Rc<RefCell<Option<Vec<Number>>>> = Rc::new(RefCell::new(None));
720    {
721        let cap = Rc::clone(&capture);
722        app.set_on_converged(Box::new(move |data, _cq, nlp, _pd| {
723            if let Some(curr) = data.borrow().curr.clone() {
724                let x = nlp.borrow().lift_x_to_full(&*curr.x);
725                *cap.borrow_mut() = Some(x);
726            }
727        }));
728    }
729
730    let max_solves = cfg.max_solves.unwrap_or(8 * cfg.n_minima);
731    // Sample ceiling for solve-gated strategies (MLSL): one round of samples
732    // per unit of solve budget. The patience-on-stall rule normally
733    // terminates first; this guarantees `--max-solves` bounds wall-clock even
734    // when the clustering filter rejects everything (pounce#103).
735    let batch = cfg.samples_per_round.unwrap_or(20).max(1);
736    let max_samples = max_solves.saturating_mul(batch);
737
738    println!(
739        "Searching for up to {} minima via `{}` (max {} solves, seed {})...",
740        cfg.n_minima,
741        cfg.method.as_str(),
742        max_solves,
743        cfg.seed
744    );
745
746    let mut driver = Driver {
747        app,
748        base: Rc::clone(base),
749        capture,
750        cfg,
751        n,
752        m,
753        nnz_h,
754        index_offset,
755        x0,
756        x_l: x_l.clone(),
757        x_u: x_u.clone(),
758        has_box,
759        l_scale: l_scale.clone(),
760        sampler: Sampler::new(cfg.seed, cfg.sobol),
761        archive: Archive::new(cfg.dedup, l_scale.clone()),
762        stagnant: 0,
763        n_solves: 0,
764        max_solves,
765        n_samples: 0,
766        max_samples,
767        psd_skipped_logged: false,
768    };
769
770    let stop = driver.run();
771    let n_solves = driver.n_solves;
772
773    // Rank the found minima by objective (best first).
774    let order = driver.archive.order_by_objective();
775    let minima: Vec<Minimum> = order
776        .iter()
777        .map(|&i| Minimum {
778            x: driver.archive.xs[i].clone(),
779            objective: driver.archive.fs[i],
780        })
781        .collect();
782    let best_obj = order.first().map(|&i| driver.archive.fs[i]);
783
784    print_table(&minima, &l_scale, stop, n_solves);
785
786    // Write the per-minimum `.sol` files: best → <stub>.sol, the rest →
787    // ranked siblings <stub>.minNNN.sol.
788    if let Some(sp) = sol_path {
789        write_sol_files(sp, &minima, m);
790    }
791
792    // JSON report: the standard single-solve report for the best minimum,
793    // plus a backward-compatible `minima` section listing all of them.
794    if let Some(json_path) = &args.json_output {
795        write_json_report(json_path, args, cfg, stop, n_solves, &minima, n, m, &info);
796    }
797
798    if best_obj.is_some() {
799        ExitCode::SUCCESS
800    } else {
801        ExitCode::from(1)
802    }
803}
804
805/// Print a ranked console table of the distinct minima found.
806fn print_table(minima: &[Minimum], l_scale: &[Number], stop: Stop, n_solves: usize) {
807    println!();
808    println!(
809        "find-minima: {} distinct minim{} in {} solves ({})",
810        minima.len(),
811        if minima.len() == 1 { "um" } else { "a" },
812        n_solves,
813        stop.as_str()
814    );
815    if minima.is_empty() {
816        println!("  (no accepted minima — try raising --max-solves or --patience)");
817        return;
818    }
819    println!("  rank        objective     dist-to-best");
820    let best = &minima[0].x;
821    for (rank, mn) in minima.iter().enumerate() {
822        let d = scaled_distance(&mn.x, best, l_scale);
823        println!("  {rank:>4}   {:>16.8e}   {:>14.6e}", mn.objective, d);
824    }
825}
826
827/// Write `.sol` files: best to `sol_path`, ranked siblings alongside.
828fn write_sol_files(sol_path: &Path, minima: &[Minimum], m: usize) {
829    let zeros = vec![0.0; m];
830    for (rank, mn) in minima.iter().enumerate() {
831        let path = if rank == 0 {
832            sol_path.to_path_buf()
833        } else {
834            sibling_sol_path(sol_path, rank)
835        };
836        let message = format!(
837            "POUNCE {} find-minima rank {rank}: Solve_Succeeded",
838            env!("CARGO_PKG_VERSION")
839        );
840        let payload = crate::nl_writer::SolutionFile {
841            message: &message,
842            x: &mn.x,
843            lambda: &zeros,
844            solve_result_num: status_to_solve_result_num(ApplicationReturnStatus::SolveSucceeded),
845            suffixes: &[],
846        };
847        match crate::nl_writer::write_sol_file(&path, &payload) {
848            Ok(_) => eprintln!("pounce: wrote {}", path.display()),
849            Err(e) => eprintln!("pounce: failed to write {}: {e}", path.display()),
850        }
851    }
852}
853
854/// `<stub>.sol` → `<stub>.minNNN.sol` for rank ≥ 1.
855fn sibling_sol_path(sol_path: &Path, rank: usize) -> PathBuf {
856    let mut stub = sol_path.to_path_buf();
857    stub.set_extension(""); // drop `.sol`
858    let base = stub.to_string_lossy().into_owned();
859    PathBuf::from(format!("{base}.min{rank:03}.sol"))
860}
861
862/// Build the JSON report (standard best-solution report + `minima` section)
863/// and write it.
864#[allow(clippy::too_many_arguments)]
865fn write_json_report(
866    json_path: &Path,
867    args: &Args,
868    cfg: &MinimaArgs,
869    stop: Stop,
870    n_solves: usize,
871    minima: &[Minimum],
872    n: usize,
873    m: usize,
874    info: &pounce_nlp::tnlp::NlpInfo,
875) {
876    let input = match &args.problem {
877        ProblemSource::Builtin(name) => InputDescriptor::Builtin { name: name.clone() },
878        ProblemSource::NlFile(p) => InputDescriptor::NlFile {
879            path: p.clone(),
880            size_bytes: std::fs::metadata(p).ok().map(|md| md.len()),
881        },
882    };
883    let mut builder = ReportBuilder::new(args.json_detail, input);
884    builder.problem.n_variables = n as Index;
885    builder.problem.n_constraints = m as Index;
886    builder.problem.n_objectives = 1;
887    builder.problem.nnz_jac_g = Some(info.nnz_jac_g);
888    builder.problem.nnz_h_lag = Some(info.nnz_h_lag);
889    if let Some(best) = minima.first() {
890        builder.solution.status = ApplicationReturnStatus::SolveSucceeded;
891        builder.solution.solve_result_num =
892            status_to_solve_result_num(ApplicationReturnStatus::SolveSucceeded);
893        builder.solution.objective = best.objective;
894        builder.solution.x = best.x.clone();
895        builder.solution.lambda = vec![0.0; m];
896    }
897    let report = builder.finish();
898
899    // Inject the `minima` section without a schema change: serialize, then
900    // splice it into the top-level object.
901    let mut value = match serde_json::to_value(&report) {
902        Ok(v) => v,
903        Err(e) => {
904            eprintln!("pounce: failed to serialize minima report: {e}");
905            return;
906        }
907    };
908    let minima_json: Vec<serde_json::Value> = minima
909        .iter()
910        .map(|mn| {
911            serde_json::json!({
912                "x": mn.x,
913                "objective": mn.objective,
914            })
915        })
916        .collect();
917    let values: Vec<Number> = minima.iter().map(|mn| mn.objective).collect();
918    if let serde_json::Value::Object(map) = &mut value {
919        map.insert(
920            "minima".to_string(),
921            serde_json::json!({
922                "method": cfg.method.as_str(),
923                "status": stop.as_str(),
924                "n_solves": n_solves,
925                "n_minima": minima.len(),
926                "minima": minima_json,
927                "values": values,
928            }),
929        );
930    }
931    match serde_json::to_string_pretty(&value) {
932        Ok(s) => match std::fs::write(json_path, s) {
933            Ok(_) => eprintln!("pounce: wrote {}", json_path.display()),
934            Err(e) => eprintln!(
935                "pounce: failed to write JSON report to {}: {e}",
936                json_path.display()
937            ),
938        },
939        Err(e) => eprintln!("pounce: failed to render minima report: {e}"),
940    }
941}
942
943/// Per-coordinate box test with a bound-magnitude-relative tolerance that
944/// absorbs the IPM's bound relaxation (`bound_relax_factor·max(1,|bound|)`,
945/// Ipopt default factor `1e-8`). A purely absolute tolerance rejects minima
946/// that legally bind large-magnitude limits (pounce#101); the relative term
947/// tracks the relaxation while staying far below any real basin spacing.
948fn coord_in_bounds(xi: Number, lo: Number, hi: Number) -> bool {
949    const ATOL: Number = 1e-9;
950    const RTOL: Number = 1e-6;
951    let tol_lo = ATOL + RTOL * lo.abs().max(1.0);
952    let tol_hi = ATOL + RTOL * hi.abs().max(1.0);
953    xi >= lo - tol_lo && xi <= hi + tol_hi
954}
955
956#[cfg(test)]
957mod tests {
958    use super::coord_in_bounds;
959
960    #[test]
961    fn accepts_interior_point() {
962        assert!(coord_in_bounds(0.0, -1.0, 1.0));
963        assert!(coord_in_bounds(250.0, 0.0, 500.0));
964    }
965
966    #[test]
967    fn accepts_bound_relaxed_point_at_large_magnitude() {
968        // A converged primal may sit ~bound_relax_factor·|bound| (≈5e-6 at a
969        // 500-unit limit) past the bound; that point is a legal minimum.
970        assert!(coord_in_bounds(500.000005, 0.0, 500.0));
971        assert!(coord_in_bounds(-500.000005, -500.0, 0.0));
972    }
973
974    #[test]
975    fn rejects_genuinely_outside_point() {
976        assert!(!coord_in_bounds(500.1, 0.0, 500.0));
977        assert!(!coord_in_bounds(-1.1, -1.0, 1.0));
978    }
979}