cpd  0.5.1
Coherent Point Drift: C++ library for point set registration
nonrigid.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 double DEFAULT_BETA = 3.0;
31 const double DEFAULT_LAMBDA = 3.0;
32 
34 struct NonrigidResult : public Result {};
35 
37 class Nonrigid : public Transform<NonrigidResult> {
38 public:
39  Nonrigid()
40  : Transform()
41  , m_lambda(DEFAULT_LAMBDA)
42  , m_beta(DEFAULT_BETA)
43  , m_linked(DEFAULT_LINKED) {}
44 
46  void init(const Matrix& fixed, const Matrix& moving);
47 
49  void modify_probabilities(Probabilities& probabilities) const;
50 
52  Nonrigid& beta(double beta) {
53  m_beta = beta;
54  return *this;
55  }
56 
58  Nonrigid& lambda(double lambda) {
59  m_lambda = lambda;
60  return *this;
61  }
62 
64  NonrigidResult compute_one(const Matrix& fixed, const Matrix& moving,
65  const Probabilities& probabilities,
66  double sigma2) const;
67 
69  Nonrigid& linked(bool linked) {
70  m_linked = linked;
71  return *this;
72  }
73 
74  virtual bool linked() const { return m_linked; }
75 
76 private:
77  Matrix m_g;
78  Matrix m_w;
79  double m_lambda;
80  double m_beta;
81  bool m_linked;
82 };
83 
85 NonrigidResult nonrigid(const Matrix& fixed, const Matrix& moving);
86 } // namespace cpd
Nonrigid & lambda(double lambda)
Sets the lambda.
Definition: nonrigid.hpp:58
Nonrigid coherent point drift.
Definition: nonrigid.hpp:37
Nonrigid & linked(bool linked)
Sets whether the scalings of the two datasets are linked.
Definition: nonrigid.hpp:69
virtual bool linked() const
Returns true if the normalization should be linked.
Definition: nonrigid.hpp:74
Nonrigid & beta(double beta)
Sets the beta.
Definition: nonrigid.hpp:52
NonrigidResult nonrigid(const Matrix &fixed, const Matrix &moving)
Runs a nonrigid 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
const double DEFAULT_LAMBDA
Default value for lambda.
Definition: nonrigid.hpp:31
Generic coherent point drift transform.
Definition: transform.hpp:75
The result of a generic transform run.
Definition: transform.hpp:52
const double DEFAULT_BETA
Default value for beta.
Definition: nonrigid.hpp:29
const bool DEFAULT_LINKED
Are the scalings of the two datasets linked by default?
Definition: transform.hpp:49
The result of a nonrigid coherent point drift run.
Definition: nonrigid.hpp:34
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