go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkRecursiveBSplineTransform.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright UMC Utrecht and contributors
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 itkRecursiveBSplineTransform_h
19#define itkRecursiveBSplineTransform_h
20
22
25#include "elxDefaultConstruct.h"
26
27namespace itk
28{
38template <typename TScalarType = double, unsigned int NDimensions = 3, unsigned int VSplineOrder = 3>
39class ITK_TEMPLATE_EXPORT RecursiveBSplineTransform
40 : public AdvancedBSplineDeformableTransform<TScalarType, NDimensions, VSplineOrder>
41{
42public:
44
48 using Pointer = SmartPointer<Self>;
49 using ConstPointer = SmartPointer<const Self>;
50
52 itkNewMacro(Self);
53
56
58 itkStaticConstMacro(SpaceDimension, unsigned int, NDimensions);
59
61 itkStaticConstMacro(SplineOrder, unsigned int, VSplineOrder);
62
64 using typename Superclass::ScalarType;
65 using typename Superclass::ParametersType;
66 using typename Superclass::ParametersValueType;
67 using typename Superclass::NumberOfParametersType;
68 using typename Superclass::DerivativeType;
69 using typename Superclass::JacobianType;
70 using typename Superclass::InputVectorType;
71 using typename Superclass::OutputVectorType;
72 using typename Superclass::InputCovariantVectorType;
73 using typename Superclass::OutputCovariantVectorType;
74 using typename Superclass::InputVnlVectorType;
75 using typename Superclass::OutputVnlVectorType;
76 using typename Superclass::InputPointType;
77 using typename Superclass::OutputPointType;
78
80 using typename Superclass::PixelType;
81 using typename Superclass::ImageType;
82 using typename Superclass::ImagePointer;
83 // using typename Superclass::CoefficientImageArray;
84
86 using typename Superclass::RegionType;
87 using typename Superclass::IndexType;
88 using typename Superclass::SizeType;
89 using typename Superclass::SpacingType;
90 using typename Superclass::DirectionType;
91 using typename Superclass::OriginType;
92 using typename Superclass::GridOffsetType;
93 using OffsetValueType = typename GridOffsetType::OffsetValueType;
94
98 using typename Superclass::SpatialHessianType;
100 using typename Superclass::InternalMatrixType;
103
105 using typename Superclass::WeightsFunctionType;
107 using typename Superclass::WeightsType;
108 using typename Superclass::ContinuousIndexType;
113
116
121 OutputPointType
122 TransformPoint(const InputPointType & point) const override;
123
125 void
126 GetJacobian(const InputPointType & inputPoint,
127 JacobianType & j,
128 NonZeroJacobianIndicesType & nonZeroJacobianIndices) const override;
129
133 void
134 EvaluateJacobianWithImageGradientProduct(const InputPointType & inputPoint,
135 const MovingImageGradientType & movingImageGradient,
136 DerivativeType & imageJacobian,
137 NonZeroJacobianIndicesType & nonZeroJacobianIndices) const override;
138
140 void
141 GetSpatialJacobian(const InputPointType & inputPoint, SpatialJacobianType & sj) const override;
142
144 void
145 GetSpatialHessian(const InputPointType & inputPoint, SpatialHessianType & sh) const override;
146
148 void
149 GetJacobianOfSpatialJacobian(const InputPointType & inputPoint,
151 NonZeroJacobianIndicesType & nonZeroJacobianIndices) const override;
152
156 void
157 GetJacobianOfSpatialJacobian(const InputPointType & inputPoint,
160 NonZeroJacobianIndicesType & nonZeroJacobianIndices) const override;
161
163 void
164 GetJacobianOfSpatialHessian(const InputPointType & inputPoint,
166 NonZeroJacobianIndicesType & nonZeroJacobianIndices) const override;
167
171 void
172 GetJacobianOfSpatialHessian(const InputPointType & inputPoint,
175 NonZeroJacobianIndicesType & nonZeroJacobianIndices) const override;
176
177protected:
179 ~RecursiveBSplineTransform() override = default;
180
181 using typename Superclass::JacobianImageType;
182 using typename Superclass::JacobianPixelType;
183
185 void
187 const RegionType & supportRegion) const override;
188
189private:
192
195
197};
198
199} // end namespace itk
200
201#ifndef ITK_MANUAL_INSTANTIATION
202# include "itkRecursiveBSplineTransform.hxx"
203#endif
204
205#endif /* itkRecursiveBSplineTransform_h */
Deformable transform using a B-spline representation.
typename WeightsFunctionType::WeightsType WeightsType
typename DerivativeWeightsFunctionType::Pointer DerivativeWeightsFunctionPointer
typename SODerivativeWeightsFunctionType::Pointer SODerivativeWeightsFunctionPointer
typename WeightsFunctionType::Pointer WeightsFunctionPointer
typename WeightsFunctionType::ContinuousIndexType ContinuousIndexType
typename SpatialJacobianType::InternalMatrixType InternalMatrixType
FixedArray< Matrix< ScalarType, InputSpaceDimension, InputSpaceDimension >, OutputSpaceDimension > SpatialHessianType
Matrix< ScalarType, OutputSpaceDimension, InputSpaceDimension > SpatialJacobianType
typename MovingImageGradientType::ValueType MovingImageGradientValueType
Returns the weights over the support region used for B-spline interpolation/reconstruction.
Returns the weights over the support region used for B-spline interpolation/reconstruction.
Returns the weights over the support region used for B-spline interpolation/reconstruction.
Returns the weights over the support region used for B-spline interpolation/reconstruction.
This helper class contains the actual implementation of the recursive B-spline transform.
A recursive implementation of the B-spline transform.
void GetJacobian(const InputPointType &inputPoint, JacobianType &j, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
typename GridOffsetType::OffsetValueType OffsetValueType
~RecursiveBSplineTransform() override=default
void GetJacobianOfSpatialJacobian(const InputPointType &inputPoint, JacobianOfSpatialJacobianType &jsj, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
void GetJacobianOfSpatialHessian(const InputPointType &inputPoint, SpatialHessianType &sh, JacobianOfSpatialHessianType &jsh, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
itkStaticConstMacro(SplineOrder, unsigned int, VSplineOrder)
void GetJacobianOfSpatialJacobian(const InputPointType &inputPoint, SpatialJacobianType &sj, JacobianOfSpatialJacobianType &jsj, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
void ComputeNonZeroJacobianIndices(NonZeroJacobianIndicesType &nonZeroJacobianIndices, const RegionType &supportRegion) const override
void GetSpatialHessian(const InputPointType &inputPoint, SpatialHessianType &sh) const override
void EvaluateJacobianWithImageGradientProduct(const InputPointType &inputPoint, const MovingImageGradientType &movingImageGradient, DerivativeType &imageJacobian, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
void GetJacobianOfSpatialHessian(const InputPointType &inputPoint, JacobianOfSpatialHessianType &jsh, NonZeroJacobianIndicesType &nonZeroJacobianIndices) const override
void GetSpatialJacobian(const InputPointType &inputPoint, SpatialJacobianType &sj) const override
OutputPointType TransformPoint(const InputPointType &point) const override
itkStaticConstMacro(SpaceDimension, unsigned int, NDimensions)
ITK_DISALLOW_COPY_AND_MOVE(RecursiveBSplineTransform)


Generated on 2024-07-17 for elastix by doxygen 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) elastix logo