gtsam/gtsam_unstable/linear/iterative.h

151 lines
3.9 KiB
C
Raw Normal View History

/* ----------------------------------------------------------------------------
2012-06-09 00:45:16 +08:00
* GTSAM Copyright 2010, Georgia Tech Research Corporation,
* Atlanta, Georgia 30332-0415
* All Rights Reserved
* Authors: Frank Dellaert, et al. (see THANKS for the full author list)
* See LICENSE for the license information
* -------------------------------------------------------------------------- */
/**
* @file iterative.h
2009-12-28 17:56:58 +08:00
* @brief Iterative methods, implementation
* @author Frank Dellaert
* @date Dec 28, 2009
2009-12-28 17:56:58 +08:00
*/
2010-01-11 06:15:34 +08:00
#pragma once
#include <gtsam/base/Matrix.h>
#include <gtsam/linear/VectorValues.h>
2012-06-09 00:45:16 +08:00
#include <gtsam/linear/ConjugateGradientSolver.h>
2009-12-31 20:55:16 +08:00
2009-12-28 17:56:58 +08:00
namespace gtsam {
2009-12-31 20:55:16 +08:00
/**
* Method of conjugate gradients (CG) template
* "System" class S needs gradient(S,v), e=S*v, v=S^e
* "Vector" class V needs dot(v,v), -v, v+v, s*v
* "Vector" class E needs dot(v,v)
* @param Ab, the "system" that needs to be solved, examples below
* @param x is the initial estimate
* @param epsilon determines the convergence criterion: norm(g)<epsilon*norm(g0)
* @param maxIterations, if 0 will be set to |x|
* @param steepest flag, if true does steepest descent, not CG
* */
template<class S, class V, class E>
V conjugateGradients(const S& Ab, V x, bool verbose, double epsilon,
size_t maxIterations, bool steepest = false);
2009-12-28 17:56:58 +08:00
2009-12-31 20:55:16 +08:00
/**
* Helper class encapsulating the combined system |Ax-b_|^2
* Needed to run Conjugate Gradients on matrices
* */
class System {
2012-06-09 10:42:45 +08:00
private:
const Matrix& A_;
const Vector& b_;
public:
System(const Matrix& A, const Vector& b) :
A_(A), b_(b) {
}
/** Access A matrix */
const Matrix& A() const { return A_; }
/** Access b vector */
const Vector& b() const { return b_; }
2010-01-31 11:33:53 +08:00
/** Apply operator A'*e */
Vector operator^(const Vector& e) const {
return A_ ^ e;
}
/**
* Print with optional string
*/
void print (const std::string& s = "System") const;
};
2009-12-28 17:56:58 +08:00
/** gradient of objective function 0.5*|Ax-b_|^2 at x = A_'*(Ax-b_) */
inline Vector gradient(const System& system, const Vector& x) {
return system.A() ^ (system.A() * x - system.b());
}
/** Apply operator A */
inline Vector operator*(const System& system, const Vector& x) {
return system.A() * x;
}
/** Apply operator A in place */
inline void multiplyInPlace(const System& system, const Vector& x, Vector& e) {
e = system.A() * x;
}
/** x += alpha* A'*e */
inline void transposeMultiplyAdd(const System& system, double alpha, const Vector& e, Vector& x) {
transposeMultiplyAdd(alpha,system.A(),e,x);
}
/**
* Method of steepest gradients, System version
*/
2010-10-28 11:26:03 +08:00
Vector steepestDescent(
const System& Ab,
const Vector& x,
const IterativeOptimizationParameters & parameters);
/**
2009-12-31 20:55:16 +08:00
* Method of conjugate gradients (CG), System version
*/
2010-10-28 11:26:03 +08:00
Vector conjugateGradientDescent(
const System& Ab,
const Vector& x,
2012-06-09 00:45:16 +08:00
const ConjugateGradientParameters & parameters);
2009-12-31 20:55:16 +08:00
/** convenience calls using matrices, will create System class internally: */
2009-12-28 17:56:58 +08:00
/**
2009-12-31 20:55:16 +08:00
* Method of steepest gradients, Matrix version
2009-12-28 17:56:58 +08:00
*/
2010-10-28 11:26:03 +08:00
Vector steepestDescent(
const Matrix& A,
const Vector& b,
const Vector& x,
2012-06-09 00:45:16 +08:00
const ConjugateGradientParameters & parameters);
2009-12-28 17:56:58 +08:00
/**
* Method of conjugate gradients (CG), Matrix version
*/
2010-10-28 11:26:03 +08:00
Vector conjugateGradientDescent(
const Matrix& A,
const Vector& b,
const Vector& x,
2012-06-09 00:45:16 +08:00
const ConjugateGradientParameters & parameters);
2009-12-28 17:56:58 +08:00
2009-12-31 20:55:16 +08:00
class GaussianFactorGraph;
/**
* Method of steepest gradients, Gaussian Factor Graph version
* */
2010-10-28 11:26:03 +08:00
VectorValues steepestDescent(
const GaussianFactorGraph& fg,
const VectorValues& x,
2012-06-09 00:45:16 +08:00
const ConjugateGradientParameters & parameters);
2009-12-31 20:55:16 +08:00
2009-12-28 17:56:58 +08:00
/**
* Method of conjugate gradients (CG), Gaussian Factor Graph version
* */
2010-10-28 11:26:03 +08:00
VectorValues conjugateGradientDescent(
const GaussianFactorGraph& fg,
const VectorValues& x,
2012-06-09 00:45:16 +08:00
const ConjugateGradientParameters & parameters);
2010-10-28 11:26:03 +08:00
2009-12-28 17:56:58 +08:00
} // namespace gtsam