cpd  0.5.1
Coherent Point Drift: C++ library for point set registration
rigid.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 const bool DEFAULT_REFLECTIONS = false;
32 
34 struct RigidResult : public Result {
40  double scale;
41 
44  Matrix matrix() const;
45 
46  void denormalize(const Normalization& normalization);
47 };
48 
52 class Rigid : public Transform<RigidResult> {
53 public:
54  Rigid()
55  : Transform()
56  , m_reflections(DEFAULT_REFLECTIONS)
57  , m_scale(DEFAULT_SCALE) {}
58 
60  Rigid& reflections(bool reflections) {
61  m_reflections = reflections;
62  return *this;
63  }
64 
66  Rigid& scale(bool scale) {
67  m_scale = scale;
68  return *this;
69  }
70 
72  RigidResult compute_one(const Matrix& fixed, const Matrix& moving,
73  const Probabilities& probabilities,
74  double sigma2) const;
75 
76  virtual bool linked() const { return !m_scale; }
77 
78 private:
79  bool m_reflections;
80  bool m_scale;
81 };
82 
84 RigidResult rigid(const Matrix& fixed, const Matrix& moving);
85 } // namespace cpd
const bool DEFAULT_REFLECTIONS
Should rigid registrations allow reflections by default?
Definition: rigid.hpp:29
RigidResult rigid(const Matrix &fixed, const Matrix &moving)
Runs a rigid registration on two matrices.
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
Rigid & scale(bool scale)
Sets whether this rigid transform allows scaling.
Definition: rigid.hpp:66
Vector translation
The translation component of the transformation.
Definition: rigid.hpp:38
virtual bool linked() const
Returns true if the normalization should be linked.
Definition: rigid.hpp:76
const bool DEFAULT_SCALE
Should rigid registrations scale the data by default?
Definition: rigid.hpp:31
Matrix rotation
The rotation component of the transformation.
Definition: rigid.hpp:36
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
Rigid coherent point drift.
Definition: rigid.hpp:52
void denormalize(const Normalization &normalization)
De-normalize this result.
Generic coherent point drift transform.
Top-level cpd namespace.
Definition: affine.hpp:26
The result of a rigid coherent point drift run.
Definition: rigid.hpp:34
Rigid & reflections(bool reflections)
Sets whether this rigid transform allows reflections.
Definition: rigid.hpp:60
Matrix matrix() const
Returns a single matrix that contains all the transformation information.
double scale
The scaling component of the transformation.
Definition: rigid.hpp:40
Probability matrices produced by comparing two data sets with a GaussTransform.
Definition: gauss_transform.hpp:32