35class CameraSet :
public std::vector<CAMERA, Eigen::aligned_allocator<CAMERA> > {
43 typedef typename CAMERA::Measurement
Z;
44 typedef typename CAMERA::MeasurementVector ZVector;
51 const ZVector& measured) {
54 size_t m = predicted.size();
55 if (measured.size() != m)
56 throw std::runtime_error(
"CameraSet::errors: size mismatch");
60 for (
size_t i = 0,
row = 0; i < m; i++,
row +=
ZDim) {
62 if(
ZDim==3 && std::isnan(bi(1))){
76 using MatrixZD = Eigen::Matrix<double, ZDim, D>;
77 using FBlocks = std::vector<MatrixZD, Eigen::aligned_allocator<MatrixZD>>;
84 virtual void print(
const std::string& s =
"")
const {
85 std::cout << s <<
"CameraSet, cameras = \n";
86 for (
size_t k = 0; k < this->size(); ++k)
92 if (this->size() != p.size())
94 bool camerasAreEqual =
true;
95 for (
size_t i = 0; i < this->size(); i++) {
96 if (this->at(i).equals(p.at(i), tol) ==
false)
97 camerasAreEqual =
false;
100 return camerasAreEqual;
109 template<
class POINT>
111 boost::optional<FBlocks&> Fs = boost::none,
112 boost::optional<Matrix&> E = boost::none)
const {
117 size_t m = this->size();
122 if (E) E->resize(
ZDim * m, N);
123 if (Fs) Fs->resize(m);
126 for (
size_t i = 0; i < m; i++) {
128 Eigen::Matrix<double, ZDim, N> Ei;
129 z.emplace_back(this->at(i).
project2(point, Fs ? &Fi : 0, E ? &Ei : 0));
130 if (Fs) (*Fs)[i] = Fi;
131 if (E) E->block<
ZDim, N>(
ZDim * i, 0) = Ei;
138 template<
class POINT>
140 boost::optional<FBlocks&> Fs = boost::none,
141 boost::optional<Matrix&> E = boost::none)
const {
155 Eigen::Matrix<double, ZDim, ND>,
156 Eigen::aligned_allocator<Eigen::Matrix<double, ZDim, ND>>>& Fs,
157 const Matrix& E,
const Eigen::Matrix<double, N, N>& P,
const Vector& b) {
159 size_t m = Fs.size();
162 size_t M1 = ND * m + 1;
163 std::vector<DenseIndex> dims(m + 1);
164 std::fill(dims.begin(), dims.end() - 1, ND);
169 for (
size_t i = 0; i < m; i++) {
171 const Eigen::Matrix<double, ZDim, ND>& Fi = Fs[i];
172 const auto FiT = Fi.transpose();
173 const Eigen::Matrix<double, ZDim, N> Ei_P =
178 - FiT * (Ei_P * (E.transpose() * b)));
182 * (Fi - Ei_P * E.block(
ZDim * i, 0,
ZDim, N).transpose() * Fi));
185 for (
size_t j = i + 1; j < m; j++) {
186 const Eigen::Matrix<double, ZDim, ND>& Fj = Fs[j];
190 * (Ei_P * E.block(
ZDim * j, 0,
ZDim, N).transpose() * Fj));
195 return augmentedHessian;
211 template <
int N,
int ND,
int NDD>
214 Eigen::Matrix<double, ZDim, ND>,
215 Eigen::aligned_allocator<Eigen::Matrix<double, ZDim, ND>>>& Fs,
216 const Matrix& E,
const Eigen::Matrix<double, N, N>& P,
const Vector& b,
218 size_t nrNonuniqueKeys = jacobianKeys.size();
219 size_t nrUniqueKeys = hessianKeys.size();
225 std::vector<DenseIndex> dims(nrUniqueKeys + 1);
226 std::fill(dims.begin(), dims.end() - 1, NDD);
231 if (nrUniqueKeys == nrNonuniqueKeys) {
238 std::vector<DenseIndex> nonuniqueDims(nrNonuniqueKeys + 1);
239 std::fill(nonuniqueDims.begin(), nonuniqueDims.end() - 1, NDD);
240 nonuniqueDims.back() = 1;
246 std::map<Key, size_t> keyToSlotMap;
247 for (
size_t k = 0; k < nrUniqueKeys; k++) {
248 keyToSlotMap[hessianKeys[k]] = k;
253 dims, Matrix::Zero(NDD * nrUniqueKeys + 1, NDD * nrUniqueKeys + 1));
258 for (
size_t i = 0; i < nrNonuniqueKeys; i++) {
259 Key key_i = jacobianKeys.at(i);
263 keyToSlotMap[key_i], nrUniqueKeys,
267 for (
size_t j = i; j < nrNonuniqueKeys; j++) {
268 Key key_j = jacobianKeys.at(j);
273 if (keyToSlotMap[key_i] != keyToSlotMap[key_j]) {
275 keyToSlotMap[key_i], keyToSlotMap[key_j],
289 nrUniqueKeys, augmentedHessian.
diagonalBlock(nrNonuniqueKeys));
291 return augmentedHessianUniqueKeys;
302 const Matrix& E,
const Eigen::Matrix<double, N, N>& P,
const Vector& b) {
303 return SchurComplement<N,D>(Fs, E, P, b);
309 const Matrix& E,
double lambda,
bool diagonalDamping =
false) {
311 Matrix EtE = E.transpose() * E;
313 if (diagonalDamping) {
314 EtE.diagonal() += lambda * EtE.diagonal();
317 EtE += lambda * Eigen::MatrixXd::Identity(n, n);
324 static Matrix
PointCov(
const Matrix& E,
const double lambda = 0.0,
325 bool diagonalDamping =
false) {
328 ComputePointCovariance<2>(P2, E, lambda, diagonalDamping);
332 ComputePointCovariance<3>(P3, E, lambda, diagonalDamping);
342 const Matrix& E,
const Vector& b,
const double lambda = 0.0,
343 bool diagonalDamping =
false) {
346 ComputePointCovariance<2>(P, E, lambda, diagonalDamping);
347 return SchurComplement<2>(Fblocks, E, P, b);
350 ComputePointCovariance<3>(P, E, lambda, diagonalDamping);
351 return SchurComplement<3>(Fblocks, E, P, b);
361 const Eigen::Matrix<double, N, N>& P,
const Vector& b,
365 assert(keys.size()==Fs.size());
366 assert(keys.size()<=allKeys.size());
369 for (
size_t slot = 0; slot < allKeys.size(); slot++)
370 KeySlotMap.insert(std::make_pair(allKeys[slot], slot));
377 size_t m = Fs.size();
378 size_t M = (augmentedHessian.
rows() - 1) /
D;
379 assert(allKeys.size()==M);
382 for (
size_t i = 0; i < m; i++) {
385 const auto FiT = Fi.transpose();
386 const Eigen::Matrix<double, 2, N> Ei_P = E.template block<ZDim, N>(
400 - FiT * (Ei_P * (E.transpose() * b)));
406 ((FiT * (Fi - Ei_P * E.template block<ZDim, N>(
ZDim * i, 0).transpose() * Fi))).eval());
409 for (
size_t j = i + 1; j < m; j++) {
419 -FiT * (Ei_P * E.template block<ZDim, N>(
ZDim * j, 0).transpose() * Fj));
430 template<
class ARCHIVE>
431 void serialize(ARCHIVE & ar,
const unsigned int ) {
439template<
class CAMERA>
442template<
class CAMERA>
445template<
class CAMERA>
449template<
class CAMERA>
#define GTSAM_MAKE_ALIGNED_OPERATOR_NEW
This marks a GTSAM object to require alignment.
Definition: types.h:277
A thin wrapper around std::map that uses boost's fast_pool_allocator.
Access to matrices via blocks of pre-defined sizes.
Concept check for values that can be used in unit tests.
Calibrated camera for which only pose is unknown.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:86
ptrdiff_t DenseIndex
The index type for Eigen objects.
Definition: types.h:75
const MATRIX::ConstRowXpr row(const MATRIX &A, size_t j)
Extracts a row view from a matrix that avoids a copy.
Definition: Matrix.h:225
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:69
A manifold defines a space in which there is a notion of a linear tangent space that can be centered ...
Definition: concepts.h:30
Give fixed size dimension of a type, fails at compile time if dynamic.
Definition: Manifold.h:164
Definition: SymmetricBlockMatrix.h:52
void setDiagonalBlock(DenseIndex I, const XprType &xpr)
Set a diagonal block. Only the upper triangular portion of xpr is evaluated.
Definition: SymmetricBlockMatrix.h:200
void setOffDiagonalBlock(DenseIndex I, DenseIndex J, const XprType &xpr)
Set an off-diagonal block. Only the upper triangular portion of xpr is evaluated.
Definition: SymmetricBlockMatrix.h:206
constBlock aboveDiagonalBlock(DenseIndex I, DenseIndex J) const
Get block above the diagonal (I, J).
Definition: SymmetricBlockMatrix.h:155
DenseIndex rows() const
Row size.
Definition: SymmetricBlockMatrix.h:119
void updateOffDiagonalBlock(DenseIndex I, DenseIndex J, const XprType &xpr)
Update an off diagonal block.
Definition: SymmetricBlockMatrix.h:233
Eigen::SelfAdjointView< Block, Eigen::Upper > diagonalBlock(DenseIndex J)
Return the J'th diagonal block as a self adjoint view.
Definition: SymmetricBlockMatrix.h:140
void updateDiagonalBlock(DenseIndex I, const XprType &xpr)
Increment the diagonal block by the values in xpr. Only reads the upper triangular part of xpr.
Definition: SymmetricBlockMatrix.h:217
Eigen::SelfAdjointView< constBlock, Eigen::Upper > selfadjointView(DenseIndex I, DenseIndex J) const
Return the square sub-matrix that contains blocks(i:j, i:j).
Definition: SymmetricBlockMatrix.h:161
A helper that implements the traits interface for GTSAM types.
Definition: Testable.h:151
A set of cameras, all with their own calibration.
Definition: CameraSet.h:35
Vector reprojectionError(const POINT &point, const ZVector &measured, boost::optional< FBlocks & > Fs=boost::none, boost::optional< Matrix & > E=boost::none) const
Calculate vector [project2(point)-z] of re-projection errors.
Definition: CameraSet.h:139
virtual void print(const std::string &s="") const
print
Definition: CameraSet.h:84
static void ComputePointCovariance(Eigen::Matrix< double, N, N > &P, const Matrix &E, double lambda, bool diagonalDamping=false)
Computes Point Covariance P, with lambda parameter.
Definition: CameraSet.h:308
static SymmetricBlockMatrix SchurComplement(const FBlocks &Fblocks, const Matrix &E, const Vector &b, const double lambda=0.0, bool diagonalDamping=false)
Do Schur complement, given Jacobian as Fs,E,P, return SymmetricBlockMatrix Dynamic version.
Definition: CameraSet.h:341
CAMERA::Measurement Z
2D measurement and noise model for each of the m views The order is kept the same as the keys that we...
Definition: CameraSet.h:43
virtual ~CameraSet()=default
Destructor.
bool equals(const CameraSet &p, double tol=1e-9) const
equals
Definition: CameraSet.h:91
static Matrix PointCov(const Matrix &E, const double lambda=0.0, bool diagonalDamping=false)
Computes Point Covariance P, with lambda parameter, dynamic version.
Definition: CameraSet.h:324
static Vector ErrorVector(const ZVector &predicted, const ZVector &measured)
Make a vector of re-projection errors.
Definition: CameraSet.h:50
static void UpdateSchurComplement(const FBlocks &Fs, const Matrix &E, const Eigen::Matrix< double, N, N > &P, const Vector &b, const KeyVector &allKeys, const KeyVector &keys, SymmetricBlockMatrix &augmentedHessian)
Applies Schur complement (exploiting block structure) to get a smart factor on cameras,...
Definition: CameraSet.h:360
static SymmetricBlockMatrix SchurComplementAndRearrangeBlocks(const std::vector< Eigen::Matrix< double, ZDim, ND >, Eigen::aligned_allocator< Eigen::Matrix< double, ZDim, ND > > > &Fs, const Matrix &E, const Eigen::Matrix< double, N, N > &P, const Vector &b, const KeyVector &jacobianKeys, const KeyVector &hessianKeys)
Do Schur complement, given Jacobian as Fs,E,P, return SymmetricBlockMatrix G = F' * F - F' * E * P * ...
Definition: CameraSet.h:212
static SymmetricBlockMatrix SchurComplement(const std::vector< Eigen::Matrix< double, ZDim, ND >, Eigen::aligned_allocator< Eigen::Matrix< double, ZDim, ND > > > &Fs, const Matrix &E, const Eigen::Matrix< double, N, N > &P, const Vector &b)
Do Schur complement, given Jacobian as Fs,E,P, return SymmetricBlockMatrix G = F' * F - F' * E * P * ...
Definition: CameraSet.h:153
static const int ZDim
Measurement dimension.
Definition: CameraSet.h:47
static SymmetricBlockMatrix SchurComplement(const FBlocks &Fs, const Matrix &E, const Eigen::Matrix< double, N, N > &P, const Vector &b)
Do Schur complement, given Jacobian as Fs,E,P, return SymmetricBlockMatrix G = F' * F - F' * E * P * ...
Definition: CameraSet.h:301
static const int D
Camera dimension.
Definition: CameraSet.h:46
friend class boost::serialization::access
Serialization function.
Definition: CameraSet.h:429
Eigen::Matrix< double, ZDim, D > MatrixZD
Definitions for blocks of F.
Definition: CameraSet.h:76
ZVector project2(const POINT &point, boost::optional< FBlocks & > Fs=boost::none, boost::optional< Matrix & > E=boost::none) const
Project a point (possibly Unit3 at infinity), with derivatives Note that F is a sparse block-diagonal...
Definition: CameraSet.h:110