argmin::core

Trait Solver

Source
pub trait Solver<O: ArgminOp>: Serialize {
    const NAME: &'static str = "UNDEFINED";

    // Required method
    fn next_iter(
        &mut self,
        op: &mut OpWrapper<O>,
        state: &IterState<O>,
    ) -> Result<ArgminIterData<O>, Error>;

    // Provided methods
    fn init(
        &mut self,
        _op: &mut OpWrapper<O>,
        _state: &IterState<O>,
    ) -> Result<Option<ArgminIterData<O>>, Error> { ... }
    fn terminate_internal(&mut self, state: &IterState<O>) -> TerminationReason { ... }
    fn terminate(&mut self, _state: &IterState<O>) -> TerminationReason { ... }
}
Expand description

Solver

Every solver needs to implement this trait.

Provided Associated Constants§

Source

const NAME: &'static str = "UNDEFINED"

Name of the solver

Required Methods§

Source

fn next_iter( &mut self, op: &mut OpWrapper<O>, state: &IterState<O>, ) -> Result<ArgminIterData<O>, Error>

Computes one iteration of the algorithm.

Provided Methods§

Source

fn init( &mut self, _op: &mut OpWrapper<O>, _state: &IterState<O>, ) -> Result<Option<ArgminIterData<O>>, Error>

Initializes the algorithm

This is executed before any iterations are performed. It can be used to perform precomputations. The default implementation corresponds to doing nothing.

Source

fn terminate_internal(&mut self, state: &IterState<O>) -> TerminationReason

Checks whether basic termination reasons apply.

Terminate if

  1. algorithm was terminated somewhere else in the Executor
  2. iteration count exceeds maximum number of iterations
  3. cost is lower than target cost

This can be overwritten in a Solver implementation; however it is not advised.

Source

fn terminate(&mut self, _state: &IterState<O>) -> TerminationReason

Checks whether the algorithm must be terminated

Dyn Compatibility§

This trait is not dyn compatible.

In older versions of Rust, dyn compatibility was called "object safety", so this trait is not object safe.

Implementors§

Source§

impl<O, B, R, F> Solver<O> for SR1TrustRegion<B, R, F>

Source§

const NAME: &'static str = "SR1 Trust Region"

Source§

impl<O, F> Solver<O> for Brent<F>
where O: ArgminOp<Param = F, Output = F, Float = F>, F: ArgminFloat,

Source§

const NAME: &'static str = "Brent"

Source§

impl<O, F> Solver<O> for GaussNewton<F>

Source§

const NAME: &'static str = "Gauss-Newton method"

Source§

impl<O, F> Solver<O> for GoldenSectionSearch<F>
where O: ArgminOp<Output = F, Param = F, Float = F>, F: ArgminFloat,

Source§

const NAME: &'static str = "Golden-section search"

Source§

impl<O, F> Solver<O> for Landweber<F>
where O: ArgminOp<Float = F>, O::Param: ArgminScaledSub<O::Param, O::Float, O::Param>, F: ArgminFloat,

Source§

const NAME: &'static str = "Landweber"

Source§

impl<O, F> Solver<O> for Newton<F>
where O: ArgminOp<Float = F>, O::Param: ArgminScaledSub<O::Param, O::Float, O::Param>, O::Hessian: ArgminInv<O::Hessian> + ArgminDot<O::Param, O::Param>, F: ArgminFloat,

Source§

const NAME: &'static str = "Newton method"

Source§

impl<O, F> Solver<O> for SimulatedAnnealing<F>
where O: ArgminOp<Output = F, Float = F>, F: ArgminFloat,

Source§

const NAME: &'static str = "Simulated Annealing"

Source§

impl<O, F> Solver<O> for CauchyPoint<F>
where O: ArgminOp<Output = F, Float = F>, O::Param: Debug + Clone + Serialize + ArgminMul<O::Float, O::Param> + ArgminWeightedDot<O::Param, F, O::Hessian> + ArgminNorm<O::Float>, O::Hessian: Clone + Serialize, F: ArgminFloat,

Source§

const NAME: &'static str = "Cauchy Point"

Source§

