File: itkBSplineTransformParametersAdaptor.h

package info (click to toggle)
insighttoolkit5 5.4.5-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 704,588 kB
  • sloc: cpp: 784,579; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,934; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 461; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (165 lines) | stat: -rw-r--r-- 6,274 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkBSplineTransformParametersAdaptor_h
#define itkBSplineTransformParametersAdaptor_h

#include "itkTransformParametersAdaptor.h"

namespace itk
{
/** \class BSplineTransformParametersAdaptor
 * \brief BSplineTransformParametersAdaptor adapts a BSplineTransform to the
 * new specified fixed parameters.
 *
 * The fixed parameters of the BSplineTransform store the following information
 * (in order as they appear in m_FixedParameters):
 *    grid size
 *    grid origin
 *    grid spacing
 *    grid direction
 *
 * During multiresolution image registration it is often desired to also increase
 * the B-spline grid resolution for greater flexibility in optimizing the
 * transform.  As defined in the base class, the user can change the resolution via
 *
     \code
     transformAdaptor->SetTransform( transform );
     transformAdaptor->SetRequiredFixedParameters( fixedParameters );
     transformAdaptor->AdaptTransformParameters();
     \endcode
 *
 * or the user can use the more intuitive API for setting the fixed parameters.
 * E.g., often the user will want to maintain the same transform domain spatial
 * extent but only increase the mesh size.  This can be done as follows:
 *
     \code
     transformAdaptor->SetTransform( transform );
     transformAdaptor->SetRequiredTransformDomainOrigin( transform->GetTransformDomainOrigin() );
     transformAdaptor->SetRequiredTransformDomainDirection( transform->GetTransformDomainDirection() );
     transformAdaptor->SetRequiredTransformDomainPhysicalDimensions( transform->GetTransformDomainPhysicalDimensions()
 ); transformAdaptor->SetRequiredTransformDomainMeshSize( newMeshSize ); transformAdaptor->AdaptTransformParameters();
     \endcode
 *
 * \author Nick Tustison
 * \author Marius Staring
 *
 * \ingroup ITKRegistrationCommon
 */
template <typename TTransform>
class ITK_TEMPLATE_EXPORT BSplineTransformParametersAdaptor : public TransformParametersAdaptor<TTransform>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(BSplineTransformParametersAdaptor);

  /** Standard class type aliases. */
  using Self = BSplineTransformParametersAdaptor;
  using Superclass = TransformParametersAdaptor<TTransform>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;

  /** New macro for creation of through a Smart Pointer. */
  itkNewMacro(Self);

  /** \see LightObject::GetNameOfClass() */
  itkOverrideGetNameOfClassMacro(BSplineTransformParametersAdaptor);

  /** Typedefs associated with the transform */
  using TransformType = TTransform;
  using TransformPointer = typename TransformType::Pointer;

  using typename Superclass::FixedParametersType;
  using typename Superclass::FixedParametersValueType;
  using typename Superclass::ParametersType;
  using typename Superclass::ParametersValueType;

  using OriginType = typename TransformType::OriginType;
  using SizeType = typename TransformType::SizeType;
  using SpacingType = typename TransformType::SpacingType;
  using IndexType = typename TransformType::IndexType;
  using MeshSizeType = typename TransformType::MeshSizeType;
  using DirectionType = typename TransformType::DirectionType;
  using PhysicalDimensionsType = typename TransformType::PhysicalDimensionsType;


  using ImageType = typename TransformType::ImageType;
  using RegionType = typename ImageType::RegionType;
  using CoefficientImageArray = typename TransformType::CoefficientImageArray;

  /** Dimension of parameters. */
  static constexpr unsigned int SpaceDimension = TransformType::SpaceDimension;

  /** Alternative method for setting the required mesh size. */
  void
  SetRequiredTransformDomainMeshSize(const MeshSizeType &);

  /** Get the required mesh size. */
  itkGetConstReferenceMacro(RequiredTransformDomainMeshSize, MeshSizeType);

  /** Alternative method for setting the required mesh size. */
  void
  SetRequiredTransformDomainPhysicalDimensions(const PhysicalDimensionsType &);

  /** Get the required physical dimensions. */
  itkGetConstReferenceMacro(RequiredTransformDomainPhysicalDimensions, PhysicalDimensionsType);

  /** Alternative method for setting the required origin. */
  void
  SetRequiredTransformDomainOrigin(const OriginType &);

  /** Get the required origin. */
  itkGetConstReferenceMacro(RequiredTransformDomainOrigin, OriginType);

  /** Alternative method for setting the required direction. */
  void
  SetRequiredTransformDomainDirection(const DirectionType &);

  /** Get the required direction. */
  itkGetConstReferenceMacro(RequiredTransformDomainDirection, DirectionType);

  void
  SetRequiredFixedParameters(const FixedParametersType) override;

  /** Initialize the transform using the specified fixed parameters */
  void
  AdaptTransformParameters() override;

protected:
  BSplineTransformParametersAdaptor();
  ~BSplineTransformParametersAdaptor() override = default;

  void
  PrintSelf(std::ostream & os, Indent indent) const override;

private:
  /** Helper function to set m_RequiredFixedParameters */
  void
  UpdateRequiredFixedParameters();

  MeshSizeType           m_RequiredTransformDomainMeshSize{};
  OriginType             m_RequiredTransformDomainOrigin{};
  DirectionType          m_RequiredTransformDomainDirection{};
  PhysicalDimensionsType m_RequiredTransformDomainPhysicalDimensions{};

}; // class BSplineTransformParametersAdaptor
} // namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkBSplineTransformParametersAdaptor.hxx"
#endif

#endif /* itkBSplineTransformParametersAdaptor_h */