ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkDisplacementFieldToBSplineImageFilter.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 __itkDisplacementFieldToBSplineImageFilter_h
19 #define __itkDisplacementFieldToBSplineImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
24 #include "itkPointSet.h"
25 #include "itkVector.h"
26 
27 namespace itk
28 {
29 
40 template <typename TInputImage, typename TOutputImage = TInputImage>
42  : public ImageToImageFilter<TInputImage, TOutputImage>
43 {
44 public:
49 
51  itkNewMacro( Self );
52 
54  itkStaticConstMacro( ImageDimension, unsigned int, TInputImage::ImageDimension );
55 
56  typedef TInputImage InputFieldType;
57  typedef TOutputImage OutputFieldType;
58 
61 
63  typedef typename OutputFieldType::PixelType PixelType;
64  typedef typename OutputFieldType::PixelType VectorType;
65  typedef typename OutputFieldType::RegionType RegionType;
66  typedef typename OutputFieldType::IndexType IndexType;
67 
68  typedef typename OutputFieldType::PointType PointType;
69  typedef typename OutputFieldType::SpacingType SpacingType;
70  typedef typename OutputFieldType::PointType OriginType;
71  typedef typename OutputFieldType::SizeType SizeType;
72  typedef typename OutputFieldType::DirectionType DirectionType;
73 
74  typedef typename VectorType::RealValueType RealType;
76 
79 
86 
88  void SetDisplacementField( const InputFieldType * field )
89  {
90  this->SetInput( field );
91  }
92 
97  {
98  return this->GetInput( 0 );
99  }
100 
105  {
107  }
108 
114  void SetConfidenceImage( const RealImageType *image )
115  {
116  this->SetNthInput( 1, const_cast<RealImageType *>( image ) );
117  }
118  void SetInput1( const RealImageType *image ) { this->SetConfidenceImage( image ); }
120 
125  {
126  return static_cast<const RealImageType*>( this->ProcessObject::GetInput( 1 ) );
127  }
128 
132  itkSetMacro( SplineOrder, unsigned int );
133 
137  itkGetConstMacro( SplineOrder, unsigned int );
138 
146  itkSetMacro( NumberOfControlPoints, ArrayType );
147 
155  itkGetConstMacro( NumberOfControlPoints, ArrayType );
156 
163  itkSetMacro( NumberOfFittingLevels, ArrayType );
164 
171  void SetNumberOfFittingLevels( unsigned int n )
172  {
173  ArrayType nlevels;
174 
175  nlevels.Fill( n );
176  this->SetNumberOfFittingLevels( nlevels );
177  }
178 
185  itkGetConstMacro( NumberOfFittingLevels, ArrayType );
186 
190  itkBooleanMacro( EstimateInverse );
191  itkSetMacro( EstimateInverse, bool );
192  itkGetConstMacro( EstimateInverse, bool );
194 
198  itkBooleanMacro( EnforceStationaryBoundary );
199  itkSetMacro( EnforceStationaryBoundary, bool );
200  itkGetConstMacro( EnforceStationaryBoundary, bool );
202 
203 protected:
204 
207 
210 
212  void PrintSelf( std::ostream& os, Indent indent ) const;
213 
215  void GenerateData();
216 
217 private:
218  DisplacementFieldToBSplineImageFilter( const Self& ); //purposely not implemented
219  void operator=( const Self& ); //purposely not implemented
220 
223  unsigned int m_SplineOrder;
226 
228 
229 
230 };
231 
232 } // end namespace itk
233 
234 #ifndef ITK_MANUAL_INSTANTIATION
235 #include "itkDisplacementFieldToBSplineImageFilter.hxx"
236 #endif
237 
238 #endif
Base class for all process objects that output image data.
void Fill(const ValueType &)
virtual void SetNumberOfFittingLevels(ArrayType _arg)
DisplacementFieldControlPointLatticeType::Pointer m_DisplacementFieldControlPointLattice
const DisplacementFieldControlPointLatticeType * GetDisplacementFieldControlPointLattice() const
BSplineScatteredDataPointSetToImageFilter< PointSetType, OutputFieldType > BSplineFilterType
virtual void SetInput(const InputImageType *image)
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
const InputImageType * GetInput(void) const
void PrintSelf(std::ostream &os, Indent indent) const
BSplineFilterType::PointDataImageType DisplacementFieldControlPointLatticeType
Image filter which provides a B-spline output approximation.
DataObject * GetInput(const DataObjectIdentifierType &key)
Base class for filters that take an image as input and produce an image as output.
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Class which takes a displacement field image and smooths it using B-splines. The inverse can also be ...
virtual void SetNthInput(DataObjectPointerArraySizeType num, DataObject *input)
Templated n-dimensional image class.
Definition: itkImage.h:75