SubgraphPreconditioner.cpp
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2 
3  * GTSAM Copyright 2010-2019, 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 
24 #include <gtsam/base/types.h>
25 #include <gtsam/base/Vector.h>
26 
27 #include <stdexcept>
28 #include <cassert>
29 
30 using std::cout;
31 using std::endl;
32 using std::vector;
33 using std::ostream;
34 
35 namespace gtsam {
36 
37 /* ************************************************************************* */
38 static Vector getSubvector(const Vector &src, const KeyInfo &keyInfo,
39  const KeyVector &keys) {
40  /* a cache of starting index and dim */
41  vector<std::pair<size_t, size_t> > cache;
42  cache.reserve(3);
43 
44  /* figure out dimension by traversing the keys */
45  size_t dim = 0;
46  for (const Key &key : keys) {
47  const KeyInfoEntry &entry = keyInfo.find(key)->second;
48  cache.emplace_back(entry.start, entry.dim);
49  dim += entry.dim;
50  }
51 
52  /* use the cache to fill the result */
53  Vector result(dim);
54  size_t index = 0;
55  for (const auto &p : cache) {
56  result.segment(index, p.second) = src.segment(p.first, p.second);
57  index += p.second;
58  }
59 
60  return result;
61 }
62 
63 /*****************************************************************************/
64 static void setSubvector(const Vector &src, const KeyInfo &keyInfo,
65  const KeyVector &keys, Vector &dst) {
66  size_t index = 0;
67  for (const Key &key : keys) {
68  const KeyInfoEntry &entry = keyInfo.find(key)->second;
69  const size_t keyDim = entry.dim;
70  dst.segment(entry.start, keyDim) = src.segment(index, keyDim);
71  index += keyDim;
72  }
73 }
74 
75 /* ************************************************************************* */
76 // Convert any non-Jacobian factors to Jacobians (e.g. Hessian -> Jacobian with
77 // Cholesky)
79  const GaussianFactorGraph &gfg) {
81  for (const auto &factor : gfg)
82  if (factor) {
83  auto jf = std::dynamic_pointer_cast<JacobianFactor>(factor);
84  if (!jf) {
85  jf = std::make_shared<JacobianFactor>(*factor);
86  }
87  result.push_back(jf);
88  }
89  return result;
90 }
91 
92 /* ************************************************************************* */
94  parameters_(p) {}
95 
96 /* ************************************************************************* */
98  const GaussianBayesNet& Rc1, const VectorValues& xbar, const SubgraphPreconditionerParameters &p) :
99  Ab2_(convertToJacobianFactors(Ab2)), Rc1_(Rc1), xbar_(xbar),
100  b2bar_(-Ab2_.gaussianErrors(xbar)), parameters_(p) {
101 }
102 
103 /* ************************************************************************* */
104 // x = xbar + inv(R1)*y
106  return xbar_ + Rc1_.backSubstitute(y);
107 }
108 
109 /* ************************************************************************* */
111  Errors e = createErrors(y);
112  VectorValues x = this->x(y);
113  Errors e2 = Ab2_.gaussianErrors(x);
114  return 0.5 * (dot(e, e) + dot(e2,e2));
115 }
116 
117 /* ************************************************************************* */
118 // gradient is y + inv(R1')*A2'*(A2*inv(R1)*y-b2bar),
120  VectorValues x = Rc1_.backSubstitute(y); /* inv(R1)*y */
121  Errors e = Ab2_ * x - b2bar_; /* (A2*inv(R1)*y-b2bar) */
123  Ab2_.transposeMultiplyAdd(1.0, e, v); /* A2'*(A2*inv(R1)*y-b2bar) */
124  return y + Rc1_.backSubstituteTranspose(v);
125 }
126 
127 /* ************************************************************************* */
128 // Apply operator A, A*y = [I;A2*inv(R1)]*y = [y; A2*inv(R1)*y]
130  Errors e = createErrors(y);
131  VectorValues x = Rc1_.backSubstitute(y); /* x=inv(R1)*y */
132  Errors e2 = Ab2_ * x; /* A2*x */
133  e.splice(e.end(), e2);
134  return e;
135 }
136 
137 /* ************************************************************************* */
138 // In-place version that overwrites e
140 
141  Errors::iterator ei = e.begin();
142  for(const auto& key_value: y) {
143  *ei = key_value.second;
144  ++ei;
145  }
146 
147  // Add A2 contribution
148  VectorValues x = Rc1_.backSubstitute(y); // x=inv(R1)*y
149  Ab2_.multiplyInPlace(x, ei); // use iterator version
150 }
151 
152 /* ************************************************************************* */
153 // Apply operator A', A'*e = [I inv(R1')*A2']*e = e1 + inv(R1')*A2'*e2
155 
156  Errors::const_iterator it = e.begin();
157  VectorValues y = zero();
158  for(auto& key_value: y) {
159  key_value.second = *it;
160  ++it;
161  }
162  transposeMultiplyAdd2(1.0, it, e.end(), y);
163  return y;
164 }
165 
166 /* ************************************************************************* */
167 // y += alpha*A'*e
169 (double alpha, const Errors& e, VectorValues& y) const {
170 
171  Errors::const_iterator it = e.begin();
172  for(auto& key_value: y) {
173  const Vector& ei = *it;
174  key_value.second += alpha * ei;
175  ++it;
176  }
177  transposeMultiplyAdd2(alpha, it, e.end(), y);
178 }
179 
180 /* ************************************************************************* */
181 // y += alpha*inv(R1')*A2'*e2
183  Errors::const_iterator it, Errors::const_iterator end, VectorValues& y) const {
184 
185  // create e2 with what's left of e
186  // TODO can we avoid creating e2 by passing iterator to transposeMultiplyAdd ?
187  Errors e2;
188  while (it != end) e2.push_back(*(it++));
189 
190  VectorValues x = VectorValues::Zero(y); // x = 0
191  Ab2_.transposeMultiplyAdd(1.0,e2,x); // x += A2'*e2
192  y += alpha * Rc1_.backSubstituteTranspose(x); // y += alpha*inv(R1')*x
193 }
194 
195 /* ************************************************************************* */
196 void SubgraphPreconditioner::print(const std::string& s) const {
197  cout << s << endl;
198  Ab2_.print();
199 }
200 
201 /*****************************************************************************/
203  assert(x.size() == y.size());
204 
205  /* back substitute */
206  for (auto it = std::make_reverse_iterator(Rc1_.end()); it != std::make_reverse_iterator(Rc1_.begin()); ++it) {
207  auto& cg = *it;
208  /* collect a subvector of x that consists of the parents of cg (S) */
209  const KeyVector parentKeys(cg->beginParents(), cg->endParents());
210  const KeyVector frontalKeys(cg->beginFrontals(), cg->endFrontals());
211  const Vector xParent = getSubvector(x, keyInfo_, parentKeys);
212  const Vector rhsFrontal = getSubvector(y, keyInfo_, frontalKeys);
213 
214  /* compute the solution for the current pivot */
215  const Vector solFrontal = cg->R().triangularView<Eigen::Upper>().solve(
216  rhsFrontal - cg->S() * xParent);
217 
218  /* assign subvector of sol to the frontal variables */
219  setSubvector(solFrontal, keyInfo_, frontalKeys, x);
220  }
221 }
222 
223 /*****************************************************************************/
225  /* copy first */
226  assert(x.size() == y.size());
227  std::copy(y.data(), y.data() + y.rows(), x.data());
228 
229  /* in place back substitute */
230  for (const auto &cg : Rc1_) {
231  const KeyVector frontalKeys(cg->beginFrontals(), cg->endFrontals());
232  const Vector rhsFrontal = getSubvector(x, keyInfo_, frontalKeys);
233  const Vector solFrontal =
234  cg->R().transpose().triangularView<Eigen::Lower>().solve(
235  rhsFrontal);
236 
237  // Check for indeterminant solution
238  if (solFrontal.hasNaN())
239  throw IndeterminantLinearSystemException(cg->keys().front());
240 
241  /* assign subvector of sol to the frontal variables */
242  setSubvector(solFrontal, keyInfo_, frontalKeys, x);
243 
244  /* substract from parent variables */
245  for (auto it = cg->beginParents(); it != cg->endParents(); it++) {
246  const KeyInfoEntry &entry = keyInfo_.find(*it)->second;
247  Eigen::Map<Vector> rhsParent(x.data() + entry.start, entry.dim, 1);
248  rhsParent -= Matrix(cg->getA(it)).transpose() * solFrontal;
249  }
250  }
251 }
252 
253 /*****************************************************************************/
254 void SubgraphPreconditioner::build(const GaussianFactorGraph &gfg, const KeyInfo &keyInfo, const std::map<Key,Vector> &lambda)
255 {
256  /* identify the subgraph structure */
258  auto subgraph = builder(gfg);
259 
260  keyInfo_ = keyInfo;
261 
262  /* build factor subgraph */
263  auto gfg_subgraph = buildFactorSubgraph(gfg, subgraph, true);
264 
265  /* factorize and cache BayesNet */
266  Rc1_ = *gfg_subgraph.eliminateSequential();
267 }
268 
269 /*****************************************************************************/
270 
271 } // nsamespace gtsam
gtsam::buildFactorSubgraph
GaussianFactorGraph buildFactorSubgraph(const GaussianFactorGraph &gfg, const Subgraph &subgraph, const bool clone)
Definition: SubgraphBuilder.cpp:407
gtsam::SubgraphPreconditioner::operator*
Errors operator*(const VectorValues &y) const
Definition: SubgraphPreconditioner.cpp:129
GaussianFactorGraph.h
Linear Factor Graph where all factors are Gaussians.
Vector.h
typedef and functions to augment Eigen's VectorXd
alpha
RealScalar alpha
Definition: level1_cplx_impl.h:147
s
RealScalar s
Definition: level1_cplx_impl.h:126
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
types.h
Typedefs for easier changing of types.
keys
const KeyVector keys
Definition: testRegularImplicitSchurFactor.cpp:40
gtsam::SubgraphPreconditioner::operator^
VectorValues operator^(const Errors &e) const
Definition: SubgraphPreconditioner.cpp:154
gtsam::SubgraphPreconditioner::SubgraphPreconditioner
SubgraphPreconditioner(const SubgraphPreconditionerParameters &p=SubgraphPreconditionerParameters())
Definition: SubgraphPreconditioner.cpp:93
gtsam::FactorGraph::print
virtual void print(const std::string &s="FactorGraph", const KeyFormatter &formatter=DefaultKeyFormatter) const
Print out graph to std::cout, with optional key formatter.
Definition: FactorGraph-inst.h:37
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::KeyInfoEntry::dim
size_t dim
Definition: IterativeSolver.h:116
gtsam::SubgraphPreconditioner::error
double error(const VectorValues &y) const
Definition: SubgraphPreconditioner.cpp:110
gtsam::SubgraphPreconditioner::b2bar_
Errors b2bar_
A2*xbar - b2.
Definition: SubgraphPreconditioner.h:63
Eigen::Upper
@ Upper
Definition: Constants.h:211
gtsam::Matrix
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
GaussianBayesNet.h
Chordal Bayes Net, the result of eliminating a factor graph.
gtsam::setSubvector
static void setSubvector(const Vector &src, const KeyInfo &keyInfo, const KeyVector &keys, Vector &dst)
Definition: SubgraphPreconditioner.cpp:64
gtsam::getSubvector
static Vector getSubvector(const Vector &src, const KeyInfo &keyInfo, const KeyVector &keys)
Definition: SubgraphPreconditioner.cpp:38
gtsam::Vector
Eigen::VectorXd Vector
Definition: Vector.h:39
gtsam::KeyVector
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:92
result
Values result
Definition: OdometryOptimize.cpp:8
gtsam::KeyInfo
Definition: IterativeSolver.h:127
pruning_fixture::factor
DecisionTreeFactor factor(D &C &B &A, "0.0 0.0 0.0 0.60658897 0.61241912 0.61241969 0.61247685 0.61247742 0.0 " "0.0 0.0 0.99995287 1.0 1.0 1.0 1.0")
gtsam::SubgraphPreconditioner::x
VectorValues x(const VectorValues &y) const
Definition: SubgraphPreconditioner.cpp:105
gtsam::GaussianFactorGraph
Definition: GaussianFactorGraph.h:73
gtsam::SubgraphPreconditioner::transposeMultiplyAdd
void transposeMultiplyAdd(double alpha, const Errors &e, VectorValues &y) const
Definition: SubgraphPreconditioner.cpp:169
gtsam::SubgraphPreconditioner::Rc1_
GaussianBayesNet Rc1_
Definition: SubgraphPreconditioner.h:61
gtsam::VectorValues
Definition: VectorValues.h:74
gtsam::SubgraphPreconditioner::gradient
VectorValues gradient(const VectorValues &y) const
Definition: SubgraphPreconditioner.cpp:119
gtsam::SubgraphPreconditioner::Ab2_
GaussianFactorGraph Ab2_
Definition: SubgraphPreconditioner.h:60
gtsam::SubgraphPreconditioner::print
void print(const std::string &s="SubgraphPreconditioner") const
Definition: SubgraphPreconditioner.cpp:196
gtsam::SubgraphPreconditioner::transposeSolve
void transposeSolve(const Vector &y, Vector &x) const override
implement x = R^{-T} y
Definition: SubgraphPreconditioner.cpp:224
gtsam::createErrors
Errors createErrors(const VectorValues &V)
Break V into pieces according to its start indices.
Definition: Errors.cpp:28
gtsam::SubgraphPreconditioner::multiplyInPlace
void multiplyInPlace(const VectorValues &y, Errors &e) const
Definition: SubgraphPreconditioner.cpp:139
gtsam::dot
double dot(const V1 &a, const V2 &b)
Definition: Vector.h:196
gtsam::VectorValues::Zero
static VectorValues Zero(const VectorValues &other)
Definition: VectorValues.cpp:77
gtsam::SubgraphPreconditioner::keyInfo_
KeyInfo keyInfo_
Definition: SubgraphPreconditioner.h:65
gtsam::SubgraphBuilder
Definition: SubgraphBuilder.h:153
gtsam::GaussianBayesNet::backSubstituteTranspose
VectorValues backSubstituteTranspose(const VectorValues &gx) const
Definition: GaussianBayesNet.cpp:145
gtsam::FastList< Vector >
Eigen::Lower
@ Lower
Definition: Constants.h:209
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
gtsam::SubgraphPreconditionerParameters::builderParams
SubgraphBuilderParameters builderParams
Definition: SubgraphPreconditioner.h:44
lambda
static double lambda[]
Definition: jv.c:524
gtsam::GaussianFactorGraph::gaussianErrors
Errors gaussianErrors(const VectorValues &x) const
Definition: GaussianFactorGraph.cpp:512
y
Scalar * y
Definition: level1_cplx_impl.h:124
key
const gtsam::Symbol key('X', 0)
JacobianFactor.h
gtsam::SubgraphPreconditionerParameters
Definition: SubgraphPreconditioner.h:40
gtsam::SubgraphPreconditioner::transposeMultiplyAdd2
void transposeMultiplyAdd2(double alpha, Errors::const_iterator begin, Errors::const_iterator end, VectorValues &y) const
Definition: SubgraphPreconditioner.cpp:182
gtsam::SubgraphPreconditioner::solve
void solve(const Vector &y, Vector &x) const override
implement x = R^{-1} y
Definition: SubgraphPreconditioner.cpp:202
gtsam::convertToJacobianFactors
static GaussianFactorGraph convertToJacobianFactors(const GaussianFactorGraph &gfg)
Definition: SubgraphPreconditioner.cpp:78
gtsam::GaussianFactorGraph::transposeMultiplyAdd
void transposeMultiplyAdd(double alpha, const Errors &e, VectorValues &x) const
Definition: GaussianFactorGraph.cpp:455
gtsam
traits
Definition: SFMdata.h:40
make_reverse_iterator
std::reverse_iterator< Iterator > make_reverse_iterator(Iterator i)
Definition: stl_iterators.cpp:16
gtsam::FactorGraph::end
const_iterator end() const
Definition: FactorGraph.h:342
SubgraphPreconditioner.h
p
float * p
Definition: Tutorial_Map_using.cpp:9
gtsam::GaussianBayesNet::backSubstitute
VectorValues backSubstitute(const VectorValues &gx) const
Definition: GaussianBayesNet.cpp:129
gtsam::VectorValues::size
size_t size() const
Definition: VectorValues.h:130
v
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
SubgraphBuilder.h
gtsam::SubgraphPreconditioner::build
void build(const GaussianFactorGraph &gfg, const KeyInfo &info, const std::map< Key, Vector > &lambda) override
build/factorize the preconditioner
Definition: SubgraphPreconditioner.cpp:254
gtsam::FactorGraph::begin
const_iterator begin() const
Definition: FactorGraph.h:339
gtsam::SubgraphPreconditioner::parameters_
SubgraphPreconditionerParameters parameters_
Definition: SubgraphPreconditioner.h:66
gtsam::IndeterminantLinearSystemException
Definition: linearExceptions.h:94
Eigen::placeholders::end
static const EIGEN_DEPRECATED end_t end
Definition: IndexedViewHelper.h:181
gtsam::KeyInfoEntry::start
size_t start
Definition: IterativeSolver.h:116
gtsam::GaussianFactorGraph::multiplyInPlace
void multiplyInPlace(const VectorValues &x, Errors &e) const
‍** In-place version e <- A*x that overwrites e. *‍/
Definition: GaussianFactorGraph.cpp:427
gtsam::Key
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:97
gtsam::SubgraphPreconditioner::xbar_
VectorValues xbar_
A1 \ b1.
Definition: SubgraphPreconditioner.h:62
gtsam::KeyInfoEntry
Definition: IterativeSolver.h:115
gtsam::GaussianBayesNet
Definition: GaussianBayesNet.h:35
gtsam::SubgraphPreconditioner::zero
VectorValues zero() const
Definition: SubgraphPreconditioner.h:104


gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:14:04