impl<O, F> Solver<O> for Dogleg<F>
where O: ArgminOp<Output = F, Float = F>, O::Param: Debug + ArgminMul<F, O::Param> + ArgminWeightedDot<O::Param, O::Float, O::Hessian> + ArgminNorm<F> + ArgminDot<O::Param, O::Float> + ArgminAdd<O::Param, O::Param> + ArgminSub<O::Param, O::Param>, O::Hessian: ArgminInv<O::Hessian> + ArgminDot<O::Param, O::Param>, F: ArgminFloat,

Source§

const NAME: &'static str = "Dogleg"

Source§

impl<O, L, F> Solver<O> for GaussNewtonLS<L, F>
where O: ArgminOp<Float = F>, O::Param: Default + Debug + ArgminScaledSub<O::Param, O::Float, O::Param> + ArgminSub<O::Param, O::Param> + ArgminMul<O::Float, O::Param>, O::Output: ArgminNorm<O::Float>, O::Jacobian: ArgminTranspose + ArgminInv<O::Jacobian> + ArgminDot<O::Jacobian, O::Jacobian> + ArgminDot<O::Output, O::Param> + ArgminDot<O::Param, O::Param>, O::Hessian: Default, L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<LineSearchOP<O>>>, F: ArgminFloat,

Source§

const NAME: &'static str = "Gauss-Newton method with Linesearch"

Source§

impl<O, L, F> Solver<O> for SteepestDescent<L>
where O: ArgminOp<Output = F, Float = F>, O::Param: Clone + Default + Serialize + ArgminSub<O::Param, O::Param> + ArgminDot<O::Param, O::Float> + ArgminScaledAdd<O::Param, O::Float, O::Param> + ArgminMul<O::Float, O::Param> + ArgminNorm<O::Float>, O::Hessian: Default, L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<O>>, F: ArgminFloat,

Source§

const NAME: &'static str = "Steepest Descent"

Source§

impl<O, L, F> Solver<O> for NewtonCG<L, F>

Source§

const NAME: &'static str = "Newton-CG"

Source§

impl<O, L, H, F> Solver<O> for BFGS<L, H, F>

Source§

const NAME: &'static str = "BFGS"

Source§

impl<O, L, H, F> Solver<O> for DFP<L, H, F>

Source§

const NAME: &'static str = "DFP"

Source§

impl<O, L, H, F> Solver<O> for SR1<L, H, F>
where O: ArgminOp<Output = F, Hessian = H, Float = F>, O::Param: Debug + Clone + Default + Serialize + ArgminSub<O::Param, O::Param> + ArgminDot<O::Param, O::Float> + ArgminDot<O::Param, O::Hessian> + ArgminNorm<O::Float> + ArgminMul<O::Float, O::Param>, O::Hessian: Debug + Clone + Default + Serialize + DeserializeOwned + ArgminSub<O::Hessian, O::Hessian> + ArgminDot<O::Param, O::Param> + ArgminDot<O::Hessian, O::Hessian> + ArgminAdd<O::Hessian, O::Hessian> + ArgminMul<F, O::Hessian>, L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<O>>, F: ArgminFloat,

Source§

const NAME: &'static str = "SR1"

Source§

impl<O, L, P, F> Solver<O> for LBFGS<L, P, F>
where O: ArgminOp<Param = P, Output = F, Float = F>, O::Param: Clone + Serialize + DeserializeOwned + Debug + Default + ArgminSub<O::Param, O::Param> + ArgminAdd<O::Param, O::Param> + ArgminDot<O::Param, O::Float> + ArgminScaledAdd<O::Param, O::Float, O::Param> + ArgminNorm<O::Float> + ArgminMul<O::Float, O::Param>, O::Hessian: Clone + Default + Serialize + DeserializeOwned, L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<O>>, F: ArgminFloat,

Source§

const NAME: &'static str = "L-BFGS"

Source§

impl<O, P, F> Solver<O> for NelderMead<P, F>
where O: ArgminOp<Output = F, Param = P, Float = F>, P: Default + Clone + Serialize + DeserializeOwned + ArgminScaledSub<O::Param, O::Float, O::Param> + ArgminSub<O::Param, O::Param> + ArgminAdd<O::Param, O::Param> + ArgminMul<O::Float, O::Param>, F: ArgminFloat + Sum<O::Float>,

Source§

const NAME: &'static str = "Nelder-Mead method"

Source§

