ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkDisplacementFieldTransform.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 itkDisplacementFieldTransform_h
19#define itkDisplacementFieldTransform_h
20
21#include "itkTransform.h"
22
23#include "itkImage.h"
24#include "itkVectorImage.h"
28
29namespace itk
30{
31
86template <typename TParametersValueType, unsigned int VDimension>
87class ITK_TEMPLATE_EXPORT DisplacementFieldTransform : public Transform<TParametersValueType, VDimension, VDimension>
88{
89public:
90 ITK_DISALLOW_COPY_AND_MOVE(DisplacementFieldTransform);
91
97
99 itkOverrideGetNameOfClassMacro(DisplacementFieldTransform);
100
102 itkNewMacro(Self);
103
106
108 using typename Superclass::ScalarType;
109
111 using typename Superclass::FixedParametersType;
113 using typename Superclass::ParametersType;
114 using typename Superclass::ParametersValueType;
115
117 using typename Superclass::JacobianType;
120
123
126
128 using typename Superclass::InputPointType;
129 using typename Superclass::OutputPointType;
130
132 using typename Superclass::InputVectorType;
133 using typename Superclass::OutputVectorType;
134
137
141
143 using typename Superclass::InputVnlVectorType;
144 using typename Superclass::OutputVnlVectorType;
145
149
153
155 using typename Superclass::DerivativeType;
156
158 static constexpr unsigned int Dimension = VDimension;
159
165
167
176
180
186 virtual void
188 virtual void
190 void SetDisplacementField(std::nullptr_t) = delete;
191 itkGetModifiableObjectMacro(DisplacementField, DisplacementFieldType);
193
197 virtual void
199 itkGetModifiableObjectMacro(InverseDisplacementField, DisplacementFieldType);
201
205 virtual void
207 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
209
213 virtual void
215 itkGetModifiableObjectMacro(InverseInterpolator, InterpolatorType);
217
219 itkGetConstReferenceMacro(DisplacementFieldSetTime, ModifiedTimeType);
220
224 TransformPoint(const InputPointType & inputPoint) const override;
225
230 TransformVector(const InputVectorType &) const override
231 {
232 itkExceptionMacro("TransformVector(Vector) unimplemented, use "
233 "TransformVector(Vector,Point)");
234 }
235
236 OutputVectorPixelType
237 TransformVector(const InputVectorPixelType &) const override
238 {
239 itkExceptionMacro("TransformVector(Vector) unimplemented, use "
240 "TransformVector(Vector,Point)");
241 }
242
243 OutputVnlVectorType
244 TransformVector(const InputVnlVectorType &) const override
245 {
246 itkExceptionMacro("TransformVector(Vector) unimplemented, use "
247 "TransformVector(Vector,Point)");
248 }
249
250
253 using Superclass::TransformDiffusionTensor3D;
254 OutputDiffusionTensor3DType
256 {
257 itkExceptionMacro("TransformDiffusionTensor(Tensor) unimplemented, use "
258 "TransformDiffusionTensor(Tensor,Point)");
259 }
260
261 OutputVectorPixelType
263 {
264 itkExceptionMacro("TransformDiffusionTensor(Tensor) unimplemented, use "
265 "TransformDiffusionTensor(Tensor,Point)");
266 }
267
268
271 using Superclass::TransformCovariantVector;
272 OutputCovariantVectorType
274 {
275 itkExceptionMacro("TransformCovariantVector(CovariantVector) "
276 "unimplemented, use TransformCovariantVector(CovariantVector,Point)");
277 }
278
279 OutputVectorPixelType
281 {
282 itkExceptionMacro("TransformCovariantVector(CovariantVector) "
283 "unimplemented, use TransformCovariantVector(CovariantVector,Point)");
284 }
285
286
289 void
290 SetParameters(const ParametersType & params) override
291 {
292 if (&(this->m_Parameters) != &params)
293 {
294 if (params.Size() != this->m_Parameters.Size())
295 {
296 itkExceptionMacro("Input parameters size (" << params.Size() << ") does not match internal size ("
297 << this->m_Parameters.Size() << ").");
298 }
299 // Copy into existing object
300 this->m_Parameters = params;
301 this->Modified();
302 }
303 }
304
314 void
316
338 void
340 {
341 j = this->m_IdentityJacobian;
342 }
343
350 virtual void
355
360 void
363
368 void
370 InverseJacobianPositionType & jacobian) const override;
372
377 virtual void
379
391 virtual void
393 JacobianPositionType & jacobian,
394 bool useSVD = false) const;
395
407 virtual void
409 JacobianPositionType & jacobian,
410 bool useSVD = false) const;
411
412 void
413 UpdateTransformParameters(const DerivativeType & update, ScalarType factor = 1.0) override;
414
417 bool
418 GetInverse(Self * inverse) const;
419
423 GetInverseTransform() const override;
424
425 virtual void
427
430 GetTransformCategory() const override
431 {
432 return Self::TransformCategoryEnum::DisplacementField;
433 }
434
435 NumberOfParametersType
437 {
438 return Dimension;
439 }
440
449 itkSetMacro(CoordinateTolerance, double);
450 itkGetConstMacro(CoordinateTolerance, double);
452
461 itkSetMacro(DirectionTolerance, double);
462 itkGetConstMacro(DirectionTolerance, double);
464
465protected:
467 ~DisplacementFieldTransform() override = default;
468 void
469 PrintSelf(std::ostream & os, Indent indent) const override;
470
474
478
482
486
487private:
496 virtual void
498 JacobianPositionType & jacobian,
499 bool doInverseJacobian) const;
500
505 virtual void
507
512 virtual void
514
517};
518
519} // end namespace itk
520
521#ifndef ITK_MANUAL_INSTANTIATION
522# include "itkDisplacementFieldTransform.hxx"
523#endif
524
525#endif // itkDisplacementFieldTransform_h
SizeValueType Size() const
Definition itkArray.h:130
A templated class holding a n-Dimensional covariant vector.
Provides local/dense/high-dimensionality transformation via a a displacement field.
CovariantVector< ScalarType, InputDiffusionTensor3DType::Dimension > InputTensorEigenVectorType
NumberOfParametersType GetNumberOfLocalParameters() const override
void ComputeInverseJacobianWithRespectToPosition(const InputPointType &point, InverseJacobianPositionType &jacobian) const override
virtual void ComputeJacobianWithRespectToParameters(const IndexType &, JacobianType &j) const
VariableLengthVector< CoordinateRepresentationType > InputVectorPixelType
typename DisplacementFieldType::SpacingType SpacingType
typename DisplacementFieldType::PixelType PixelType
virtual void GetInverseJacobianOfForwardFieldWithRespectToPosition(const InputPointType &point, JacobianPositionType &jacobian, bool useSVD=false) const
virtual void SetInverseInterpolator(InterpolatorType *interpolator)
OutputVectorPixelType TransformDiffusionTensor(const InputVectorPixelType &) const
typename DisplacementFieldType::IndexType IndexType
void ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &j) const override
void SetParameters(const ParametersType &params) override
virtual void ComputeJacobianWithRespectToPosition(const IndexType &index, JacobianPositionType &jacobian) const
OutputPointType TransformPoint(const InputPointType &inputPoint) const override
vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension > JacobianPositionType
Point< CoordinateRepresentationType, VInputDimension > InputPointType
typename DisplacementFieldType::Pointer DisplacementFieldPointer
OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const override
vnl_vector_fixed< CoordinateRepresentationType, VInputDimension > InputVnlVectorType
CovariantVector< ScalarType, OutputDiffusionTensor3DType::Dimension > OutputTensorEigenVectorType
typename DisplacementFieldType::RegionType RegionType
OutputVnlVectorType TransformVector(const InputVnlVectorType &) const override
OutputDiffusionTensor3DType TransformDiffusionTensor(const InputDiffusionTensor3DType &) const
virtual void SetFixedParametersFromDisplacementField() const
typename DisplacementFieldType::DirectionType DirectionType
virtual void SetInterpolator(InterpolatorType *interpolator)
void ComputeJacobianWithRespectToPosition(const InputPointType &point, JacobianPositionType &jacobian) const override
VectorInterpolateImageFunction< DisplacementFieldType, ScalarType > InterpolatorType
~DisplacementFieldTransform() override=default
OutputVectorPixelType TransformCovariantVector(const InputVectorPixelType &) const override
virtual void SetDisplacementField(DisplacementFieldType *field)
void SetFixedParameters(const FixedParametersType &) override
DiffusionTensor3D< CoordinateRepresentationType > InputDiffusionTensor3DType
virtual void GetInverseJacobianOfForwardFieldWithRespectToPosition(const IndexType &index, JacobianPositionType &jacobian, bool useSVD=false) const
virtual void SetInverseDisplacementField(DisplacementFieldType *inverseField)
void SetDisplacementField(std::nullptr_t)=delete
typename DisplacementFieldType::PointType PointType
virtual void ComputeJacobianWithRespectToPositionInternal(const IndexType &index, JacobianPositionType &jacobian, bool doInverseJacobian) const
typename DisplacementFieldType::ConstPointer DisplacementFieldConstPointer
ImageVectorOptimizerParametersHelper< ScalarType, OutputVectorType::Dimension, Dimension > OptimizerParametersHelperType
Transform< TParametersValueType, VDimension, VDimension > Superclass
OutputVectorType TransformVector(const InputVectorType &) const override
void UpdateTransformParameters(const DerivativeType &update, ScalarType factor=1.0) override
virtual void SetDisplacementField(VectorImageDisplacementFieldType *field)
void PrintSelf(std::ostream &os, Indent indent) const override
Image< OutputVectorType, Dimension > DisplacementFieldType
typename DisplacementFieldType::SizeType SizeType
OutputVectorPixelType TransformVector(const InputVectorPixelType &) const override
TransformCategoryEnum GetTransformCategory() const override
Point< CoordinateRepresentationType, VOutputDimension > OutputPointType
bool GetInverse(Self *inverse) const
Vector< CoordinateRepresentationType, VInputDimension > InputVectorType
VectorImage< TParametersValueType, Dimension > VectorImageDisplacementFieldType
CovariantVector< CoordinateRepresentationType, VInputDimension > InputCovariantVectorType
Vector< CoordinateRepresentationType, VOutputDimension > OutputVectorType
virtual void VerifyFixedParametersInformation()
vnl_matrix_fixed< ParametersValueType, VInputDimension, VOutputDimension > InverseJacobianPositionType
InverseTransformBasePointer GetInverseTransform() const override
Class to hold and manage parameters of type Image<Vector<...>,...>, used in Transforms,...
Templated n-dimensional image class.
Definition itkImage.h:89
SmartPointer< const Self > ConstPointer
Definition itkImage.h:97
Point< PointValueType, VImageDimension > PointType
ImageRegion< VImageDimension > RegionType
Vector< SpacingValueType, VImageDimension > SpacingType
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void Modified() const
Implements transparent reference counting.
OptimizerParameters< ParametersValueType > ParametersType
VariableLengthVector< TParametersValueType > OutputVectorPixelType
DiffusionTensor3D< TParametersValueType > InputDiffusionTensor3DType
ParametersType m_Parameters
virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType &pnt, InverseJacobianPositionType &jacobian) const
virtual void ComputeJacobianWithRespectToPosition(const InputPointType &x, JacobianPositionType &jacobian) const
IdentifierType NumberOfParametersType
TParametersValueType ParametersValueType
Point< TParametersValueType, VOutputDimension > OutputPointType
OptimizerParameters< FixedParametersValueType > FixedParametersType
VariableLengthVector< TParametersValueType > InputVectorPixelType
CovariantVector< TParametersValueType, VOutputDimension > OutputCovariantVectorType
virtual OutputVectorType TransformVector(const InputVectorType &) const
TransformBaseTemplateEnums::TransformCategory TransformCategoryEnum
double FixedParametersValueType
vnl_vector_fixed< TParametersValueType, VInputDimension > InputVnlVectorType
vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension > JacobianPositionType
vnl_matrix_fixed< ParametersValueType, VInputDimension, VOutputDimension > InverseJacobianPositionType
DiffusionTensor3D< TParametersValueType > OutputDiffusionTensor3DType
Array2D< ParametersValueType > JacobianType
Point< TParametersValueType, VInputDimension > InputPointType
typename InverseTransformBaseType::Pointer InverseTransformBasePointer
CovariantVector< TParametersValueType, VInputDimension > InputCovariantVectorType
Vector< TParametersValueType, VInputDimension > InputVectorType
ParametersValueType ScalarType
Vector< TParametersValueType, VOutputDimension > OutputVectorType
vnl_vector_fixed< TParametersValueType, VOutputDimension > OutputVnlVectorType
Templated n-dimensional vector image class.
Base class for all vector image interpolators.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType ModifiedTimeType