ITK  4.13.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 TParametersValueType, unsigned int NDimensions>
68 class ITK_TEMPLATE_EXPORT TimeVaryingBSplineVelocityFieldTransform :
69  public VelocityFieldTransform<TParametersValueType, NDimensions>
70 {
71 public:
72 
78 
81 
83  itkNewMacro( Self );
84 
86  typedef typename Superclass::InverseTransformBasePointer InverseTransformBasePointer;
87 
89  typedef typename Superclass::InterpolatorType InterpolatorType;
90  typedef typename Superclass::VelocityFieldInterpolatorType VelocityFieldIntegratorType;
91 
93  typedef typename Superclass::DisplacementFieldType DisplacementFieldType;
94  typedef typename Superclass::VelocityFieldType VelocityFieldType;
95 
97  typedef typename Superclass::ScalarType ScalarType;
98 
100  typedef typename Superclass::ParametersType ParametersType;
102  typedef typename Superclass::FixedParametersType FixedParametersType;
104  typedef typename Superclass::NumberOfParametersType NumberOfParametersType;
105 
107  typedef typename Superclass::DerivativeType DerivativeType;
108 
110  itkStaticConstMacro( Dimension, unsigned int, NDimensions );
111 
113  itkStaticConstMacro( VelocityFieldDimension, unsigned int, NDimensions + 1 );
114 
117  typedef typename VelocityFieldType::SpacingType VelocityFieldSpacingType;
119 
121  typedef typename VelocityFieldType::Pointer TimeVaryingVelocityFieldControlPointLatticePointer;
122 
128 
131  {
132  return this->GetModifiableVelocityField();
133  }
134 
137  {
138  this->SetVelocityField( fieldLattice );
139  }
140 
148  virtual void UpdateTransformParameters( const DerivativeType & update, ScalarType factor = 1.0 ) ITK_OVERRIDE;
149 
151  virtual void IntegrateVelocityField() ITK_OVERRIDE;
152 
154  itkSetMacro( VelocityFieldOrigin, VelocityFieldPointType );
155  itkGetConstMacro( VelocityFieldOrigin, VelocityFieldPointType );
157 
159  itkSetMacro( VelocityFieldSpacing, VelocityFieldSpacingType );
160  itkGetConstMacro( VelocityFieldSpacing, VelocityFieldSpacingType );
162 
164  itkSetMacro( VelocityFieldSize, VelocityFieldSizeType );
165  itkGetConstMacro( VelocityFieldSize, VelocityFieldSizeType );
167 
169  itkSetMacro( VelocityFieldDirection, VelocityFieldDirectionType );
170  itkGetConstMacro( VelocityFieldDirection, VelocityFieldDirectionType );
172 
174  itkSetMacro( SplineOrder, unsigned int );
175  itkGetConstMacro( SplineOrder, unsigned int );
177 
178 protected:
180  virtual ~TimeVaryingBSplineVelocityFieldTransform() ITK_OVERRIDE;
181  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
182 
183 private:
184  ITK_DISALLOW_COPY_AND_ASSIGN(TimeVaryingBSplineVelocityFieldTransform);
185 
186  unsigned int m_SplineOrder;
187  bool m_TemporalPeriodicity;
188 
189  VelocityFieldPointType m_VelocityFieldOrigin;
190  VelocityFieldSpacingType m_VelocityFieldSpacing;
191  VelocityFieldDirectionType m_VelocityFieldDirection;
192  VelocityFieldSizeType m_VelocityFieldSize;
193 };
194 
195 } // end namespace itk
196 
197 #ifndef ITK_MANUAL_INSTANTIATION
198 # include "itkTimeVaryingBSplineVelocityFieldTransform.hxx"
199 #endif
200 
201 #endif // itkTimeVaryingBSplineVelocityFieldTransform_h
Light weight base class for most itk classes.
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
virtual void SetTimeVaryingVelocityFieldControlPointLattice(VelocityFieldType *fieldLattice)
Integrate a time-varying velocity field represented by a B-spline control point lattice.
VelocityFieldTransform< TParametersValueType, NDimensions > Superclass
TPixel PixelType
Definition: itkImage.h:89
TParametersValueType ValueType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Provides local/dense/high-dimensionality transformation via a a velocity field.