impl<O, P, F> Solver<O> for ParticleSwarm<P, F>
where O: ArgminOp<Output = F, Param = P, Float = F>, O::Param: Position<F> + DeserializeOwned + Serialize, O::Hessian: Clone + Default, F: ArgminFloat,

Source§

const NAME: &'static str = "Particle Swarm Optimization"

Source§

impl<O, P, L, B, F> Solver<O> for NonlinearConjugateGradient<P, L, B, F>
where O: ArgminOp<Param = P, Output = F, Float = F>, P: Clone + Default + Serialize + DeserializeOwned + ArgminSub<O::Param, O::Param> + ArgminDot<O::Param, O::Float> + ArgminScaledAdd<O::Param, O::Float, O::Param> + ArgminAdd<O::Param, O::Param> + ArgminMul<F, O::Param> + ArgminNorm<O::Float>, O::Hessian: Default, L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<O>>, B: ArgminNLCGBetaUpdate<O::Param, O::Float>, F: ArgminFloat,

Source§

const NAME: &'static str = "Nonlinear Conjugate Gradient"

Source§

impl<O, P, L, F> Solver<O> for BacktrackingLineSearch<P, L, F>
where P: Clone + Default + Serialize + DeserializeOwned + ArgminSub<P, P> + ArgminDot<P, F> + ArgminScaledAdd<P, F, P>, O: ArgminOp<Param = P, Output = F, Float = F>, L: LineSearchCondition<P, F>, F: ArgminFloat,

Source§

const NAME: &'static str = "Backtracking Line search"

Source§

impl<O, R, F> Solver<O> for TrustRegion<R, F>
where O: ArgminOp<Output = F, Float = F>, O::Param: Default + Clone + Debug + Serialize + ArgminMul<F, O::Param> + ArgminWeightedDot<O::Param, F, O::Hessian> + ArgminNorm<F> + ArgminDot<O::Param, F> + ArgminAdd<O::Param, O::Param> + ArgminSub<O::Param, O::Param> + ArgminZeroLike, O::Hessian: Default + Clone + Debug + Serialize + ArgminDot<O::Param, O::Param>, R: ArgminTrustRegion<F> + Solver<OpWrapper<O>>, F: ArgminFloat,

Source§

const NAME: &'static str = "Trust region"

Source§

impl<P, O, F> Solver<O> for HagerZhangLineSearch<P, F>
where O: ArgminOp<Param = P, Output = F, Float = F>, P: Clone + Default + Serialize + DeserializeOwned + ArgminSub<P, P> + ArgminDot<P, F> + ArgminScaledAdd<P, F, P>, F: ArgminFloat,

Source§

const NAME: &'static str = "Hager-Zhang Line search"

Source§

impl<P, O, F> Solver<O> for MoreThuenteLineSearch<P, F>
where O: ArgminOp<Param = P, Output = F, Float = F>, P: Clone + Serialize + DeserializeOwned + ArgminSub<P, P> + ArgminDot<P, F> + ArgminScaledAdd<P, F, P>, F: ArgminFloat,

Source§

const NAME: &'static str = "More-Thuente Line search"

Source§

impl<P, O, F> Solver<O> for Steihaug<P, F>
where O: ArgminOp<Param = P, Output = F, Float = F>, P: Clone + Serialize + DeserializeOwned + Default + ArgminMul<F, P> + ArgminWeightedDot<P, F, O::Hessian> + ArgminNorm<F> + ArgminDot<P, F> + ArgminAdd<P, P> + ArgminSub<P, P> + ArgminZeroLike, O::Hessian: ArgminDot<P, P>, F: ArgminFloat,

Source§

const NAME: &'static str = "Steihaug"

Source§

impl<P, O, S, F> Solver<O> for ConjugateGradient<P, S>
where O: ArgminOp<Param = P, Output = P, Float = F>, P: Clone + Serialize + DeserializeOwned + ArgminDot<O::Param, S> + ArgminSub<O::Param, O::Param> + ArgminScaledAdd<O::Param, S, O::Param> + ArgminAdd<O::Param, O::Param> + ArgminConj + ArgminMul<O::Float, O::Param>, S: Debug + ArgminDiv<S, S> + ArgminNorm<O::Float> + ArgminConj, F: ArgminFloat,

Source§

const NAME: &'static str = "Conjugate Gradient"