pub struct NonlinearConjugateGradient<P, L, B, F> { /* private fields */ }
Expand description

§Non-linear Conjugate Gradient method

A generalization of the conjugate gradient method for nonlinear optimization problems.

Requires an initial parameter vector.

§Requirements on the optimization problem

The optimization problem is required to implement CostFunction and Gradient.

§Reference

Jorge Nocedal and Stephen J. Wright (2006). Numerical Optimization. Springer. ISBN 0-387-30303-0.

Implementations§

source§

impl<P, L, B, F> NonlinearConjugateGradient<P, L, B, F>
where F: ArgminFloat,

source

pub fn new(linesearch: L, beta_method: B) -> Self

Construct a new instance of NonlinearConjugateGradient.

Takes a LineSearch and a NLCGBetaUpdate as input.

§Example
let nlcg: NonlinearConjugateGradient<Vec<f64>, _, _, f64> =
    NonlinearConjugateGradient::new(linesearch, beta_method);
source

pub fn restart_iters(self, iters: u64) -> Self

Specify the number of iterations after which a restart should be performed.

This allows the algorithm to “forget” previous information which may not be helpful anymore.

§Example
let nlcg = nlcg.restart_iters(100);
source

pub fn restart_orthogonality(self, v: F) -> Self

Set the value for the orthogonality measure.

Setting this parameter leads to a restart of the algorithm (setting beta = 0) after consecutive search directions stop being orthogonal. In other words, if this condition is met:

|\nabla f_k^T * \nabla f_{k-1}| / | \nabla f_k |^2 >= v

A typical value for v is 0.1.

§Example
let nlcg = nlcg.restart_orthogonality(0.1);

Trait Implementations§

source§

impl<P: Clone, L: Clone, B: Clone, F: Clone> Clone for NonlinearConjugateGradient<P, L, B, F>

source§

fn clone(&self) -> NonlinearConjugateGradient<P, L, B, F>

Returns a copy of the value. Read more
1.0.0 · source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
source§

impl<'de, P, L, B, F> Deserialize<'de> for NonlinearConjugateGradient<P, L, B, F>
where P: Deserialize<'de>, L: Deserialize<'de>, B: Deserialize<'de>, F: Deserialize<'de>,

source§

fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>
where __D: Deserializer<'de>,

Deserialize this value from the given Serde deserializer. Read more
source§

impl<P, L, B, F> Serialize for NonlinearConjugateGradient<P, L, B, F>
where P: Serialize, L: Serialize, B: Serialize, F: Serialize,

source§

fn serialize<__S>(&self, __serializer: __S) -> Result<__S::Ok, __S::Error>
where __S: Serializer,

Serialize this value into the given Serde serializer. Read more
source§

impl<O, P, G, L, B, F> Solver<O, IterState<P, G, (), (), (), F>> for NonlinearConjugateGradient<P, L, B, F>
where O: CostFunction<Param = P, Output = F> + Gradient<Param = P, Gradient = G>, P: Clone + ArgminAdd<P, P> + ArgminMul<F, P>, G: Clone + ArgminMul<F, P> + ArgminDot<G, F> + ArgminL2Norm<F>, L: Clone + LineSearch<P, F> + Solver<O, IterState<P, G, (), (), (), F>>, B: NLCGBetaUpdate<G, P, F>, F: ArgminFloat,

source§

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

Name of the solver. Mainly used in Observers.
source§

fn init( &mut self, problem: &mut Problem<O>, state: IterState<P, G, (), (), (), F> ) -> Result<(IterState<P, G, (), (), (), F>, Option<KV>), Error>

Initializes the algorithm. Read more
source§

fn next_iter( &mut self, problem: &mut Problem<O>, state: IterState<P, G, (), (), (), F> ) -> Result<(IterState<P, G, (), (), (), F>, Option<KV>), Error>

Computes a single iteration of the algorithm and has access to the optimization problem definition and the internal state of the solver. Returns an updated state and optionally a KV which holds key-value pairs used in Observers.
source§

fn terminate_internal(&mut self, state: &I) -> TerminationStatus

Checks whether basic termination reasons apply. Read more
source§

fn terminate(&mut self, _state: &I) -> TerminationStatus

Used to implement stopping criteria, in particular criteria which are not covered by (terminate_internal. Read more

Auto Trait Implementations§

§

impl<P, L, B, F> RefUnwindSafe for NonlinearConjugateGradient<P, L, B, F>

§

impl<P, L, B, F> Send for NonlinearConjugateGradient<P, L, B, F>
where B: Send, F: Send, L: Send, P: Send,

§

impl<P, L, B, F> Sync for NonlinearConjugateGradient<P, L, B, F>
where B: Sync, F: Sync, L: Sync, P: Sync,

§

impl<P, L, B, F> Unpin for NonlinearConjugateGradient<P, L, B, F>
where B: Unpin, F: Unpin, L: Unpin, P: Unpin,

§

impl<P, L, B, F> UnwindSafe for NonlinearConjugateGradient<P, L, B, F>

Blanket Implementations§

source§

impl<T> Any for T
where T: 'static + ?Sized,

source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
source§

impl<T> Borrow<T> for T
where T: ?Sized,

source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
source§

impl<T> From<T> for T

source§

fn from(t: T) -> T

Returns the argument unchanged.

source§

impl<T, U> Into<U> for T
where U: From<T>,

source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

source§

impl<T> ToOwned for T
where T: Clone,

§

type Owned = T

The resulting type after obtaining ownership.
source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

§

type Error = Infallible

The type returned in the event of a conversion error.
source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
§

impl<V, T> VZip<V> for T
where V: MultiLane<T>,

§

fn vzip(self) -> V

source§

impl<T> DeserializeOwned for T
where T: for<'de> Deserialize<'de>,

source§

impl<T> SendAlias for T

source§

impl<T> SyncAlias for T