LinearizedFactor.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 <iostream>
20 
21 namespace gtsam {
22 
23 /* ************************************************************************* */
25  const GaussianFactor::shared_ptr& gaussian, const Values& lin_points)
26 : NonlinearFactor(gaussian->keys())
27 {
28  // Extract the keys and linearization points
29  for(const Key& key: gaussian->keys()) {
30  // extract linearization point
31  assert(lin_points.exists(key));
32  this->lin_points_.insert(key, lin_points.at(key));
33  }
34 }
35 
36 /* ************************************************************************* */
37 // LinearizedJacobianFactor
38 /* ************************************************************************* */
40 }
41 
42 /* ************************************************************************* */
44  const JacobianFactor::shared_ptr& jacobian, const Values& lin_points)
45 : Base(jacobian, lin_points) {
46 
47  // Create the dims array
48  size_t *dims = (size_t *)alloca(sizeof(size_t) * (jacobian->size() + 1));
49  size_t index = 0;
50  for(JacobianFactor::const_iterator iter = jacobian->begin(); iter != jacobian->end(); ++iter) {
51  dims[index++] = jacobian->getDim(iter);
52  }
53  dims[index] = 1;
54 
55  // Update the BlockInfo accessor
56  Ab_ = VerticalBlockMatrix(dims, dims+jacobian->size()+1, jacobian->rows());
57  // Get the Ab matrix from the Jacobian factor, with any covariance baked in
58  Ab_.matrix() = jacobian->augmentedJacobian();
59 }
60 
61 /* ************************************************************************* */
62 void LinearizedJacobianFactor::print(const std::string& s, const KeyFormatter& keyFormatter) const {
63 
64  std::cout << s << std::endl;
65 
66  std::cout << "Nonlinear Keys: ";
67  for(const Key& key: this->keys())
68  std::cout << keyFormatter(key) << " ";
69  std::cout << std::endl;
70 
71  for(const_iterator key=begin(); key!=end(); ++key)
72  std::cout << boost::format("A[%1%]=\n")%keyFormatter(*key) << A(*key) << std::endl;
73  std::cout << "b=\n" << b() << std::endl;
74 
75  lin_points_.print("Linearization Point: ");
76 }
77 
78 /* ************************************************************************* */
80 
81  const This *e = dynamic_cast<const This*> (&expected);
82  if (e) {
83 
84  Matrix thisMatrix = this->Ab_.range(0, Ab_.nBlocks());
85  Matrix rhsMatrix = e->Ab_.range(0, Ab_.nBlocks());
86 
87  return Base::equals(expected, tol)
88  && lin_points_.equals(e->lin_points_, tol)
89  && equal_with_abs_tol(thisMatrix, rhsMatrix, tol);
90  } else {
91  return false;
92  }
93 }
94 
95 /* ************************************************************************* */
97  Vector errorVector = error_vector(c);
98  return 0.5 * errorVector.dot(errorVector);
99 }
100 
101 /* ************************************************************************* */
102 boost::shared_ptr<GaussianFactor>
104 
105  // Create the 'terms' data structure for the Jacobian constructor
106  std::vector<std::pair<Key, Matrix> > terms;
107  for(Key key: keys()) {
108  terms.push_back(std::make_pair(key, this->A(key)));
109  }
110 
111  // compute rhs
112  Vector b = -error_vector(c);
113 
114  return boost::shared_ptr<GaussianFactor>(new JacobianFactor(terms, b, noiseModel::Unit::Create(dim())));
115 }
116 
117 /* ************************************************************************* */
119 
120  Vector errorVector = -b();
121 
122  for(Key key: this->keys()) {
123  const Value& newPt = c.at(key);
124  const Value& linPt = lin_points_.at(key);
125  Vector d = linPt.localCoordinates_(newPt);
126  const constABlock A = this->A(key);
127  errorVector += A * d;
128  }
129 
130  return errorVector;
131 }
132 
133 /* ************************************************************************* */
134 // LinearizedHessianFactor
135 /* ************************************************************************* */
137 }
138 
139 /* ************************************************************************* */
141  const HessianFactor::shared_ptr& hessian, const Values& lin_points)
142  : Base(hessian, lin_points), info_(hessian->info()) {}
143 
144 /* ************************************************************************* */
145 void LinearizedHessianFactor::print(const std::string& s, const KeyFormatter& keyFormatter) const {
146 
147  std::cout << s << std::endl;
148 
149  std::cout << "Nonlinear Keys: ";
150  for(const Key& key: this->keys())
151  std::cout << keyFormatter(key) << " ";
152  std::cout << std::endl;
153 
154  gtsam::print(Matrix(info_.selfadjointView()), "Ab^T * Ab: ");
155 
156  lin_points_.print("Linearization Point: ");
157 }
158 
159 /* ************************************************************************* */
161 
162  const This *e = dynamic_cast<const This*> (&expected);
163  if (e) {
164 
165  Matrix thisMatrix = this->info_.selfadjointView();
166  thisMatrix(thisMatrix.rows()-1, thisMatrix.cols()-1) = 0.0;
167  Matrix rhsMatrix = e->info_.selfadjointView();
168  rhsMatrix(rhsMatrix.rows()-1, rhsMatrix.cols()-1) = 0.0;
169 
170  return Base::equals(expected, tol)
171  && lin_points_.equals(e->lin_points_, tol)
172  && equal_with_abs_tol(thisMatrix, rhsMatrix, tol);
173  } else {
174  return false;
175  }
176 }
177 
178 /* ************************************************************************* */
180 
181  // Construct an error vector in key-order from the Values
182  Vector dx = Vector::Zero(dim());
183  size_t index = 0;
184  for(unsigned int i = 0; i < this->size(); ++i){
185  Key key = this->keys()[i];
186  const Value& newPt = c.at(key);
187  const Value& linPt = lin_points_.at(key);
188  dx.segment(index, linPt.dim()) = linPt.localCoordinates_(newPt);
189  index += linPt.dim();
190  }
191 
192  // error 0.5*(f - 2*x'*g + x'*G*x)
193  double f = constantTerm();
194  double xtg = dx.dot(linearTerm());
195  double xGx = dx.transpose() * squaredTerm() * dx;
196 
197  return 0.5 * (f - 2.0 * xtg + xGx);
198 }
199 
200 /* ************************************************************************* */
201 boost::shared_ptr<GaussianFactor>
203 
204  // Construct an error vector in key-order from the Values
205  Vector dx = Vector::Zero(dim());
206  size_t index = 0;
207  for(unsigned int i = 0; i < this->size(); ++i){
208  Key key = this->keys()[i];
209  const Value& newPt = c.at(key);
210  const Value& linPt = lin_points_.at(key);
211  dx.segment(index, linPt.dim()) = linPt.localCoordinates_(newPt);
212  index += linPt.dim();
213  }
214 
215  // f2 = f1 - 2*dx'*g1 + dx'*G1*dx
216  //newInfo(this->size(), this->size())(0,0) += -2*dx.dot(linearTerm()) + dx.transpose() * squaredTerm().selfadjointView<Eigen::Upper>() * dx;
217  double f = constantTerm() - 2*dx.dot(linearTerm()) + dx.transpose() * squaredTerm() * dx;
218 
219  // g2 = g1 - G1*dx
220  //newInfo.rangeColumn(0, this->size(), this->size(), 0) -= squaredTerm().selfadjointView<Eigen::Upper>() * dx;
221  Vector g = linearTerm() - squaredTerm() * dx;
222  std::vector<Vector> gs;
223  std::size_t offset = 0;
224  for(DenseIndex i = 0; i < info_.nBlocks()-1; ++i) {
225  const std::size_t dim = info_.getDim(i);
226  gs.push_back(g.segment(offset, dim));
227  offset += dim;
228  }
229 
230  // G2 = G1
231  // Do Nothing
232  std::vector<Matrix> Gs;
233  for(DenseIndex i = 0; i < info_.nBlocks()-1; ++i) {
234  Gs.push_back(info_.diagonalBlock(i));
235  for(DenseIndex j = i + 1; j < info_.nBlocks()-1; ++j) {
236  Gs.push_back(info_.aboveDiagonalBlock(i, j));
237  }
238  }
239 
240  // Create a Hessian Factor from the modified info matrix
241  //return boost::shared_ptr<GaussianFactor>(new HessianFactor(js, newInfo));
242  return boost::shared_ptr<GaussianFactor>(new HessianFactor(keys(), Gs, gs, f));
243 }
244 
245 } // \namespace aspn
void print(const Matrix &A, const string &s, ostream &stream)
Definition: Matrix.cpp:155
Vector error_vector(const Values &c) const
bool exists(Key j) const
Definition: Values.cpp:104
void insert(Key j, const Value &val)
Definition: Values.cpp:140
Matrix expected
Definition: testMatrix.cpp:974
bool equals(const NonlinearFactor &expected, double tol=1e-9) const override
const Matrix & matrix() const
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:43
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
const_iterator end() const
Definition: Factor.h:130
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 y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate offset
iterator iter(handle obj)
Definition: pytypes.h:1547
size_t dim() const override
virtual Vector localCoordinates_(const Value &value) const =0
bool equals(const Values &other, double tol=1e-9) const
Definition: Values.cpp:88
ptrdiff_t DenseIndex
The index type for Eigen objects.
Definition: types.h:67
void g(const string &key, int i)
Definition: testBTree.cpp:43
static shared_ptr Create(size_t dim)
Definition: NoiseModel.h:608
else if n * info
boost::shared_ptr< GaussianFactor > linearize(const Values &c) const override
DenseIndex nBlocks() const
Block count.
double error(const Values &c) const override
Eigen::VectorXd Vector
Definition: Vector.h:38
const ValueType at(Key j) const
Definition: Values-inl.h:342
const_iterator begin() const
Definition: Factor.h:127
size_t dim() const override
constColumn linearTerm() const
boost::shared_ptr< This > shared_ptr
A shared_ptr to this class.
A dummy factor that allows a linear factor to act as a nonlinear factor.
Eigen::SelfAdjointView< constBlock, Eigen::Upper > squaredTerm() const
std::function< std::string(Key)> KeyFormatter
Typedef for a function to format a key, i.e. to convert it to a string.
Definition: Key.h:35
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
Array< double, 1, 3 > e(1./3., 0.5, 2.)
RealScalar s
const constBVector b() const
constBlock aboveDiagonalBlock(DenseIndex I, DenseIndex J) const
Get block above the diagonal (I, J).
virtual bool equals(const NonlinearFactor &f, double tol=1e-9) const
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
void print(const std::string &s="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const override
virtual size_t dim() const =0
bool equal_with_abs_tol(const Eigen::DenseBase< MATRIX > &A, const Eigen::DenseBase< MATRIX > &B, double tol=1e-9)
Definition: base/Matrix.h:84
Block range(DenseIndex startBlock, DenseIndex endBlock)
traits
Definition: chartTesting.h:28
const constABlock A() const
void print(const std::string &str="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const
Definition: Values.cpp:77
bool equals(const NonlinearFactor &expected, double tol=1e-9) const override
size_t size() const
Definition: Factor.h:135
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
const KeyVector & keys() const
Access the factor&#39;s involved variable keys.
Definition: Factor.h:124
boost::shared_ptr< GaussianFactor > linearize(const Values &c) const override
double error(const Values &c) const override
Eigen::SelfAdjointView< constBlock, Eigen::Upper > selfadjointView(DenseIndex I, DenseIndex J) const
Return the square sub-matrix that contains blocks(i:j, i:j).
A Gaussian factor using the canonical parameters (information form)
const G double tol
Definition: Group.h:83
DenseIndex getDim(DenseIndex block) const
Number of dimensions for variable on this diagonal block.
Eigen::SelfAdjointView< Block, Eigen::Upper > diagonalBlock(DenseIndex J)
Return the J&#39;th diagonal block as a self adjoint view.
KeyVector::const_iterator const_iterator
Const iterator over keys.
Definition: Factor.h:67
const KeyVector keys
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
DenseIndex nBlocks() const
Block count.
void print(const std::string &s="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const override
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:61
std::ptrdiff_t j


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:42:31