ITK  5.2.0
Insight Toolkit
itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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 itkBSplineSmoothingOnUpdateDisplacementFieldTransform_h
19 #define itkBSplineSmoothingOnUpdateDisplacementFieldTransform_h
20 
22 
24 #include "itkPointSet.h"
25 
26 namespace itk
27 {
28 
51 template <typename TParametersValueType, unsigned int NDimensions>
53  : public DisplacementFieldTransform<TParametersValueType, NDimensions>
54 {
55 public:
56  ITK_DISALLOW_COPY_AND_MOVE(BSplineSmoothingOnUpdateDisplacementFieldTransform);
57 
63 
66 
68  itkNewMacro(Self);
69 
71  static constexpr unsigned int Dimension = NDimensions;
72 
74  using ScalarType = typename Superclass::ScalarType;
75  using DerivativeType = typename Superclass::DerivativeType;
76  using DerivativeValueType = typename DerivativeType::ValueType;
77  using DisplacementFieldType = typename Superclass::DisplacementFieldType;
78  using DisplacementFieldPointer = typename Superclass::DisplacementFieldPointer;
79  using DisplacementFieldConstPointer = typename Superclass::DisplacementFieldConstPointer;
80 
82 
87  using DisplacementVectorType = typename DisplacementFieldType::PixelType;
89  using SplineOrderType = unsigned int;
94  using ArrayValueType = typename ArrayType::ValueType;
95 
106  void
107  UpdateTransformParameters(const DerivativeType & update, ScalarType factor = 1.0) override;
108 
112  itkSetMacro(SplineOrder, SplineOrderType);
113 
117  itkGetConstMacro(SplineOrder, SplineOrderType);
118 
126  itkSetMacro(NumberOfControlPointsForTheUpdateField, ArrayType);
127 
135  itkGetConstMacro(NumberOfControlPointsForTheUpdateField, ArrayType);
136 
143  void
144  SetMeshSizeForTheUpdateField(const ArrayType &);
145 
153  itkSetMacro(NumberOfControlPointsForTheTotalField, ArrayType);
154 
162  itkGetConstMacro(NumberOfControlPointsForTheTotalField, ArrayType);
163 
170  void
171  SetMeshSizeForTheTotalField(const ArrayType &);
172 
176  itkBooleanMacro(EnforceStationaryBoundary);
177  itkSetMacro(EnforceStationaryBoundary, bool);
178  itkGetConstMacro(EnforceStationaryBoundary, bool);
180 
181 protected:
184 
185  void
186  PrintSelf(std::ostream & os, Indent indent) const override;
187 
189  typename LightObject::Pointer
190  InternalClone() const override;
191 
196  BSplineSmoothDisplacementField(const DisplacementFieldType *, const ArrayType &);
197 
198 private:
199  SplineOrderType m_SplineOrder{ 3 };
200  bool m_EnforceStationaryBoundary{ true };
203 };
204 
205 } // end namespace itk
206 
207 #ifndef ITK_MANUAL_INSTANTIATION
208 # include "itkBSplineSmoothingOnUpdateDisplacementFieldTransform.hxx"
209 #endif
210 
211 #endif // itkBSplineSmoothingOnUpdateDisplacementFieldTransform_h
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::SplineOrderType
unsigned int SplineOrderType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:89
itkDisplacementFieldTransform.h
itk::PointSet
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform
Representation of a smooth deformation field with B-splines.
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:52
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::DerivativeType
typename Superclass::DerivativeType DerivativeType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:75
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::ScalarType
typename Superclass::ScalarType ScalarType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:74
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::ArrayValueType
typename ArrayType::ValueType ArrayValueType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:94
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::m_NumberOfControlPointsForTheTotalField
ArrayType m_NumberOfControlPointsForTheTotalField
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:202
itk::DisplacementFieldToBSplineImageFilter
Class which takes a dense displacement field image and/or a set of points with associated displacemen...
Definition: itkDisplacementFieldToBSplineImageFilter.h:44
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::DisplacementFieldTransform
Provides local/dense/high-dimensionality transformation via a a displacement field.
Definition: itkDisplacementFieldTransform.h:86
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::DisplacementVectorType
typename DisplacementFieldType::PixelType DisplacementVectorType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:87
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::WeightsContainerType
typename BSplineFilterType::WeightsContainerType WeightsContainerType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:91
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::m_NumberOfControlPointsForTheUpdateField
ArrayType m_NumberOfControlPointsForTheUpdateField
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:201
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::DisplacementFieldControlPointLatticeType
DisplacementFieldType DisplacementFieldControlPointLatticeType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:92
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::DisplacementFieldPointer
typename Superclass::DisplacementFieldPointer DisplacementFieldPointer
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:78
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::DisplacementFieldConstPointer
typename Superclass::DisplacementFieldConstPointer DisplacementFieldConstPointer
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:79
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::TransformPointer
typename Transform< TParametersValueType, NDimensions, NDimensions >::Pointer TransformPointer
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:81
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::DisplacementFieldType
typename Superclass::DisplacementFieldType DisplacementFieldType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:77
itk::DisplacementFieldToBSplineImageFilter::WeightsContainerType
typename BSplineFilterType::WeightsContainerType WeightsContainerType
Definition: itkDisplacementFieldToBSplineImageFilter.h:90
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::ArrayType
typename BSplineFilterType::ArrayType ArrayType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:93
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itkDisplacementFieldToBSplineImageFilter.h
itk::BSplineSmoothingOnUpdateDisplacementFieldTransform::DerivativeValueType
typename DerivativeType::ValueType DerivativeValueType
Definition: itkBSplineSmoothingOnUpdateDisplacementFieldTransform.h:76
itk::DisplacementFieldToBSplineImageFilter::ArrayType
typename BSplineFilterType::ArrayType ArrayType
Definition: itkDisplacementFieldToBSplineImageFilter.h:92
itk::Transform
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
itkPointSet.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44