ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
62 using typename Superclass::ParametersType;
63
65 using typename Superclass::JacobianType;
68
71
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;
137
154 void
155 SetParametersByValue(const ParametersType & parameters) override;
156
165 void
167
169 const ParametersType &
170 GetParameters() const override;
171
173 const FixedParametersType &
174 GetFixedParameters() const override;
175
181
193 virtual void
195
199 {
200 return this->m_CoefficientImages;
201 }
202
203 using typename Superclass::DerivativeType;
204
215 void
216 UpdateTransformParameters(const DerivativeType & update, TParametersValueType factor = 1.0) override;
217
220
226
229 TransformPoint(const InputPointType & point) const override;
230
233
236
239
242
251 virtual void
252 TransformPoint(const InputPointType & inputPoint,
253 OutputPointType & outputPoint,
254 WeightsType & weights,
255 ParameterIndexArrayType & indices,
256 bool & inside) const = 0;
257
258#if !defined(ITK_LEGACY_REMOVE)
260 itkLegacyMacro(unsigned long GetNumberOfWeights() const)
261 {
262 return m_WeightsFunction->GetNumberOfWeights();
263 }
264#endif
265
268 using Superclass::TransformVector;
269 OutputVectorType
270 TransformVector(const InputVectorType &) const override
271 {
272 itkExceptionMacro("Method not applicable for deformable transform.");
273 }
274
277 OutputVnlVectorType
278 TransformVector(const InputVnlVectorType &) const override
279 {
280 itkExceptionMacro("Method not applicable for deformable transform. ");
281 }
282
285 using Superclass::TransformCovariantVector;
286 OutputCovariantVectorType
288 {
289 itkExceptionMacro("Method not applicable for deformable transform. ");
290 }
291
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
349 itkSetObjectMacro(WeightsFunction, WeightsFunctionType);
350 itkGetModifiableObjectMacro(WeightsFunction, WeightsFunctionType);
352
354 void
356
357protected:
359 void
361
363 virtual void
365
367 virtual void
369
371 virtual void
373
375 virtual void
377
379 virtual void
381
383 virtual bool
385
386 // NOTE: There is a natural duality between the
387 // two representations of of the coefficients
388 // whereby the m_InternalParametersBuffer is
389 // needed to fit into the optimization framework
390 // and the m_CoefficientImages is needed for
391 // the Jacobian computations. This implementation
392 // is an attempt to remove as much redundancy as possible
393 // and share as much information between the two
394 // instances as possible.
395 //
401
404
407
408private:
409 static CoefficientImageArray
411}; // class BSplineBaseTransform
412} // namespace itk
413
414#ifndef ITK_MANUAL_INSTANTIATION
415# include "itkBSplineBaseTransform.hxx"
416#endif
417
418#endif /* itkBSplineBaseTransform_h */
FixedArray< unsigned long, NumberOfWeights > ParameterIndexArrayType
void SetParameters(const ParametersType &parameters) override
typename ImageType::PixelType PixelType
virtual void SetFixedParametersGridDirectionFromTransformDomainInformation() const =0
Vector< TParametersValueType, Self::SpaceDimension > OutputVectorType
typename ImageType::Pointer ImagePointer
CovariantVector< TParametersValueType, Self::SpaceDimension > InputCovariantVectorType
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
Image< ParametersValueType, Self::SpaceDimension > ImageType
vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension > JacobianPositionType
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
Vector< TParametersValueType, Self::SpaceDimension > InputVectorType
Point< TParametersValueType, Self::SpaceDimension > InputPointType
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
SmartPointer< const Self > ConstPointer
void SetParametersByValue(const ParametersType &parameters) override
typename ImageType::SpacingType PhysicalDimensionsType
Transform< TParametersValueType, VDimension, VDimension > Superclass
virtual void SetCoefficientImageInformationFromFixedParameters()=0
static CoefficientImageArray ArrayOfImagePointerGeneratorHelper()
virtual void SetFixedParametersGridSpacingFromTransformDomainInformation() const =0
Point< TParametersValueType, Self::SpaceDimension > OutputPointType
OutputVectorType TransformVector(const InputVectorType &) const override
BSplineInterpolationWeightFunction< ScalarType, Self::SpaceDimension, Self::SplineOrder > WeightsFunctionType
typename RegionType::IndexType IndexType
void ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &) const override=0
FixedArray< ImagePointer, VDimension > CoefficientImageArray
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
ImageRegion< Self::SpaceDimension > RegionType
typename ParametersType::ValueType ParametersValueType
void SetFixedParameters(const FixedParametersType &parameters) override=0
CovariantVector< TParametersValueType, Self::SpaceDimension > OutputCovariantVectorType
typename ImageType::DirectionType DirectionType
Returns the weights over the support region used for B-spline interpolation/reconstruction.
A templated class holding a n-Dimensional covariant vector.
Simulate a standard C array with copy semantics.
An image region represents a structured region of data.
Templated n-dimensional image class.
Definition itkImage.h:89
Point< PointValueType, VImageDimension > PointType
Vector< SpacingValueType, VImageDimension > SpacingType
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
A templated class holding a geometric point in n-Dimensional space.
Definition itkPoint.h:54
Implements transparent reference counting.
OptimizerParameters< ParametersValueType > ParametersType
OptimizerParameters< FixedParametersValueType > FixedParametersType
TransformBaseTemplateEnums::TransformCategory TransformCategoryEnum
vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension > JacobianPositionType
vnl_matrix_fixed< ParametersValueType, VInputDimension, VOutputDimension > InverseJacobianPositionType
A templated class holding a n-Dimensional vector.
Definition itkVector.h:63
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....