Trait argmin::core::Solver[][src]

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

    fn next_iter(
        &mut self,
        op: &mut OpWrapper<O>,
        state: &IterState<O>
    ) -> Result<ArgminIterData<O>, Error>; 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 { ... } }

Solver

Every solver needs to implement this trait.

Associated Constants

const NAME: &'static str[src]

Name of the solver

Loading content...

Required methods

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

Computes one iteration of the algorithm.

Loading content...

Provided methods

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

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.

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

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.

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

Checks whether the algorithm must be terminated

Loading content...

Implementors

impl<O, B, R, F> Solver<O> for SR1TrustRegion<B, R, F> where
    O: ArgminOp<Output = F, Hessian = B, Float = F>,
    O::Param: Debug + Clone + Default + Serialize + ArgminSub<O::Param, O::Param> + ArgminAdd<O::Param, O::Param> + ArgminDot<O::Param, O::Float> + ArgminDot<O::Param, O::Hessian> + ArgminNorm<O::Float> + ArgminZeroLike + ArgminMul<F, O::Param>,
    B: 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>,
    R: ArgminTrustRegion<F> + Solver<OpWrapper<O>>,
    F: ArgminFloat + ArgminNorm<O::Float>, 
[src]

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

impl<O, F> Solver<O> for GaussNewton<F> where
    O: ArgminOp<Float = F>,
    O::Param: 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<O::Jacobian> + ArgminInv<O::Jacobian> + ArgminDot<O::Jacobian, O::Jacobian> + ArgminDot<O::Output, O::Param> + ArgminDot<O::Param, O::Param>,
    F: ArgminFloat
[src]

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

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

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
[src]

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

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

Perform one iteration of SA algorithm

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
[src]

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
[src]

impl<O, L, F> Solver<O> for GaussNewtonLS<L, F> where
    O: ArgminOp<Float = F>,
    O::Param: 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<O::Jacobian> + ArgminInv<O::Jacobian> + ArgminDot<O::Jacobian, O::Jacobian> + ArgminDot<O::Output, O::Param> + ArgminDot<O::Param, O::Param>,
    L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<LineSearchOP<O>>>,
    F: ArgminFloat
[src]

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
[src]

impl<O, L, F> Solver<O> for NewtonCG<L, F> where
    O: ArgminOp<Output = F, Float = F>,
    O::Param: Send + Sync + Clone + Serialize + Default + ArgminSub<O::Param, O::Param> + ArgminAdd<O::Param, O::Param> + ArgminDot<O::Param, O::Float> + ArgminScaledAdd<O::Param, O::Float, O::Param> + ArgminMul<F, O::Param> + ArgminConj + ArgminZeroLike + ArgminNorm<O::Float>,
    O::Hessian: Send + Sync + Default + Clone + Serialize + ArgminInv<O::Hessian> + ArgminDot<O::Param, O::Param>,
    L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<O>>,
    F: ArgminFloat + Default + ArgminDiv<O::Float, O::Float> + ArgminNorm<O::Float> + ArgminConj
[src]

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

impl<O, L, H, F> Solver<O> for DFP<L, H, F> where
    O: ArgminOp<Output = F, Hessian = H, Float = F>,
    O::Param: Clone + Default + Serialize + ArgminSub<O::Param, O::Param> + ArgminDot<O::Param, O::Float> + ArgminDot<O::Param, O::Hessian> + ArgminScaledAdd<O::Param, F, O::Param> + ArgminNorm<O::Float> + ArgminMul<O::Float, O::Param> + ArgminTranspose<O::Param>,
    H: 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> + ArgminTranspose<O::Hessian> + ArgminEye,
    L: Clone + ArgminLineSearch<O::Param, O::Float> + Solver<OpWrapper<O>>,
    F: ArgminFloat
[src]

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>,
    H: 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
[src]

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
[src]

impl<O, P, F> Solver<O> for NelderMead<P, F> where
    O: ArgminOp<Output = F, Param = P, Float = F>,
    P: 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>, 
[src]

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

fn next_iter(
    &mut self,
    _op: &mut OpWrapper<O>,
    _state: &IterState<O>
) -> Result<ArgminIterData<O>, Error>
[src]

Perform one iteration of algorithm

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
[src]

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
[src]

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
[src]

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
[src]

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
[src]

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
[src]

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
[src]

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

Perform one iteration of CG algorithm

Loading content...