testPCGSolver.cpp
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2 
3  * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4  * Atlanta, Georgia 30332-0415
5  * All Rights Reserved
6  * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7 
8  * See LICENSE for the license information
9 
10  * -------------------------------------------------------------------------- */
11 
19 #include <tests/smallExample.h>
22 #include <gtsam/linear/PCGSolver.h>
24 #include <gtsam/inference/Symbol.h>
25 #include <gtsam/base/Matrix.h>
26 
28 
29 #include <iostream>
30 #include <fstream>
31 
32 using namespace std;
33 using namespace gtsam;
34 
35 const double tol = 1e-3;
36 
39 
40 /* ************************************************************************* */
41 // Test cholesky decomposition
43  Matrix R = (Matrix(3,3) <<
44  1., -1., -1.,
45  0., 2., -1.,
46  0., 0., 1.).finished();
47  Matrix AtA = R.transpose() * R;
48 
49  Vector Rvector = (Vector(9) << 1., -1., -1.,
50  0., 2., -1.,
51  0., 0., 1.).finished();
52 // Vector Rvector = (Vector(6) << 1., -1., -1.,
53 // 2., -1.,
54 // 1.).finished();
55 
56  Vector b = Vector3(1., 2., 3.);
57 
58  Vector x = Vector3(6.5, 2.5, 3.) ;
59 
60  /* test cholesky */
61  Matrix Rhat = AtA.llt().matrixL().transpose();
62  EXPECT(assert_equal(R, Rhat, 1e-5));
63 
64  /* test backward substitution */
65  Vector xhat = Rhat.triangularView<Eigen::Upper>().solve(b);
66  EXPECT(assert_equal(x, xhat, 1e-5));
67 
68  /* test in-place back substitution */
69  xhat = b;
70  Rhat.triangularView<Eigen::Upper>().solveInPlace(xhat);
71  EXPECT(assert_equal(x, xhat, 1e-5));
72 
73  /* test triangular matrix map */
74  Eigen::Map<Eigen::MatrixXd> Radapter(Rvector.data(), 3, 3);
75  xhat = Radapter.transpose().triangularView<Eigen::Upper>().solve(b);
76  EXPECT(assert_equal(x, xhat, 1e-5));
77 
78 }
79 
80 /* ************************************************************************* */
81 // Test GaussianFactorGraphSystem::multiply and getb
83 {
84  // Create a Gaussian Factor Graph
85  GaussianFactorGraph simpleGFG;
86  SharedDiagonal unit2 = noiseModel::Diagonal::Sigmas(Vector2(0.5, 0.3));
87  simpleGFG.emplace_shared<JacobianFactor>(2, (Matrix(2,2)<< 10, 0, 0, 10).finished(), (Vector(2) << -1, -1).finished(), unit2);
88  simpleGFG.emplace_shared<JacobianFactor>(2, (Matrix(2,2)<< -10, 0, 0, -10).finished(), 0, (Matrix(2,2)<< 10, 0, 0, 10).finished(), (Vector(2) << 2, -1).finished(), unit2);
89  simpleGFG.emplace_shared<JacobianFactor>(2, (Matrix(2,2)<< -5, 0, 0, -5).finished(), 1, (Matrix(2,2)<< 5, 0, 0, 5).finished(), (Vector(2) << 0, 1).finished(), unit2);
90  simpleGFG.emplace_shared<JacobianFactor>(0, (Matrix(2,2)<< -5, 0, 0, -5).finished(), 1, (Matrix(2,2)<< 5, 0, 0, 5).finished(), (Vector(2) << -1, 1.5).finished(), unit2);
91  simpleGFG.emplace_shared<JacobianFactor>(0, (Matrix(2,2)<< 1, 0, 0, 1).finished(), (Vector(2) << 0, 0).finished(), unit2);
92  simpleGFG.emplace_shared<JacobianFactor>(1, (Matrix(2,2)<< 1, 0, 0, 1).finished(), (Vector(2) << 0, 0).finished(), unit2);
93  simpleGFG.emplace_shared<JacobianFactor>(2, (Matrix(2,2)<< 1, 0, 0, 1).finished(), (Vector(2) << 0, 0).finished(), unit2);
94 
95  // Create a dummy-preconditioner and a GaussianFactorGraphSystem
96  DummyPreconditioner dummyPreconditioner;
97  KeyInfo keyInfo(simpleGFG);
98  std::map<Key,Vector> lambda;
99  dummyPreconditioner.build(simpleGFG, keyInfo, lambda);
100  GaussianFactorGraphSystem gfgs(simpleGFG, dummyPreconditioner, keyInfo, lambda);
101 
102  // Prepare container for each variable
103  Vector initial, residual, preconditionedResidual, p, actualAp;
104  initial = (Vector(6) << 0., 0., 0., 0., 0., 0.).finished();
105 
106  // Calculate values using GaussianFactorGraphSystem same as inside of PCGSolver
107  gfgs.residual(initial, residual); /* r = b-Ax */
108  gfgs.leftPrecondition(residual, preconditionedResidual); /* pr = L^{-1} (b-Ax) */
109  gfgs.rightPrecondition(preconditionedResidual, p); /* p = L^{-T} pr */
110  gfgs.multiply(p, actualAp); /* A p */
111 
112  // Expected value of Ap for the first iteration of this example problem
113  Vector expectedAp = (Vector(6) << 100400, -249074.074, -2080, 148148.148, -146480, 37962.963).finished();
114  EXPECT(assert_equal(expectedAp, actualAp, 1e-3));
115 
116  // Expected value of getb
117  Vector expectedb = (Vector(6) << 100.0, -194.444, -20.0, 138.889, -120.0, -55.556).finished();
118  Vector actualb;
119  gfgs.getb(actualb);
120  EXPECT(assert_equal(expectedb, actualb, 1e-3));
121 }
122 
123 /* ************************************************************************* */
124 // Test Dummy Preconditioner
127  params.linearSolverType = LevenbergMarquardtParams::Iterative;
128  auto pcg = std::make_shared<PCGSolverParameters>(
129  std::make_shared<DummyPreconditionerParameters>());
130  params.iterativeParams = pcg;
131 
133 
134  Point2 x0(10, 10);
135  Values c0;
136  c0.insert(X(1), x0);
137 
138  Values actualPCG = LevenbergMarquardtOptimizer(fg, c0, params).optimize();
139 
140  DOUBLES_EQUAL(0, fg.error(actualPCG), tol);
141 }
142 
143 /* ************************************************************************* */
144 // Test Block-Jacobi Precondioner
145 TEST(PCGSolver, blockjacobi) {
147  params.linearSolverType = LevenbergMarquardtParams::Iterative;
148  auto pcg = std::make_shared<PCGSolverParameters>(
149  std::make_shared<BlockJacobiPreconditionerParameters>());
150  params.iterativeParams = pcg;
151 
153 
154  Point2 x0(10, 10);
155  Values c0;
156  c0.insert(X(1), x0);
157 
158  Values actualPCG = LevenbergMarquardtOptimizer(fg, c0, params).optimize();
159 
160  DOUBLES_EQUAL(0, fg.error(actualPCG), tol);
161 }
162 
163 /* ************************************************************************* */
164 // Test Incremental Subgraph PCG Solver
165 TEST(PCGSolver, subgraph) {
167  params.linearSolverType = LevenbergMarquardtParams::Iterative;
168  auto pcg = std::make_shared<PCGSolverParameters>(
169  std::make_shared<SubgraphPreconditionerParameters>());
170  params.iterativeParams = pcg;
171 
173 
174  Point2 x0(10, 10);
175  Values c0;
176  c0.insert(X(1), x0);
177 
178  Values actualPCG = LevenbergMarquardtOptimizer(fg, c0, params).optimize();
179 
180  DOUBLES_EQUAL(0, fg.error(actualPCG), tol);
181 }
182 
183 /* ************************************************************************* */
184 int main() {
185  TestResult tr;
186  return TestRegistry::runAllTests(tr);
187 }
188 
TestRegistry::runAllTests
static int runAllTests(TestResult &result)
Definition: TestRegistry.cpp:27
main
int main()
Definition: testPCGSolver.cpp:184
gtsam::GaussianFactorGraphSystem::rightPrecondition
void rightPrecondition(const Vector &x, Vector &y) const
Definition: PCGSolver.cpp:126
gtsam::NonlinearOptimizer::optimize
virtual const Values & optimize()
Definition: NonlinearOptimizer.h:98
GaussianFactorGraph.h
Linear Factor Graph where all factors are Gaussians.
gtsam::GaussianFactorGraphSystem::residual
void residual(const Vector &x, Vector &r) const
Definition: PCGSolver.cpp:77
llt
EIGEN_DONT_INLINE void llt(const Mat &A, const Mat &B, Mat &C)
Definition: llt.cpp:5
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
gtsam::GaussianFactorGraphSystem::getb
void getb(Vector &b) const
Definition: PCGSolver.cpp:107
EXPECT
#define EXPECT(condition)
Definition: Test.h:150
gtsam::Vector2
Eigen::Vector2d Vector2
Definition: Vector.h:42
TestHarness.h
b
Scalar * b
Definition: benchVecAdd.cpp:17
LevenbergMarquardtOptimizer.h
A nonlinear optimizer that uses the Levenberg-Marquardt trust-region scheme.
Matrix.h
typedef and functions to augment Eigen's MatrixXd
initial
Values initial
Definition: OdometryOptimize.cpp:2
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
gtsam::JacobianFactor
Definition: JacobianFactor.h:91
Eigen::Upper
@ Upper
Definition: Constants.h:211
gtsam::Matrix
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
X
#define X
Definition: icosphere.cpp:20
gtsam::GaussianFactorGraphSystem::multiply
void multiply(const Vector &x, Vector &y) const
Definition: PCGSolver.cpp:90
gtsam::Vector3
Eigen::Vector3d Vector3
Definition: Vector.h:43
vanilla::params
static const SmartProjectionParams params
Definition: smartFactorScenarios.h:69
gtsam::Vector
Eigen::VectorXd Vector
Definition: Vector.h:38
pybind_wrapper_test_script.dummy
dummy
Definition: pybind_wrapper_test_script.py:42
gtsam::KeyInfo
Definition: IterativeSolver.h:127
gtsam::GaussianFactorGraph
Definition: GaussianFactorGraph.h:73
gtsam::NonlinearFactorGraph::error
double error(const Values &values) const
Definition: NonlinearFactorGraph.cpp:170
ceres::Matrix
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > Matrix
Definition: gtsam/3rdparty/ceres/eigen.h:42
gtsam::PCGSolver
Definition: PCGSolver.h:55
DOUBLES_EQUAL
#define DOUBLES_EQUAL(expected, actual, threshold)
Definition: Test.h:141
gtsam::NonlinearFactorGraph
Definition: NonlinearFactorGraph.h:55
gtsam::SharedDiagonal
noiseModel::Diagonal::shared_ptr SharedDiagonal
Definition: NoiseModel.h:764
PCGSolver.h
x0
static Symbol x0('x', 0)
gtsam::LevenbergMarquardtOptimizer
Definition: LevenbergMarquardtOptimizer.h:35
L
MatrixXd L
Definition: LLT_example.cpp:6
Symbol.h
gtsam::Point2
Vector2 Point2
Definition: Point2.h:32
gtsam::Rot2::transpose
Matrix2 transpose() const
Definition: Rot2.cpp:92
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
tol
const double tol
Definition: testPCGSolver.cpp:35
lambda
static double lambda[]
Definition: jv.c:524
TestResult
Definition: TestResult.h:26
gtsam::DummyPreconditioner
Definition: Preconditioner.h:103
gtsam::GaussianFactorGraphSystem::leftPrecondition
void leftPrecondition(const Vector &x, Vector &y) const
Definition: PCGSolver.cpp:118
gtsam
traits
Definition: SFMdata.h:40
gtsam::GaussianFactorGraphSystem
Definition: PCGSolver.h:82
gtsam::Values
Definition: Values.h:65
std
Definition: BFloat16.h:88
SubgraphPreconditioner.h
p
float * p
Definition: Tutorial_Map_using.cpp:9
gtsam::Values::insert
void insert(Key j, const Value &val)
Definition: Values.cpp:155
gtsam::assert_equal
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:40
gtsam::LevenbergMarquardtParams
Definition: LevenbergMarquardtParams.h:35
gtsam::example::createReallyNonlinearFactorGraph
NonlinearFactorGraph createReallyNonlinearFactorGraph()
Definition: smallExample.h:378
initial
Definition: testScenarioRunner.cpp:148
TEST
TEST(PCGSolver, llt)
Definition: testPCGSolver.cpp:42
smallExample.h
Create small example with two poses and one landmark.
ceres::Vector
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
Definition: gtsam/3rdparty/ceres/eigen.h:38
R
Rot2 R(Rot2::fromAngle(0.1))
gtsam::FactorGraph::emplace_shared
IsDerived< DERIVEDFACTOR > emplace_shared(Args &&... args)
Emplace a shared pointer to factor of given type.
Definition: FactorGraph.h:153


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:08:06