ITK  4.3.0
Insight Segmentation and Registration Toolkit
itkAffineTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef __itkAffineTransform_h
19 #define __itkAffineTransform_h
20 
21 #include <iostream>
22 
24 
25 namespace itk
26 {
98 
99 template<
100  class TScalarType = double, // Data type for scalars
101  // (e.g. float or double)
102  unsigned int NDimensions = 3 >
103 // Number of dimensions in the input space
105  public MatrixOffsetTransformBase< TScalarType, NDimensions, NDimensions >
106 {
107 public:
110  typedef MatrixOffsetTransformBase< TScalarType,
111  NDimensions,
112  NDimensions > Superclass;
113 
116 
119 
121  itkNewMacro(Self);
122 
124  itkStaticConstMacro(InputSpaceDimension, unsigned int, NDimensions);
125  itkStaticConstMacro(OutputSpaceDimension, unsigned int, NDimensions);
126  itkStaticConstMacro(SpaceDimension, unsigned int, NDimensions);
127  itkStaticConstMacro( ParametersDimension, unsigned int,
128  NDimensions *( NDimensions + 1 ) );
130 
148 
152  typedef typename InverseTransformBaseType::Pointer InverseTransformBasePointer;
153 
160  void Translate(const OutputVectorType & offset, bool pre = 0);
161 
173  void Scale(const OutputVectorType & factor, bool pre = 0);
174 
175  void Scale(const TScalarType & factor, bool pre = 0);
176 
192  void Rotate(int axis1, int axis2, TScalarType angle, bool pre = 0);
194 
208  void Rotate2D(TScalarType angle, bool pre = 0);
209 
223  void Rotate3D(const OutputVectorType & axis, TScalarType angle, bool pre = 0);
224 
236  void Shear(int axis1, int axis2, TScalarType coef, bool pre = 0);
237 
239  bool GetInverse(Self *inverse) const;
240 
243 
252  itkLegacyMacro(InputPointType BackTransform(const OutputPointType & point) const);
253  itkLegacyMacro(InputVectorType BackTransform(const OutputVectorType & vector) const);
254  itkLegacyMacro(InputVnlVectorType BackTransform( const OutputVnlVectorType & vector) const);
255  itkLegacyMacro(InputCovariantVectorType BackTransform( const OutputCovariantVectorType & vector) const);
257 
267  itkLegacyMacro(InputPointType BackTransformPoint(const OutputPointType & point) const);
268 
280  ScalarType Metric(const Self *other) const;
281 
285  ScalarType Metric(void) const;
286 
287 protected:
295  AffineTransform(const MatrixType & matrix,
296  const OutputVectorType & offset);
297  AffineTransform(unsigned int paramDims);
298  AffineTransform();
300 
302  virtual ~AffineTransform();
303 
305  void PrintSelf(std::ostream & s, Indent indent) const;
306 
307 private:
308 
309  AffineTransform(const Self & other);
310  const Self & operator=(const Self &);
311 }; //class AffineTransform
312 
313 #if !defined(ITK_LEGACY_REMOVE)
314 
315 template< class TScalarType, unsigned int NDimensions >
316 inline
318 AffineTransform< TScalarType, NDimensions >::BackTransform(const OutputVectorType & vect) const
319 {
320  itkWarningMacro(
321  << "BackTransform(): This method is slated to be removed "
322  << "from ITK. Instead, please use GetInverse() to generate an inverse "
323  << "transform and then perform the transform using that inverted transform.");
324  return this->GetInverseMatrix() * vect;
325 }
327 
329 template< class TScalarType, unsigned int NDimensions >
330 inline
331 typename AffineTransform< TScalarType, NDimensions >::InputVnlVectorType
332 AffineTransform< TScalarType, NDimensions >::BackTransform(const OutputVnlVectorType & vect) const
333 {
334  itkWarningMacro(
335  << "BackTransform(): This method is slated to be removed "
336  << "from ITK. Instead, please use GetInverse() to generate an inverse "
337  << "transform and then perform the transform using that inverted transform.");
338  return this->GetInverseMatrix() * vect;
339 }
341 
343 template< class TScalarType, unsigned int NDimensions >
344 inline
345 typename AffineTransform< TScalarType, NDimensions >::InputCovariantVectorType
346 AffineTransform< TScalarType, NDimensions >::BackTransform(const OutputCovariantVectorType & vec) const
347 {
348  itkWarningMacro(
349  << "BackTransform(): This method is slated to be removed "
350  << "from ITK. Instead, please use GetInverse() to generate an inverse "
351  << "transform and then perform the transform using that inverted transform.");
352 
353  InputCovariantVectorType result; // Converted vector
354 
355  for ( unsigned int i = 0; i < NDimensions; i++ )
356  {
358  for ( unsigned int j = 0; j < NDimensions; j++ )
359  {
360  result[i] += this->GetMatrix()[j][i] * vec[j]; // Direct matrix transposed
361  }
362  }
363  return result;
364 }
365 
367 template< class TScalarType, unsigned int NDimensions >
368 inline
369 typename AffineTransform< TScalarType, NDimensions >::InputPointType
370 AffineTransform< TScalarType, NDimensions >::BackTransformPoint(const OutputPointType & point) const
371 {
372  return this->BackTransform(point);
373 }
374 
376 template< class TScalarType, unsigned int NDimensions >
377 inline
378 typename AffineTransform< TScalarType, NDimensions >::InputPointType
379 AffineTransform< TScalarType, NDimensions >::BackTransform(const OutputPointType & point) const
380 {
381  itkWarningMacro(
382  << "BackTransform(): This method is slated to be removed "
383  << "from ITK. Instead, please use GetInverse() to generate an inverse "
384  << "transform and then perform the transform using that inverted transform.");
385  InputPointType result; // Converted point
386  ScalarType temp[NDimensions];
387  unsigned int i, j;
389 
390  for ( j = 0; j < NDimensions; j++ )
391  {
392  temp[j] = point[j] - this->GetOffset()[j];
393  }
394 
395  for ( i = 0; i < NDimensions; i++ )
396  {
397  result[i] = 0.0;
398  for ( j = 0; j < NDimensions; j++ )
399  {
400  result[i] += this->GetInverseMatrix()[i][j] * temp[j];
401  }
402  }
403  return result;
404 }
405 #endif
406 } // namespace itk
407 
408 #ifndef ITK_MANUAL_INSTANTIATION
409 #include "itkAffineTransform.hxx"
410 #endif
411 
412 #endif /* __itkAffineTransform_h */
413