ITK  5.4.0
Insight Toolkit
itkBSplineBaseTransform.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 itkBSplineBaseTransform_h
19#define itkBSplineBaseTransform_h
20
21#include <iostream>
22#include "itkTransform.h"
23#include "itkImage.h"
25
26namespace itk
27{
33template <typename TParametersValueType = double, unsigned int VDimension = 3, unsigned int VSplineOrder = 3>
34class ITK_TEMPLATE_EXPORT BSplineBaseTransform : public Transform<TParametersValueType, VDimension, VDimension>
35{
36public:
37 ITK_DISALLOW_COPY_AND_MOVE(BSplineBaseTransform);
38
44
46 itkOverrideGetNameOfClassMacro(BSplineBaseTransform);
47
49 static constexpr unsigned int SpaceDimension = VDimension;
50
52 static constexpr unsigned int SplineOrder = VSplineOrder;
53
56
58 using typename Superclass::ScalarType;
59
61 using typename Superclass::FixedParametersType;
62 using typename Superclass::ParametersType;
63
65 using typename Superclass::JacobianType;
66 using typename Superclass::JacobianPositionType;
67 using typename Superclass::InverseJacobianPositionType;
68
70 using typename Superclass::TransformCategoryEnum;
71
73 using typename Superclass::NumberOfParametersType;
74
78
82
84 using InputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;
85 using OutputVnlVectorType = vnl_vector_fixed<TParametersValueType, SpaceDimension>;
86
90
110 void
111 SetParameters(const ParametersType & parameters) override;
112
135 void
136 SetFixedParameters(const FixedParametersType & parameters) override = 0;
155 void
156 SetParametersByValue(const ParametersType & parameters) override;
157
166 void
168
170 const ParametersType &
171 GetParameters() const override;
172
174 const FixedParametersType &
175 GetFixedParameters() const override;
176
182
194 virtual void
196
200 {
201 return this->m_CoefficientImages;
202 }
203
204 using typename Superclass::DerivativeType;
205
216 void
217 UpdateTransformParameters(const DerivativeType & update, TParametersValueType factor = 1.0) override;
218
221
227
230 TransformPoint(const InputPointType & point) const override;
231
234
237
239 static constexpr unsigned int NumberOfWeights{ WeightsFunctionType::NumberOfWeights };
240
243
252 virtual void
253 TransformPoint(const InputPointType & inputPoint,
254 OutputPointType & outputPoint,
255 WeightsType & weights,
256 ParameterIndexArrayType & indices,
257 bool & inside) const = 0;
258
259#if !defined(ITK_LEGACY_REMOVE)
261 itkLegacyMacro(unsigned long GetNumberOfWeights() const) { return m_WeightsFunction->GetNumberOfWeights(); }
262#endif
263
266 using Superclass::TransformVector;
267 OutputVectorType
268 TransformVector(const InputVectorType &) const override
269 {
270 itkExceptionMacro("Method not applicable for deformable transform.");
271 }
276 OutputVnlVectorType
277 TransformVector(const InputVnlVectorType &) const override
278 {
279 itkExceptionMacro("Method not applicable for deformable transform. ");
280 }
281
284 using Superclass::TransformCovariantVector;
285 OutputCovariantVectorType
287 {
288 itkExceptionMacro("Method not applicable for deformable transform. ");
289 }
293 void
295 WeightsType &,
297
298 void
300
301 void
303 {
304 itkExceptionMacro("ComputeJacobianWithRespectToPosition not yet implemented "
305 "for "
306 << this->GetNameOfClass());
307 }
308 using Superclass::ComputeJacobianWithRespectToPosition;
309
311 NumberOfParametersType
312 GetNumberOfParameters() const override = 0;
313
317
319 GetTransformCategory() const override
320 {
321 return Self::TransformCategoryEnum::BSpline;
322 }
323
324 unsigned int
326
329
331
335 {
336 return this->GetNumberOfParameters();
337 }
338
339protected:
341 void
342 PrintSelf(std::ostream & os, Indent indent) const override;
343
345 ~BSplineBaseTransform() override = default;
346
348 itkSetObjectMacro(WeightsFunction, WeightsFunctionType);
349 itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType);
353 void
355
356protected:
358 void
360
362 virtual void
364
366 virtual void
368
370 virtual void
372
374 virtual void
376
378 virtual void
380
382 virtual bool
384
385 // NOTE: There is a natural duality between the
386 // two representations of of the coefficients
387 // whereby the m_InternalParametersBuffer is
388 // needed to fit into the optimization framework
389 // and the m_CoefficientImages is needed for
390 // the Jacobian computations. This implementation
391 // is an attempt to remove as much redundancy as possible
392 // and share as much information between the two
393 // instances as possible.
394 //
399 CoefficientImageArray m_CoefficientImages{ Self::ArrayOfImagePointerGeneratorHelper() };
400
402 ParametersType m_InternalParametersBuffer{};
403
406
407private:
408 static CoefficientImageArray
410}; // class BSplineBaseTransform
411} // namespace itk
412
413#ifndef ITK_MANUAL_INSTANTIATION
414# include "itkBSplineBaseTransform.hxx"
415#endif
416
417#endif /* itkBSplineBaseTransform_h */
Array2D class representing a 2D array.
Definition: itkArray2D.h:43
Array class with size defined at construction time.
Definition: itkArray.h:48
A base class with common elements of BSplineTransform and BSplineDeformableTransform.
void SetParameters(const ParametersType &parameters) override
typename ImageType::PixelType PixelType
virtual void SetFixedParametersGridDirectionFromTransformDomainInformation() const =0
typename ImageType::Pointer ImagePointer
typename ImageType::PointType OriginType
typename WeightsFunctionType::ContinuousIndexType ContinuousIndexType
virtual NumberOfParametersType GetNumberOfParametersPerDimension() const =0
virtual void TransformPoint(const InputPointType &inputPoint, OutputPointType &outputPoint, WeightsType &weights, ParameterIndexArrayType &indices, bool &inside) const =0
NumberOfParametersType GetNumberOfLocalParameters() const override
OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const override
void SetFixedParametersFromTransformDomainInformation() const
virtual bool InsideValidRegion(ContinuousIndexType &) const =0
const FixedParametersType & GetFixedParameters() const override
vnl_vector_fixed< TParametersValueType, SpaceDimension > OutputVnlVectorType
virtual void SetFixedParametersGridSizeFromTransformDomainInformation() const =0
const ParametersType & GetParameters() const override
typename ImageType::SpacingType SpacingType
const CoefficientImageArray GetCoefficientImages() const
virtual void SetFixedParametersGridOriginFromTransformDomainInformation() const =0
NumberOfParametersType GetNumberOfParameters() const override=0
virtual void SetCoefficientImages(const CoefficientImageArray &images)=0
void ComputeJacobianWithRespectToPosition(const InputPointType &, JacobianPositionType &) const override
typename RegionType::SizeType SizeType
vnl_vector_fixed< TParametersValueType, SpaceDimension > InputVnlVectorType
void ComputeJacobianFromBSplineWeightsWithRespectToPosition(const InputPointType &, WeightsType &, ParameterIndexArrayType &) const
void PrintSelf(std::ostream &os, Indent indent) const override
TransformCategoryEnum GetTransformCategory() const override
void SetParametersByValue(const ParametersType &parameters) override
typename ImageType::SpacingType PhysicalDimensionsType
virtual void SetCoefficientImageInformationFromFixedParameters()=0
static CoefficientImageArray ArrayOfImagePointerGeneratorHelper()
virtual void SetFixedParametersGridSpacingFromTransformDomainInformation() const =0
OutputVectorType TransformVector(const InputVectorType &) const override
typename RegionType::IndexType IndexType
void ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &) const override=0
OutputVnlVectorType TransformVector(const InputVnlVectorType &) const override
void UpdateTransformParameters(const DerivativeType &update, TParametersValueType factor=1.0) override
OutputPointType TransformPoint(const InputPointType &point) const override
typename WeightsFunctionType::WeightsType WeightsType
unsigned int GetNumberOfAffectedWeights() const
~BSplineBaseTransform() override=default
typename ParametersType::ValueType ParametersValueType
void SetFixedParameters(const FixedParametersType &parameters) override=0
typename ImageType::DirectionType DirectionType
Returns the weights over the support region used for B-spline interpolation/reconstruction.
A templated class holding a point in n-Dimensional image space.
A templated class holding a n-Dimensional covariant vector.
An image region represents a structured region of data.
Templated n-dimensional image class.
Definition: itkImage.h:89
TPixel PixelType
Definition: itkImage.h:108
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Class to hold and manage different parameter types used during optimization.
TParametersValueType ValueType
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:54
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:84
vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension > JacobianPositionType
Definition: itkTransform.h:145
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:63
const char * GetNameOfClass() const override
static Pointer New()
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar images