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,
impl<P, L, B, F> NonlinearConjugateGradient<P, L, B, F>where
F: ArgminFloat,
Sourcepub fn new(linesearch: L, beta_method: B) -> Self
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);Sourcepub fn restart_iters(self, iters: u64) -> Self
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);Sourcepub fn restart_orthogonality(self, v: F) -> Self
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>
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>
fn clone(&self) -> NonlinearConjugateGradient<P, L, B, F>
1.0.0 · Source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
source. Read moreSource§impl<'de, P, L, B, F> Deserialize<'de> for NonlinearConjugateGradient<P, L, B, F>
impl<'de, P, L, B, F> Deserialize<'de> for NonlinearConjugateGradient<P, L, B, F>
Source§fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error>where
__D: Deserializer<'de>,
Source§impl<P, L, B, F> Serialize for NonlinearConjugateGradient<P, L, B, F>
impl<P, L, B, F> Serialize for NonlinearConjugateGradient<P, L, B, F>
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,
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§fn init(
&mut self,
problem: &mut Problem<O>,
state: IterState<P, G, (), (), (), F>,
) -> Result<(IterState<P, G, (), (), (), F>, Option<KV>), Error>
fn init( &mut self, problem: &mut Problem<O>, state: IterState<P, G, (), (), (), F>, ) -> Result<(IterState<P, G, (), (), (), F>, Option<KV>), Error>
Source§fn next_iter(
&mut self,
problem: &mut Problem<O>,
state: IterState<P, G, (), (), (), F>,
) -> Result<(IterState<P, G, (), (), (), F>, Option<KV>), Error>
fn next_iter( &mut self, problem: &mut Problem<O>, state: IterState<P, G, (), (), (), F>, ) -> Result<(IterState<P, G, (), (), (), F>, Option<KV>), Error>
state and optionally a KV which holds key-value pairs used in
Observers.Source§fn terminate_internal(&mut self, state: &I) -> TerminationStatus
fn terminate_internal(&mut self, state: &I) -> TerminationStatus
Source§fn terminate(&mut self, _state: &I) -> TerminationStatus
fn terminate(&mut self, _state: &I) -> TerminationStatus
terminate_internal. Read moreAuto Trait Implementations§
impl<P, L, B, F> Freeze for NonlinearConjugateGradient<P, L, B, F>
impl<P, L, B, F> RefUnwindSafe for NonlinearConjugateGradient<P, L, B, F>
impl<P, L, B, F> Send for NonlinearConjugateGradient<P, L, B, F>
impl<P, L, B, F> Sync for NonlinearConjugateGradient<P, L, B, F>
impl<P, L, B, F> Unpin for NonlinearConjugateGradient<P, L, B, F>
impl<P, L, B, F> UnwindSafe for NonlinearConjugateGradient<P, L, B, F>
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self from the equivalent element of its
superset. Read more§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self is actually part of its subset T (and can be converted to it).§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset but without any property checks. Always succeeds.§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self to the equivalent element of its superset.