LinearFreeFormTransformation4D.h
1 /*
2  * Medical Image Registration ToolKit (MIRTK)
3  *
4  * Copyright 2008-2015 Imperial College London
5  * Copyright 2008-2013 Daniel Rueckert, Julia Schnabel
6  * Copyright 2013-2015 Andreas Schuh
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef MIRTK_LinearFreeFormTransformation4D_H
22 #define MIRTK_LinearFreeFormTransformation4D_H
23 
24 #include "mirtk/FreeFormTransformation4D.h"
25 #include "mirtk/LinearInterpolateImageFunction4D.h"
26 
27 #include <cmath>
28 
29 
30 namespace mirtk {
31 
32 
33 // Forward declaration
34 class BSplineFreeFormTransformation4D;
35 
36 
37 /**
38  * Free-form transformation with linear 4D interpolation kernel
39  */
41 {
42  mirtkTransformationMacro(LinearFreeFormTransformation4D);
43 
44  // ---------------------------------------------------------------------------
45  // Types
46 
48 
49  // ---------------------------------------------------------------------------
50  // Data members
51 
52 protected:
53 
54  /// Interpolates control point values at arbitrary lattice locations
55  Interpolator _FFD;
56 
57 public:
58 
59  // ---------------------------------------------------------------------------
60  // Construction/Destruction
61 
62  /// Default constructor
64 
65  /// Construct free-form transformation for given image domain and lattice spacing
66  LinearFreeFormTransformation4D(double, double, double, double,
67  double, double, double, double,
68  double, double, double, double,
69  double *, double *, double *);
70 
71  /// Construct free-form transformation for given image domain and lattice spacing
73  double = -1, double = -1, double = -1, double = -1);
74 
75  /// Construct free-form transformation for given target image and lattice spacing
77  double, double, double, double);
78 
79  /// Construct free-form transformation from given B-spline FFD
81 
82  /// Copy Constructor
84 
85  /// Destructor
87 
88  // ---------------------------------------------------------------------------
89  // Approximation/Interpolation
90 
91  /// Approximate displacements: This function takes a set of points and a set
92  /// of displacements and finds !new! parameters such that the resulting
93  /// transformation approximates the displacements as good as possible.
94  virtual void ApproximateDOFs(const double *, const double *, const double *, const double *,
95  const double *, const double *, const double *, int);
96 
97  /// Finds gradient of approximation error: This function takes a set of points
98  /// and a set of errors. It finds a gradient w.r.t. the transformation parameters
99  /// which minimizes the L2 norm of the approximation error and adds it to the
100  /// input gradient with the given weight.
101  virtual void ApproximateDOFsGradient(const double *, const double *, const double *, const double *,
102  const double *, const double *, const double *, int,
103  double *, double = 1.0) const;
104 
105  /// Interpolates displacements: This function takes a set of displacements defined
106  /// at the control points and finds a FFD which interpolates these displacements.
107  virtual void Interpolate(const double *, const double *, const double *);
108 
109  // ---------------------------------------------------------------------------
110  // Lattice
111 
112  // Import other overloads
114 
115  /// Gets the temporal bounding box for a control point. The last parameter
116  /// specifies what fraction of the bounding box to return. The default
117  /// is 1 which equals 100% of the bounding box.
118  virtual void BoundingBox(int, double &, double &, double = 1) const;
119 
120  /// Gets the spatial bounding box for a control point. The last parameter
121  /// specifies what fraction of the bounding box to return. The default
122  /// is 1 which equals 100% of the bounding box.
123  virtual void BoundingBox(int, double &, double &, double &,
124  double &, double &, double &, double = 1) const;
125 
126  /// Gets the spatio-temporal bounding box for a control point. The last parameter
127  /// specifies what fraction of the bounding box to return. The default
128  /// is 1 which equals 100% of the bounding box.
129  virtual void BoundingBox(int, double &, double &, double &, double &,
130  double &, double &, double &, double &, double = 1) const;
131 
132  // ---------------------------------------------------------------------------
133  // Evaluation
134 
135  /// Evaluates the free-form transformation at a point in lattice coordinates
136  void Evaluate(double &, double &, double &, double) const;
137 
138  /// Evaluates the free-form transformation at a point in lattice coordinates
139  void EvaluateInside(double &, double &, double &, double) const;
140 
141  /// Calculates the Jacobian of the local transformation w.r.t world coordinates
142  void EvaluateJacobian(Matrix &, double, double, double, double) const;
143 
144  // ---------------------------------------------------------------------------
145  // Point transformation
146 
147  /// Transforms a single point using the local transformation only
148  virtual void LocalTransform(double &, double &, double &, double, double = 1) const;
149 
150  /// Transforms a single point using the inverse of the local transformation only
151  virtual bool LocalInverse(double &, double &, double &, double, double = 1) const;
152 
153  // ---------------------------------------------------------------------------
154  // Derivatives
155 
156  // Do not hide base class methods
159 
160  /// Calculates the Jacobian of the local transformation w.r.t world coordinates
161  virtual void LocalJacobian(Matrix &, double, double, double, double, double = NaN) const;
162 
163  /// Calculates the Jacobian of the transformation w.r.t the transformation parameters
164  virtual void JacobianDOFs(double [3], int, int, int, int, double, double, double, double) const;
165 
166  // ---------------------------------------------------------------------------
167  // Properties
169 
170  /// Size of support region of the used kernel
171  virtual int KernelSize() const;
172 
173 protected:
174 
175  /// Calculate the bending energy of the transformation at control points
176  double BendingEnergy(int i, int j, int k, int l) const;
177 
178 public:
179 
180  /// Approximates the bending energy on the control point lattice
181  virtual double BendingEnergy(bool = false, bool = true) const;
182 
183  /// Approximates the gradient of the bending energy on the control point
184  /// lattice w.r.t the transformation parameters and adds it with the given weight
185  virtual void BendingEnergyGradient(double *, double = NaN, bool = false, bool = true, bool = true) const;
186 
187  // ---------------------------------------------------------------------------
188  // I/O
189 
190  // Do not hide methods of base class
192 
193  /// Prints the parameters of the transformation
194  virtual void Print(ostream &, Indent = 0) const;
195 
196  /// Whether this transformation can read a file of specified type (i.e. format)
197  virtual bool CanRead(TransformationType) const;
198 
199 };
200 
201 ////////////////////////////////////////////////////////////////////////////////
202 // Inline definitions
203 ////////////////////////////////////////////////////////////////////////////////
204 
205 // =============================================================================
206 // Evaluation
207 // =============================================================================
208 
209 // -----------------------------------------------------------------------------
211 ::Evaluate(double &x, double &y, double &z, double t) const
212 {
213  Vector d = _FFD(x, y, z, t);
214  x = d._x, y = d._y, z = d._z;
215 }
216 
217 // -----------------------------------------------------------------------------
219 ::EvaluateInside(double &x, double &y, double &z, double t) const
220 {
221  Vector d = _FFD.GetInside(x, y, z, t);
222  x = d._x, y = d._y, z = d._z;
223 }
224 
225 // =============================================================================
226 // Point transformation
227 // =============================================================================
228 
229 // -----------------------------------------------------------------------------
231 ::LocalTransform(double &x, double &y, double &z, double t, double) const
232 {
233  // Convert to lattice coordinates
234  double dx = x, dy = y, dz = z;
235  this->WorldToLattice(dx, dy, dz);
236  t = this->TimeToLattice(t);
237  // Evaluate displacement
238  this->Evaluate(dx, dy, dz, t);
239  // Transform point
240  x += dx, y += dy, z += dz;
241 }
242 
243 // -----------------------------------------------------------------------------
245 ::LocalInverse(double &x, double &y, double &z, double t, double) const
246 {
247  // Convert to lattice coordinates
248  double dx = x, dy = y, dz = z;
249  this->WorldToLattice(dx, dy, dz);
250  t = this->TimeToLattice(t);
251  // Evaluate displacement
252  this->Evaluate(dx, dy, dz, t);
253  // Transform point
254  x -= dx, y -= dy, z -= dz;
255  return true;
256 }
257 
258 // ===========================================================================
259 // Derivatives
260 // ===========================================================================
261 
262 // ---------------------------------------------------------------------------
264 ::LocalJacobian(Matrix &jac, double x, double y, double z, double t, double) const
265 {
266  // Convert to lattice coordinates
267  this->WorldToLattice(x, y, z);
268  t = this->TimeToLattice(t);
269  // Compute 1st order derivatives
270  this->EvaluateJacobian(jac, x, y, z, t);
271  // Convert derivatives to world coordinates
272  JacobianToWorld(jac);
273  // Add derivatives of "x" term in T(x) = x + this->Evaluate(x)
274  jac(0, 0) += 1.0;
275  jac(1, 1) += 1.0;
276  jac(2, 2) += 1.0;
277 }
278 
279 // -----------------------------------------------------------------------------
281 ::JacobianDOFs(double jac[3], int ci, int cj, int ck, int cl, double x, double y, double z, double t) const
282 {
283  // Convert point to lattice coordinates
284  this->WorldToLattice(x, y, z);
285  t = this->TimeToLattice(t);
286  // Calculate tensor product of 1D linear interpolation kernels
287  const double dx = abs(ci - x);
288  const double dy = abs(cj - y);
289  const double dz = abs(ck - z);
290  const double dt = abs(cl - t);
291  if (dx < 1.0 && dy < 1.0 && dz < 1.0 && dt < 1.0) {
292  jac[0] = jac[1] = jac[2] = (1.0 - dx) * (1.0 - dy) * (1.0 - dz) * (1.0 - dt);
293  } else {
294  jac[0] = jac[1] = jac[3] = .0;
295  }
296 }
297 
298 
299 } // namespace mirtk
300 
301 #endif // MIRTK_LinearFreeFormTransformation4D_H
virtual void BendingEnergyGradient(double *, double=NaN, bool=false, bool=true, bool=true) const
double BendingEnergy(int i, int j, int k, int l) const
Calculate the bending energy of the transformation at control points.
void Evaluate(double &, double &, double &, double) const
Evaluates the free-form transformation at a point in lattice coordinates.
virtual double TimeToLattice(double) const
Transforms time (in ms) to temporal lattice coordinate.
Interpolator _FFD
Interpolates control point values at arbitrary lattice locations.
void JacobianToWorld(double &, double &) const
virtual void Print(ostream &, Indent=0) const
Prints the parameters of the transformation.
virtual void ApproximateDOFs(const double *, const double *, const double *, const double *, const double *, const double *, const double *, int)
virtual double BendingEnergy(double, double, double, double=0, double=NaN, bool=true) const
Calculates the bending of the transformation.
virtual void Print(ostream &, Indent=0) const
Prints the parameters of the transformation.
virtual void ApproximateDOFsGradient(const double *, const double *, const double *, const double *, const double *, const double *, const double *, int, double *, double=1.0) const
virtual void LocalJacobian(Matrix &, double, double, double, double, double=NaN) const
Calculates the Jacobian of the local transformation w.r.t world coordinates.
void EvaluateJacobian(Matrix &, double, double, double, double) const
Calculates the Jacobian of the local transformation w.r.t world coordinates.
void EvaluateInside(double &, double &, double &, double) const
Evaluates the free-form transformation at a point in lattice coordinates.
virtual void LocalTransform(double &, double &, double &, double, double=1) const
Transforms a single point using the local transformation only.
virtual bool LocalInverse(double &, double &, double &, double, double=1) const
Transforms a single point using the inverse of the local transformation only.
Definition: IOConfig.h:41
virtual VoxelType GetInside(double, double, double, double) const
virtual ~LinearFreeFormTransformation4D()
Destructor.
void BoundingBox(double &, double &) const
Gets the temporal bounding box of the free-form deformation (in ms)
virtual void JacobianDOFs(Matrix &, int, int, int, int, double, double, double, double, double=NaN) const
Calculates the Jacobian of the transformation w.r.t a control point.
T _z
The z component.
Definition: Vector3D.h:60
T _y
The y component.
Definition: Vector3D.h:59
LinearFreeFormTransformation4D()
Default constructor.
virtual bool CanRead(TransformationType) const
Whether this transformation can read a file of specified type (i.e. format)
virtual void BoundingBox(int, double &, double &, double=1) const
T _x
The x component.
Definition: Vector3D.h:58
virtual void LocalJacobian(Matrix &, double, double, double, double=0, double=NaN) const
Calculates the Jacobian of the local transformation w.r.t world coordinates.
virtual void Interpolate(const double *, const double *, const double *)
virtual void JacobianDOFs(double [3], int, int, int, int, double, double, double, double) const
Calculates the Jacobian of the transformation w.r.t the transformation parameters.
virtual int KernelSize() const
Size of support region of the used kernel.
virtual void WorldToLattice(double &, double &) const
Transforms world coordinates (in mm) to lattice coordinates.