Struct argmin::solver::conjugategradient::NonlinearConjugateGradient
source · [−]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
sourceimpl<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
Specifiy 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
sourceimpl<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>
sourcefn clone(&self) -> NonlinearConjugateGradient<P, L, B, F>
fn clone(&self) -> NonlinearConjugateGradient<P, L, B, F>
Returns a copy of the value. Read more
1.0.0 · sourcefn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
Performs copy-assignment from source
. Read more
sourceimpl<'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>,
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>,
sourcefn 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>,
Deserialize this value from the given Serde deserializer. Read more
sourceimpl<P, L, B, F> Serialize for NonlinearConjugateGradient<P, L, B, F> where
P: Serialize,
L: Serialize,
B: Serialize,
F: Serialize,
impl<P, L, B, F> Serialize for NonlinearConjugateGradient<P, L, B, F> where
P: Serialize,
L: Serialize,
B: Serialize,
F: Serialize,
sourceimpl<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 + SerializeAlias + DeserializeOwnedAlias + ArgminAdd<P, P> + ArgminMul<F, P>,
G: Clone + SerializeAlias + DeserializeOwnedAlias + ArgminMul<F, P> + ArgminDot<G, F> + ArgminNorm<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 + SerializeAlias + DeserializeOwnedAlias + ArgminAdd<P, P> + ArgminMul<F, P>,
G: Clone + SerializeAlias + DeserializeOwnedAlias + ArgminMul<F, P> + ArgminDot<G, F> + ArgminNorm<F>,
L: Clone + LineSearch<P, F> + Solver<O, IterState<P, G, (), (), F>>,
B: NLCGBetaUpdate<G, P, F>,
F: ArgminFloat,
sourceconst NAME: &'static str = "Nonlinear Conjugate Gradient"
const NAME: &'static str = "Nonlinear Conjugate Gradient"
Name of the solver. Mainly used in Observers.
sourcefn 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>
Initializes the algorithm. Read more
sourcefn 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>
sourcefn terminate_internal(&mut self, state: &I) -> TerminationReason
fn terminate_internal(&mut self, state: &I) -> TerminationReason
Checks whether basic termination reasons apply. Read more
sourcefn terminate(&mut self, _state: &I) -> TerminationReason
fn terminate(&mut self, _state: &I) -> TerminationReason
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> where
B: RefUnwindSafe,
F: RefUnwindSafe,
L: RefUnwindSafe,
P: RefUnwindSafe,
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> where
B: UnwindSafe,
F: UnwindSafe,
L: UnwindSafe,
P: UnwindSafe,
Blanket Implementations
sourceimpl<T> BorrowMut<T> for T where
T: ?Sized,
impl<T> BorrowMut<T> for T where
T: ?Sized,
const: unstable · sourcefn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more