ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkTimeVaryingBSplineVelocityFieldTransform.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 __itkTimeVaryingBSplineVelocityFieldTransform_h
19 #define __itkTimeVaryingBSplineVelocityFieldTransform_h
20 
22 
23 namespace itk
24 {
25 
67 template<typename TScalar, unsigned int NDimensions>
69  public VelocityFieldTransform<TScalar, NDimensions>
70 {
71 public:
72 
78 
81 
83  itkNewMacro( Self );
84 
87 
91 
95 
98 
103 
106 
108  itkStaticConstMacro( Dimension, unsigned int, NDimensions );
109 
111  itkStaticConstMacro( VelocityFieldDimension, unsigned int, NDimensions + 1 );
112 
113  typedef typename VelocityFieldType::PointType VelocityFieldPointType;
114  typedef typename VelocityFieldType::SizeType VelocityFieldSizeType;
115  typedef typename VelocityFieldType::SpacingType VelocityFieldSpacingType;
116  typedef typename VelocityFieldType::DirectionType VelocityFieldDirectionType;
117 
119  typedef typename VelocityFieldType::Pointer TimeVaryingVelocityFieldControlPointLatticePointer;
120 
126 
128  typename VelocityFieldType::Pointer GetTimeVaryingVelocityFieldControlPointLattice()
129  {
130  return this->GetModifiableVelocityField();
131  }
132 
135  {
136  this->SetVelocityField( fieldLattice );
137  }
138 
146  virtual void UpdateTransformParameters( const DerivativeType & update, ScalarType factor = 1.0 );
147 
149  virtual void IntegrateVelocityField();
150 
152  itkSetMacro( VelocityFieldOrigin, VelocityFieldPointType );
153  itkGetConstMacro( VelocityFieldOrigin, VelocityFieldPointType );
155 
157  itkSetMacro( VelocityFieldSpacing, VelocityFieldSpacingType );
158  itkGetConstMacro( VelocityFieldSpacing, VelocityFieldSpacingType );
160 
162  itkSetMacro( VelocityFieldSize, VelocityFieldSizeType );
163  itkGetConstMacro( VelocityFieldSize, VelocityFieldSizeType );
165 
167  itkSetMacro( VelocityFieldDirection, VelocityFieldDirectionType );
168  itkGetConstMacro( VelocityFieldDirection, VelocityFieldDirectionType );
170 
172  itkSetMacro( SplineOrder, unsigned int );
173  itkGetConstMacro( SplineOrder, unsigned int );
175 
176 protected:
179  void PrintSelf( std::ostream& os, Indent indent ) const;
180 
181 private:
182  TimeVaryingBSplineVelocityFieldTransform( const Self& ); //purposely not implementen
183  void operator=( const Self& ); //purposely not implemented
184 
185  unsigned int m_SplineOrder;
187 
192 };
193 
194 } // end namespace itk
195 
196 #ifndef ITK_MANUAL_INSTANTIATION
197 # include "itkTimeVaryingBSplineVelocityFieldTransform.hxx"
198 #endif
199 
200 #endif // __itkTimeVaryingBSplineVelocityFieldTransform_h
void PrintSelf(std::ostream &os, Indent indent) const
Light weight base class for most itk classes.
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
Superclass::NumberOfParametersType NumberOfParametersType
virtual void SetVelocityField(VelocityFieldType *)
virtual void UpdateTransformParameters(const DerivativeType &update, ScalarType factor=1.0)
Integrate a time-varying velocity field represented by a B-spline control point lattice.
TPixel PixelType
Definition: itkImage.h:89
Image< OutputVectorType, VelocityFieldDimension > VelocityFieldType
VectorInterpolateImageFunction< VelocityFieldType, ScalarType > VelocityFieldInterpolatorType
Superclass::DerivativeType DerivativeType
Superclass::ParametersType ParametersType
Superclass::InverseTransformBasePointer InverseTransformBasePointer
Control indentation during Print() invocation.
Definition: itkIndent.h:49
VectorInterpolateImageFunction< DisplacementFieldType, ScalarType > InterpolatorType
virtual VelocityFieldType * GetModifiableVelocityField()
Superclass::DisplacementFieldType DisplacementFieldType
virtual void SetTimeVaryingVelocityFieldControlPointLattice(VelocityFieldType *fieldLattice)
Provides local/dense/high-dimensionality transformation via a a velocity field.