cpd  0.5.1
Coherent Point Drift: C++ library for point set registration
affine.hpp
Go to the documentation of this file.
1 // cpd - Coherent Point Drift
2 // Copyright (C) 2017 Pete Gadomski <pete.gadomski@gmail.com>
3 //
4 // This program is free software; you can redistribute it and/or modify
5 // it under the terms of the GNU General Public License as published by
6 // the Free Software Foundation; either version 2 of the License, or
7 // (at your option) any later version.
8 //
9 // This program is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 // GNU General Public License for more details.
13 //
14 // You should have received a copy of the GNU General Public License along
15 // with this program; if not, write to the Free Software Foundation, Inc.,
16 // 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
17 
21 
22 #pragma once
23 
24 #include <cpd/transform.hpp>
25 
26 namespace cpd {
27 
29 struct AffineResult : public Result {
32 
35 
37  Matrix matrix() const;
38 
40  void denormalize(const Normalization& normalization);
41 };
42 
44 class Affine : public Transform<AffineResult> {
45 public:
46  Affine()
47  : Transform()
48  , m_linked(DEFAULT_LINKED) {}
49 
51  AffineResult compute_one(const Matrix& fixed, const Matrix& moving,
52  const Probabilities& probabilities,
53  double sigma2) const;
54 
56  Affine& linked(bool linked) {
57  m_linked = linked;
58  return *this;
59  }
60 
61  virtual bool linked() const { return m_linked; }
62 
63 private:
64  bool m_linked;
65 };
66 
68 AffineResult affine(const Matrix& fixed, const Matrix& moving);
69 } // namespace cpd
Matrix transform
The affine transformation.
Definition: affine.hpp:31
Affine coherent point drift.
Definition: affine.hpp:44
Affine & linked(bool linked)
Sets whether the scalings of the two datasets are linked.
Definition: affine.hpp:56
Matrix matrix() const
Returns the transform and the translation as one matrix.
The result of a affine coherent point drift run.
Definition: affine.hpp:29
double sigma2
The final sigma2 value.
Definition: transform.hpp:56
Eigen::MatrixXd Matrix
Our base matrix class.
Definition: matrix.hpp:29
Eigen::VectorXd Vector
Typedef for our specific type of vector.
Definition: matrix.hpp:32
AffineResult affine(const Matrix &fixed, const Matrix &moving)
Runs a affine registration on two matrices.
virtual bool linked() const
Returns true if the normalization should be linked.
Definition: affine.hpp:61
Vector translation
The translation vector.
Definition: affine.hpp:34
void denormalize(const Normalization &normalization)
Denormalize this result.
Generic coherent point drift transform.
Definition: transform.hpp:75
The result of a generic transform run.
Definition: transform.hpp:52
const bool DEFAULT_LINKED
Are the scalings of the two datasets linked by default?
Definition: transform.hpp:49
The results of normalizing data to a unit cube (or whatever dimensionality).
Definition: normalization.hpp:30
Generic coherent point drift transform.
Top-level cpd namespace.
Definition: affine.hpp:26
Probability matrices produced by comparing two data sets with a GaussTransform.
Definition: gauss_transform.hpp:32