ITK  6.0.0
Insight Toolkit
itkBSplineTransform.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 itkBSplineTransform_h
19 #define itkBSplineTransform_h
20 
22 
23 namespace itk
24 {
102 template <typename TParametersValueType = double, unsigned int VDimension = 3, unsigned int VSplineOrder = 3>
103 class ITK_TEMPLATE_EXPORT BSplineTransform : public BSplineBaseTransform<TParametersValueType, VDimension, VSplineOrder>
104 {
105 public:
106  ITK_DISALLOW_COPY_AND_MOVE(BSplineTransform);
107 
113 
115  itkNewMacro(Self);
116 
118  itkOverrideGetNameOfClassMacro(BSplineTransform);
119 
121  static constexpr unsigned int SpaceDimension = VDimension;
122 
124  static constexpr unsigned int SplineOrder = VSplineOrder;
125 
127  using typename Superclass::ScalarType;
128 
130  using typename Superclass::ParametersType;
131  using typename Superclass::ParametersValueType;
132  using typename Superclass::FixedParametersType;
133  using typename Superclass::FixedParametersValueType;
134 
136  using typename Superclass::JacobianType;
137  using typename Superclass::JacobianPositionType;
139 
141  using typename Superclass::NumberOfParametersType;
142 
144  using typename Superclass::InputVectorType;
145  using typename Superclass::OutputVectorType;
146 
150 
152  using typename Superclass::InputVnlVectorType;
153  using typename Superclass::OutputVnlVectorType;
154 
156  using typename Superclass::InputPointType;
157  using typename Superclass::OutputPointType;
158 
159 
160  std::string
161  GetTransformTypeAsString() const override;
162 
185  void
186  SetFixedParameters(const FixedParametersType & passedParameters) override;
190  using typename Superclass::ImageType;
191  using typename Superclass::ImagePointer;
192  using typename Superclass::CoefficientImageArray;
193 
205  void
206  SetCoefficientImages(const CoefficientImageArray & images) override;
207 
209  using typename Superclass::RegionType;
210 
211  using typename Superclass::IndexType;
212  using typename Superclass::SizeType;
213  using typename Superclass::SpacingType;
214  using typename Superclass::DirectionType;
215  using typename Superclass::OriginType;
216 
218  using typename Superclass::WeightsFunctionType;
219 
220  using typename Superclass::WeightsType;
221  using typename Superclass::ContinuousIndexType;
222 
224  using typename Superclass::ParameterIndexArrayType;
225 
234  using Superclass::TransformPoint;
235  void
236  TransformPoint(const InputPointType & point,
237  OutputPointType & outputPoint,
238  WeightsType & weights,
239  ParameterIndexArrayType & indices,
240  bool & inside) const override;
244  void
245  ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &) const override;
246 
248  NumberOfParametersType
249  GetNumberOfParameters() const override;
250 
252  NumberOfParametersType
253  GetNumberOfParametersPerDimension() const override;
254 
255  using PhysicalDimensionsType = typename Superclass::SpacingType;
256  using typename Superclass::PixelType;
257 
258  using typename Superclass::MeshSizeType;
259 
261  virtual void
262  SetTransformDomainOrigin(const OriginType &);
263 
265  virtual OriginType
266  GetTransformDomainOrigin() const;
267 
269  virtual void
270  SetTransformDomainPhysicalDimensions(const PhysicalDimensionsType &);
271 
273  virtual PhysicalDimensionsType
274  GetTransformDomainPhysicalDimensions() const;
275 
277  virtual void
278  SetTransformDomainDirection(const DirectionType &);
279 
281  virtual DirectionType
282  GetTransformDomainDirection() const;
283 
285  virtual void
286  SetTransformDomainMeshSize(const MeshSizeType &);
287 
289  virtual MeshSizeType
290  GetTransformDomainMeshSize() const;
291 
292 protected:
294  void
295  PrintSelf(std::ostream & os, Indent indent) const override;
296 
298  ~BSplineTransform() override = default;
299 
300 private:
302  void
303  SetCoefficientImageInformationFromFixedParameters() override;
304 
306  void
308  {}
309  void
311  {}
312  void
314  {}
315  void
317  {}
321  bool
322  InsideValidRegion(ContinuousIndexType &) const override;
323 
324  void
325  SetFixedParametersFromCoefficientImageInformation();
326 
327  void
328  SetFixedParametersFromTransformDomainInformation(const OriginType & meshOrigin,
329  const PhysicalDimensionsType & meshPhysical,
330  const DirectionType & meshDirection,
331  const MeshSizeType & meshSize);
332 
333 }; // class BSplineTransform
334 } // namespace itk
335 
336 #ifndef ITK_MANUAL_INSTANTIATION
337 # include "itkBSplineTransform.hxx"
338 #endif
339 
340 #endif /* itkBSplineTransform_h */
itk::BSplineBaseTransform::MeshSizeType
SizeType MeshSizeType
Definition: itkBSplineBaseTransform.h:333
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::BSplineBaseTransform::PhysicalDimensionsType
typename ImageType::SpacingType PhysicalDimensionsType
Definition: itkBSplineBaseTransform.h:330
itk::Transform< TParametersValueType, VDimension, VDimension >::OutputVnlVectorType
vnl_vector_fixed< TParametersValueType, VOutputDimension > OutputVnlVectorType
Definition: itkTransform.h:157
itk::BSplineBaseTransform::WeightsType
typename WeightsFunctionType::WeightsType WeightsType
Definition: itkBSplineBaseTransform.h:235
itk::BSplineBaseTransform::DirectionType
typename ImageType::DirectionType DirectionType
Definition: itkBSplineBaseTransform.h:225
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BSplineBaseTransform::OriginType
typename ImageType::PointType OriginType
Definition: itkBSplineBaseTransform.h:226
itk::BSplineTransform
Deformable transform using a BSpline representation.
Definition: itkBSplineTransform.h:103
itk::BSplineTransform::SetFixedParametersGridOriginFromTransformDomainInformation
void SetFixedParametersGridOriginFromTransformDomainInformation() const override
Definition: itkBSplineTransform.h:310
itk::point
*par Constraints *The filter requires an 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
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::Transform< TParametersValueType, VDimension, VDimension >::ScalarType
ParametersValueType ScalarType
Definition: itkTransform.h:127
itk::FixedArray< ImagePointer, VDimension >
itk::BSplineTransform::SetFixedParametersGridSizeFromTransformDomainInformation
void SetFixedParametersGridSizeFromTransformDomainInformation() const override
Definition: itkBSplineTransform.h:307
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::Transform< TParametersValueType, VDimension, VDimension >::JacobianPositionType
vnl_matrix_fixed< ParametersValueType, VOutputDimension, VInputDimension > JacobianPositionType
Definition: itkTransform.h:131
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itkBSplineBaseTransform.h
itk::BSplineTransform::SetFixedParametersGridSpacingFromTransformDomainInformation
void SetFixedParametersGridSpacingFromTransformDomainInformation() const override
Definition: itkBSplineTransform.h:313
itk::Transform< TParametersValueType, VDimension, VDimension >::InverseJacobianPositionType
vnl_matrix_fixed< ParametersValueType, VInputDimension, VOutputDimension > InverseJacobianPositionType
Definition: itkTransform.h:132
itk::BSplineBaseTransform
A base class with common elements of BSplineTransform and BSplineDeformableTransform.
Definition: itkBSplineBaseTransform.h:34
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Transform< TParametersValueType, VDimension, VDimension >
itk::BSplineTransform::SetFixedParametersGridDirectionFromTransformDomainInformation
void SetFixedParametersGridDirectionFromTransformDomainInformation() const override
Definition: itkBSplineTransform.h:316
itk::Array2D
Array2D class representing a 2D array.
Definition: itkArray2D.h:42
itk::images
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar images
Definition: itkVectorGradientMagnitudeImageFilter.h:107
itk::Transform< TParametersValueType, VDimension, VDimension >::InputVnlVectorType
vnl_vector_fixed< TParametersValueType, VInputDimension > InputVnlVectorType
Definition: itkTransform.h:156