1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157
use crate::{
algebra::{
abstr::Real,
linear::{matrix::General, vector::vector::Vector},
},
optimization::{Optim, OptimResult},
};
#[cfg(feature = "serde")]
use serde::{Deserialize, Serialize};
use std::clone::Clone;
/// Conjugate Gradient method
///
/// The conjugate gradient method is a solver for systems of linear equations
/// with a symmetric and positive-definite matrix. Ax = b where A is a symmetric
/// and positive-definite matrix
///
/// input: $A \in \mathbb{R}^{n \times n}$ and $ b \in \mathbb{R}^{n} $ and
/// initial approximation $x_{0} \in \mathbb{R}^{n} $
///
/// output: $ x_k $
///
/// 1. $ d_{0} = r_{0} := b - Ax_{0} $ and set $ k := 0 $
/// 2. $ \alpha_{k} := \frac{\lvert \lvert r_{k} \rvert
/// \rvert_{2}^{2}}{d_{k}^{T}Ad_{k}} $ <br> $ x_{k+1} := x_{k} +
/// \alpha_{j}d_{k} $ <br> $ r_{k+1} := r_{k} - \alpha_{k}Ad_{k} $ <br>
/// $ \beta_{k} := \frac{\lvert \lvert r_{k+1} \rvert
/// \rvert_{2}^{2}}{\lvert \lvert r_{k} \rvert \rvert_{2}^{2}} $ <br>
/// $ d_{k+1} := r_{k+1} + \beta_{k}d_{k} $ <br>
/// 3. if $ \lvert \lvert r_{k+ 1} \rvert \rvert_{2} > \epsilon $ increase
/// $k:= k + 1$ and goto 2.
///
/// # Example
///
/// ```
/// use mathru::*;
/// use mathru::algebra::linear::{vector::Vector, matrix::General};
/// use mathru::optimization::{Optim, ConjugateGradient};
///
/// struct LinearEquation
/// { ///
/// a: General<f64>,
/// b: Vector<f64>,
/// }
///
/// //Ax = b
/// impl LinearEquation
/// {
/// pub fn new() -> LinearEquation
/// {
/// LinearEquation
/// {
/// a: matrix![1.0, 3.0; 3.0, 5.0],
/// b: vector![-7.0; 7.0]
/// }
/// }
/// }
///
/// impl Optim<f64> for LinearEquation
/// {
/// // A
/// fn jacobian(&self, _input: &Vector<f64>) -> General<f64>
/// {
/// return self.a.clone();
/// }
///
/// // f = b-Ax
/// fn eval(&self, x: &Vector<f64>) -> Vector<f64>
/// {
/// return self.b.clone() - self.a.clone() * x.clone()
/// }
///
/// //Computes the Hessian at the given value x
/// fn hessian(&self, _x: &Vector<f64>) -> General<f64>
/// {
/// unimplemented!();
/// }
/// }
///
/// //create optimizer instance
/// let optim: ConjugateGradient<f64> = ConjugateGradient::new(10, 0.01);
///
/// let leq: LinearEquation = LinearEquation::new();
///
/// // Initial approximation
/// let x_0: Vector<f64> = vector![1.0; 1.0];
///
/// // Minimize function
/// let x_min: Vector<f64> = optim.minimize(&leq, &x_0).arg();
/// ```
#[cfg_attr(feature = "serde", derive(Serialize, Deserialize))]
#[derive(Clone, Copy, Debug)]
pub struct ConjugateGradient<T> {
iters: u64,
epsilon: T,
}
impl<T> ConjugateGradient<T> {
/// Creates an instance of ConjugateGradient method
///
/// # Arguments
///
/// * 'iters': Number of iterations
/// * 'epsilon':
pub fn new(iters: u64, epsilon: T) -> ConjugateGradient<T> {
ConjugateGradient { iters, epsilon }
}
}
impl<T> ConjugateGradient<T>
where
T: Real,
{
/// Minimize function func
///
/// # Arguments
///
/// * 'func0': Function to be minimized
/// * 'x_0': Initial guess for the minimum
///
/// # Return
///
/// local minimum
pub fn minimize<F>(&self, func: &F, x_0: &Vector<T>) -> OptimResult<Vector<T>>
where
F: Optim<T>,
{
let mut x_n: Vector<T> = x_0.clone();
let mut d_n: Vector<T> = func.eval(&x_n);
let mut r_n: Vector<T> = d_n.clone();
for _i in 0..self.iters {
let jacobian: General<T> = func.jacobian(&x_n);
let temp: Vector<T> = d_n.clone().transpose() * jacobian.clone();
let alpha_n: T = r_n.dotp(&r_n) / temp.transpose().dotp(&d_n);
x_n += d_n.clone() * alpha_n;
let r_n_1: Vector<T> = r_n.clone() - jacobian * d_n.clone() * alpha_n;
let beta_n: T = r_n_1.dotp(&r_n_1) / r_n.dotp(&r_n);
d_n = r_n_1.clone() + d_n * beta_n;
if r_n_1.dotp(&r_n_1).pow(T::from_f64(0.5)) < self.epsilon {
break;
}
r_n = r_n_1
}
OptimResult::new(x_n)
}
}