go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
Common/Transforms/itkAdvancedTranslationTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright UMC Utrecht and contributors
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 /*=========================================================================
19 
20  Program: Insight Segmentation & Registration Toolkit
21  Module: $RCSfile: itkAdvancedTranslationTransform.h,v $
22  Language: C++
23  Date: $Date: 2007-07-15 16:38:25 $
24  Version: $Revision: 1.36 $
25 
26  Copyright (c) Insight Software Consortium. All rights reserved.
27  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
28 
29  This software is distributed WITHOUT ANY WARRANTY; without even
30  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
31  PURPOSE. See the above copyright notices for more information.
32 
33 =========================================================================*/
34 #ifndef __itkAdvancedTranslationTransform_h
35 #define __itkAdvancedTranslationTransform_h
36 
37 #include <iostream>
38 #include "itkAdvancedTransform.h"
39 #include "itkMacro.h"
40 #include "itkMatrix.h"
41 
42 namespace itk
43 {
44 
52 template<
53 class TScalarType = double, // Data type for scalars (float or double)
54 unsigned int NDimensions = 3 >
55 // Number of dimensions
56 class ITK_EXPORT AdvancedTranslationTransform :
57  public AdvancedTransform< TScalarType, NDimensions, NDimensions >
58 {
59 public:
60 
64  typedef SmartPointer< Self > Pointer;
65  typedef SmartPointer< const Self > ConstPointer;
66 
68  itkNewMacro( Self );
69 
72 
74  itkStaticConstMacro( SpaceDimension, unsigned int, NDimensions );
75  itkStaticConstMacro( ParametersDimension, unsigned int, NDimensions );
76 
79 
84  typedef typename Superclass::TransformCategoryEnum TransformCategoryEnum;
85 
88 
90  typedef Vector< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > InputVectorType;
91  typedef Vector< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > OutputVectorType;
92 
94  typedef CovariantVector< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > InputCovariantVectorType;
95  typedef CovariantVector< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > OutputCovariantVectorType;
96 
98  typedef vnl_vector_fixed< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > InputVnlVectorType;
99  typedef vnl_vector_fixed< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > OutputVnlVectorType;
100 
102  typedef Point< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > InputPointType;
103  typedef Point< TScalarType, itkGetStaticConstMacro( SpaceDimension ) > OutputPointType;
104 
106  typedef typename Superclass
107  ::NonZeroJacobianIndicesType NonZeroJacobianIndicesType;
109  typedef typename Superclass
110  ::JacobianOfSpatialJacobianType JacobianOfSpatialJacobianType;
112  typedef typename Superclass
113  ::JacobianOfSpatialHessianType JacobianOfSpatialHessianType;
115 
119  const OutputVectorType & GetOffset( void ) const
120  { return m_Offset; }
121 
124  void SetParameters( const ParametersType & parameters ) override;
125 
127  const ParametersType & GetParameters( void ) const override;
128 
132  void SetOffset( const OutputVectorType & offset )
133  { m_Offset = offset; return; }
134 
136  void Compose( const Self * other, bool pre = 0 );
137 
142  void Translate( const OutputVectorType & offset, bool pre = 0 );
143 
148  OutputPointType TransformPoint( const InputPointType & point ) const override;
149 
150  OutputVectorType TransformVector( const InputVectorType & vector ) const override;
151 
152  OutputVnlVectorType TransformVector( const InputVnlVectorType & vector ) const override;
153 
155  const InputCovariantVectorType & vector ) const override;
156 
160  inline InputPointType BackTransform( const OutputPointType & point ) const;
161 
162  inline InputVectorType BackTransform( const OutputVectorType & vector ) const;
163 
164  inline InputVnlVectorType BackTransform( const OutputVnlVectorType & vector ) const;
165 
167  const OutputCovariantVectorType & vector ) const;
168 
173  bool GetInverse( Self * inverse ) const;
174 
177  const InputPointType &,
178  JacobianType &,
179  NonZeroJacobianIndicesType & ) const override;
180 
183  const InputPointType &,
184  SpatialJacobianType & ) const override;
185 
188  const InputPointType &,
189  SpatialHessianType & ) const override;
190 
193  const InputPointType &,
195  NonZeroJacobianIndicesType & ) const override;
196 
199  const InputPointType &,
202  NonZeroJacobianIndicesType & ) const override;
203 
206  const InputPointType &,
208  NonZeroJacobianIndicesType & ) const override;
209 
214  const InputPointType & ipp,
215  SpatialHessianType & sh,
217  NonZeroJacobianIndicesType & nonZeroJacobianIndices ) const override;
218 
220  void SetIdentity( void );
221 
224  { return NDimensions; }
225 
231  bool IsLinear() const override { return true; }
232 
237  {
238  return TransformCategoryEnum::Linear;
239  }
240 
241 
245  void SetFixedParameters( const FixedParametersType & ) override
246  { /* purposely blank */ }
247 
251  const FixedParametersType & GetFixedParameters( void ) const override
252  {
253  this->m_FixedParameters.SetSize( 0 );
254  return this->m_FixedParameters;
255  }
256 
257 
258 protected:
259 
263  void PrintSelf( std::ostream & os, Indent indent ) const override;
264 
265 private:
266 
267  AdvancedTranslationTransform( const Self & ); //purposely not implemented
268  void operator=( const Self & ); //purposely not implemented
269 
270  OutputVectorType m_Offset; // Offset of the transformation
271 
278 
279 };
280 
281 //class AdvancedTranslationTransform
282 
283 // Back transform a point
284 template< class TScalarType, unsigned int NDimensions >
285 inline
288 {
289  return point - m_Offset;
290 }
291 
292 
293 // Back transform a vector
294 template< class TScalarType, unsigned int NDimensions >
295 inline
298 {
299  return vect;
300 }
301 
302 
303 // Back transform a vnl_vector
304 template< class TScalarType, unsigned int NDimensions >
305 inline
308 {
309  return vect;
310 }
311 
312 
313 // Back Transform a CovariantVector
314 template< class TScalarType, unsigned int NDimensions >
315 inline
318 {
319  return vect;
320 }
321 
322 
323 } // namespace itk
324 
325 #ifndef ITK_MANUAL_INSTANTIATION
326 #include "itkAdvancedTranslationTransform.hxx"
327 #endif
328 
329 #endif /* __itkAdvancedTranslationTransform_h */
Transform maps points, vectors and covariant vectors from an input space to an output space.
Superclass::OutputVectorType OutputVectorType
Superclass::OutputVnlVectorType OutputVnlVectorType
Superclass::InputVectorType InputVectorType
Superclass::OutputPointType OutputPointType
Superclass ::OutputCovariantVectorType OutputCovariantVectorType
Superclass::InputPointType InputPointType
FixedArray< Matrix< ScalarType, InputSpaceDimension, InputSpaceDimension >, OutputSpaceDimension > SpatialHessianType
Superclass ::InputCovariantVectorType InputCovariantVectorType
Superclass::InputVnlVectorType InputVnlVectorType
Matrix< ScalarType, OutputSpaceDimension, InputSpaceDimension > SpatialJacobianType
Translation transformation of a vector space (e.g. space coordinates)
void GetSpatialJacobian(const InputPointType &, SpatialJacobianType &) const override
void GetSpatialHessian(const InputPointType &, SpatialHessianType &) const override
NumberOfParametersType GetNumberOfParameters(void) const override
vnl_vector_fixed< TScalarType, itkGetStaticConstMacro(SpaceDimension) > InputVnlVectorType
void SetFixedParameters(const FixedParametersType &) override
void SetParameters(const ParametersType &parameters) override
OutputVectorType TransformVector(const InputVectorType &vector) const override
InputCovariantVectorType BackTransform(const OutputCovariantVectorType &vector) const
itkStaticConstMacro(SpaceDimension, unsigned int, NDimensions)
CovariantVector< TScalarType, itkGetStaticConstMacro(SpaceDimension) > InputCovariantVectorType
AdvancedTransform< TScalarType, NDimensions, NDimensions > Superclass
Superclass ::NonZeroJacobianIndicesType NonZeroJacobianIndicesType
Superclass ::JacobianOfSpatialHessianType JacobianOfSpatialHessianType
OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &vector) const override
void GetJacobianOfSpatialHessian(const InputPointType &ipp, SpatialHessianType &sh, JacobianOfSpatialHessianType &jsh, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
void PrintSelf(std::ostream &os, Indent indent) const override
Vector< TScalarType, itkGetStaticConstMacro(SpaceDimension) > InputVectorType
vnl_vector_fixed< TScalarType, itkGetStaticConstMacro(SpaceDimension) > OutputVnlVectorType
CovariantVector< TScalarType, itkGetStaticConstMacro(SpaceDimension) > OutputCovariantVectorType
void GetJacobianOfSpatialHessian(const InputPointType &, JacobianOfSpatialHessianType &, NonZeroJacobianIndicesType &) const override
InputPointType BackTransform(const OutputPointType &point) const
InputVnlVectorType BackTransform(const OutputVnlVectorType &vector) const
void Compose(const Self *other, bool pre=0)
Point< TScalarType, itkGetStaticConstMacro(SpaceDimension) > InputPointType
Superclass ::JacobianOfSpatialJacobianType JacobianOfSpatialJacobianType
void Translate(const OutputVectorType &offset, bool pre=0)
OutputVnlVectorType TransformVector(const InputVnlVectorType &vector) const override
const ParametersType & GetParameters(void) const override
Point< TScalarType, itkGetStaticConstMacro(SpaceDimension) > OutputPointType
AdvancedTranslationTransform(const Self &)
void GetJacobian(const InputPointType &, JacobianType &, NonZeroJacobianIndicesType &) const override
void GetJacobianOfSpatialJacobian(const InputPointType &, JacobianOfSpatialJacobianType &, NonZeroJacobianIndicesType &) const override
InputVectorType BackTransform(const OutputVectorType &vector) const
Vector< TScalarType, itkGetStaticConstMacro(SpaceDimension) > OutputVectorType
const FixedParametersType & GetFixedParameters(void) const override
OutputPointType TransformPoint(const InputPointType &point) const override
bool GetInverse(Self *inverse) const
itkStaticConstMacro(ParametersDimension, unsigned int, NDimensions)
void GetJacobianOfSpatialJacobian(const InputPointType &, SpatialJacobianType &, JacobianOfSpatialJacobianType &, NonZeroJacobianIndicesType &) const override


Generated on OURCE_DATE_EPOCH for elastix by doxygen 1.9.1 elastix logo