ITK  6.0.0
Insight Toolkit
itkBSplineExponentialDiffeomorphicTransform.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  * https://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 itkBSplineExponentialDiffeomorphicTransform_h
19 #define itkBSplineExponentialDiffeomorphicTransform_h
20 
23 #include "itkPointSet.h"
24 
25 namespace itk
26 {
27 
55 template <typename TParametersValueType, unsigned int VDimension>
57  : public ConstantVelocityFieldTransform<TParametersValueType, VDimension>
58 {
59 public:
60  ITK_DISALLOW_COPY_AND_MOVE(BSplineExponentialDiffeomorphicTransform);
61 
67 
69  itkOverrideGetNameOfClassMacro(BSplineExponentialDiffeomorphicTransform);
70 
72  itkNewMacro(Self);
73 
75  static constexpr unsigned int ConstantVelocityFieldDimension = VDimension;
76 
78  static constexpr unsigned int Dimension = VDimension;
79 
81  using typename Superclass::ScalarType;
82  using typename Superclass::DerivativeType;
84 
85  using typename Superclass::ParametersType;
86  using typename Superclass::ParametersValueType;
87  using typename Superclass::FixedParametersType;
88  using typename Superclass::FixedParametersValueType;
89 
90  using typename Superclass::DisplacementFieldType;
94 
96 
102  using SplineOrderType = unsigned int;
106  using ArrayValueType = typename ArrayType::ValueType;
107 
113  void
114  UpdateTransformParameters(const DerivativeType & update, ScalarType factor = 1.0) override;
115 
120  BSplineSmoothConstantVelocityField(const ConstantVelocityFieldType *, const ArrayType &);
121 
125  itkSetMacro(SplineOrder, SplineOrderType);
126  itkGetConstMacro(SplineOrder, SplineOrderType);
136  itkSetMacro(NumberOfControlPointsForTheConstantVelocityField, ArrayType);
137  itkGetConstMacro(NumberOfControlPointsForTheConstantVelocityField, ArrayType);
147  itkSetMacro(NumberOfControlPointsForTheUpdateField, ArrayType);
148  itkGetConstMacro(NumberOfControlPointsForTheUpdateField, ArrayType);
157  void
158  SetMeshSizeForTheConstantVelocityField(const ArrayType &);
159 
166  void
167  SetMeshSizeForTheUpdateField(const ArrayType &);
168 
169 protected:
171  ~BSplineExponentialDiffeomorphicTransform() override = default;
172 
173  void
174  PrintSelf(std::ostream &, Indent) const override;
175 
176 private:
177  ArrayType m_NumberOfControlPointsForTheConstantVelocityField{};
178  ArrayType m_NumberOfControlPointsForTheUpdateField{};
179 
180  SplineOrderType m_SplineOrder{ 3 };
181 };
182 
183 } // end namespace itk
184 
185 #ifndef ITK_MANUAL_INSTANTIATION
186 # include "itkBSplineExponentialDiffeomorphicTransform.hxx"
187 #endif
188 
189 #endif // itkBSplineExponentialDiffeomorphicTransform_h
itk::BSplineExponentialDiffeomorphicTransform::WeightsContainerType
typename BSplineFilterType::WeightsContainerType WeightsContainerType
Definition: itkBSplineExponentialDiffeomorphicTransform.h:104
itk::PointSet
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:81
itk::BSplineExponentialDiffeomorphicTransform::SplineOrderType
unsigned int SplineOrderType
Definition: itkBSplineExponentialDiffeomorphicTransform.h:102
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::BSplineExponentialDiffeomorphicTransform::DisplacementVectorType
typename DisplacementFieldType::PixelType DisplacementVectorType
Definition: itkBSplineExponentialDiffeomorphicTransform.h:95
itk::BSplineExponentialDiffeomorphicTransform::DerivativeValueType
typename DerivativeType::ValueType DerivativeValueType
Definition: itkBSplineExponentialDiffeomorphicTransform.h:83
itkConstantVelocityFieldTransform.h
itk::BSplineExponentialDiffeomorphicTransform
Exponential transform using B-splines as the smoothing kernel.
Definition: itkBSplineExponentialDiffeomorphicTransform.h:56
itk::Transform< TParametersValueType, VDimension, VDimension >::ScalarType
ParametersValueType ScalarType
Definition: itkTransform.h:127
itk::BSplineExponentialDiffeomorphicTransform::ArrayValueType
typename ArrayType::ValueType ArrayValueType
Definition: itkBSplineExponentialDiffeomorphicTransform.h:106
itk::ConstantVelocityFieldTransform::ConstantVelocityFieldPointer
typename ConstantVelocityFieldType::Pointer ConstantVelocityFieldPointer
Definition: itkConstantVelocityFieldTransform.h:98
itk::Image::PixelType
TPixel PixelType
Definition: itkImage.h:108
itk::ConstantVelocityFieldTransform::DisplacementFieldPointer
typename DisplacementFieldType::Pointer DisplacementFieldPointer
Definition: itkConstantVelocityFieldTransform.h:94
itk::BSplineExponentialDiffeomorphicTransform::ArrayType
typename BSplineFilterType::ArrayType ArrayType
Definition: itkBSplineExponentialDiffeomorphicTransform.h:105
itk::DisplacementFieldToBSplineImageFilter::WeightsContainerType
typename BSplineFilterType::WeightsContainerType WeightsContainerType
Definition: itkDisplacementFieldToBSplineImageFilter.h:93
itk::Array::ValueType
TValue ValueType
Definition: itkArray.h:52
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itkDisplacementFieldToBSplineImageFilter.h
itk::Array
Array class with size defined at construction time.
Definition: itkArray.h:47
itk::DisplacementFieldToBSplineImageFilter::ArrayType
typename BSplineFilterType::ArrayType ArrayType
Definition: itkDisplacementFieldToBSplineImageFilter.h:95
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itkPointSet.h
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::ConstantVelocityFieldTransform
Provides local/dense/high-dimensionality transformation via a a constant velocity field.
Definition: itkConstantVelocityFieldTransform.h:36