27class GaussianFactorGraph;
31struct PreconditionerParameters;
39 typedef boost::shared_ptr<PCGSolverParameters> shared_ptr;
44 void print(std::ostream &os)
const override;
48 return *preconditioner_;
52 void print(
const std::string &s)
const;
54 boost::shared_ptr<PreconditionerParameters> preconditioner_;
56 void setPreconditionerParams(
const boost::shared_ptr<PreconditionerParameters> preconditioner);
65 typedef boost::shared_ptr<PCGSolver> shared_ptr;
70 boost::shared_ptr<Preconditioner> preconditioner_;
78 using IterativeSolver::optimize;
81 const KeyInfo &keyInfo,
const std::map<Key, Vector> &lambda,
94 const std::map<Key, Vector> &lambda);
99 const std::map<Key, Vector> &lambda_;
101 void residual(
const Vector &x, Vector &r)
const;
102 void multiply(
const Vector &x, Vector& y)
const;
103 void leftPrecondition(
const Vector &x, Vector &y)
const;
104 void rightPrecondition(
const Vector &x, Vector &y)
const;
105 inline void scal(
const double alpha, Vector &x)
const {
108 inline double dot(
const Vector &x,
const Vector &y)
const {
111 inline void axpy(
const double alpha,
const Vector &x, Vector &y)
const {
115 void getb(Vector &b)
const;
123 const std::map<Key, size_t> & dimensions);
Implementation of Conjugate Gradient solver for a linear system.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
VectorValues buildVectorValues(const Vector &v, const Ordering &ordering, const map< Key, size_t > &dimensions)
Create VectorValues from a Vector.
Definition: PCGSolver.cpp:140
Point3 optimize(const NonlinearFactorGraph &graph, const Values &values, Key landmarkKey)
Optimize for triangulation.
Definition: triangulation.cpp:73
void print(const Matrix &A, const string &s, ostream &stream)
print without optional string, must specify cout yourself
Definition: Matrix.cpp:155
void GTSAM_DEPRECATED scal(double alpha, Vector &x)
BLAS Level 1 scal: x <- alpha*x.
Definition: Vector.h:210
double dot(const V1 &a, const V2 &b)
Dot product.
Definition: Vector.h:194
void GTSAM_DEPRECATED axpy(double alpha, const V1 &x, V2 &y)
BLAS Level 1 axpy: y <- alpha*x + y.
Definition: Vector.h:217
Definition: Ordering.h:34
parameters for the conjugate gradient method
Definition: ConjugateGradientSolver.h:29
A Linear Factor Graph is a factor graph where all factors are Gaussian, i.e.
Definition: GaussianFactorGraph.h:69
parameters for iterative linear solvers
Definition: IterativeSolver.h:44
Base class for Iterative Solvers like SubgraphSolver.
Definition: IterativeSolver.h:86
Handy data structure for iterative solvers.
Definition: IterativeSolver.h:126
Parameters for PCG.
Definition: PCGSolver.h:36
A virtual base class for the preconditioned conjugate gradient solver.
Definition: PCGSolver.h:62
System class needed for calling preconditionedConjugateGradient.
Definition: PCGSolver.h:89
Definition: Preconditioner.h:24
Definition: Preconditioner.h:64
This class represents a collection of vector-valued variables associated each with a unique integer i...
Definition: VectorValues.h:74