ITK  5.2.0
Insight Toolkit
itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.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  * http://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 itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor_h
19 #define itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor_h
20 
22 
23 namespace itk
24 {
65 template <typename TTransform>
67  : public TransformParametersAdaptor<TTransform>
68 {
69 public:
71 
77 
79  itkNewMacro(Self);
80 
83 
85  using TransformType = TTransform;
86  using TransformPointer = typename TransformType::Pointer;
87  using ParametersType = typename TransformType::ParametersType;
88  using ParametersValueType = typename TransformType::ParametersValueType;
89  using FixedParametersType = typename TransformType::FixedParametersType;
90  using FixedParametersValueType = typename TransformType::FixedParametersValueType;
91 
93  typename TransformType::TimeVaryingVelocityFieldControlPointLatticeType;
95  typename TimeVaryingVelocityFieldControlPointLatticeType::Pointer;
98  using VectorType = typename TimeVaryingVelocityFieldControlPointLatticeType::PixelType;
100  using SpacingType = typename TimeVaryingVelocityFieldControlPointLatticeType::SpacingType;
105 
107  static constexpr unsigned int TotalDimension = TransformType::Dimension + 1;
108 
110  itkSetMacro(SplineOrder, SizeValueType);
111 
113  itkGetConstMacro(SplineOrder, SizeValueType);
114 
116  void
117  SetRequiredTransformDomainMeshSize(const MeshSizeType &);
118 
120  itkGetConstReferenceMacro(RequiredTransformDomainMeshSize, MeshSizeType);
121 
123  void
124  SetRequiredTransformDomainSize(const SizeType &);
125 
127  itkGetConstReferenceMacro(RequiredTransformDomainSize, SizeType);
128 
130  void
131  SetRequiredTransformDomainSpacing(const SpacingType &);
132 
134  itkGetConstReferenceMacro(RequiredTransformDomainSpacing, SpacingType);
135 
137  void
138  SetRequiredTransformDomainOrigin(const OriginType &);
139 
141  itkGetConstReferenceMacro(RequiredTransformDomainOrigin, OriginType);
142 
144  void
145  SetRequiredTransformDomainDirection(const DirectionType &);
146 
148  itkGetConstReferenceMacro(RequiredTransformDomainDirection, DirectionType);
149 
151  const OriginType
153  {
154  OriginType requiredLatticeOrigin;
155  for (SizeValueType i = 0; i < TotalDimension; i++)
156  {
157  requiredLatticeOrigin[i] = this->m_RequiredFixedParameters[TotalDimension + i];
158  }
159  return requiredLatticeOrigin;
160  }
162 
164  const SpacingType
166  {
167  SpacingType requiredLatticeSpacing;
168  for (SizeValueType i = 0; i < TotalDimension; i++)
169  {
170  FixedParametersValueType domainPhysicalDimensions =
171  static_cast<FixedParametersValueType>(this->m_RequiredTransformDomainSize[i] - 1.0) *
172  this->m_RequiredTransformDomainSpacing[i];
173  requiredLatticeSpacing[i] =
174  domainPhysicalDimensions / static_cast<FixedParametersValueType>(this->m_RequiredTransformDomainMeshSize[i]);
175  }
176  return requiredLatticeSpacing;
177  }
179 
181  const SizeType
183  {
184  SizeType requiredLatticeSize;
185  for (SizeValueType i = 0; i < TotalDimension; i++)
186  {
187  requiredLatticeSize[i] = static_cast<SizeValueType>(this->m_RequiredFixedParameters[i]);
188  }
189  return requiredLatticeSize;
190  }
192 
194  const DirectionType
196  {
197  return this->m_RequiredTransformDomainDirection;
198  }
199 
201  void
202  AdaptTransformParameters() override;
203 
204  void
205  SetRequiredFixedParameters(const FixedParametersType) override;
206 
207 protected:
210 
211  void
212  PrintSelf(std::ostream & os, Indent indent) const override;
213 
214 private:
216  void
217  UpdateRequiredFixedParameters();
218 
224 
226 
227 }; // class TimeVaryingBSplineVelocityFieldTransformParametersAdaptor
228 } // namespace itk
229 
230 #ifndef ITK_MANUAL_INSTANTIATION
231 # include "itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.hxx"
232 #endif
233 
234 #endif /* itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor_h */
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::m_SplineOrder
SizeValueType m_SplineOrder
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:225
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::SizeType
typename TimeVaryingVelocityFieldControlPointLatticeType::SizeType SizeType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:101
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::SpacingType
typename TimeVaryingVelocityFieldControlPointLatticeType::SpacingType SpacingType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:100
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::TransformPointer
typename TransformType::Pointer TransformPointer
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:86
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::ParametersValueType
typename TransformType::ParametersValueType ParametersValueType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:88
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::TransformType
TTransform TransformType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:85
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::TimeVaryingVelocityFieldControlPointLatticePointer
typename TimeVaryingVelocityFieldControlPointLatticeType::Pointer TimeVaryingVelocityFieldControlPointLatticePointer
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:95
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::m_RequiredTransformDomainSpacing
SpacingType m_RequiredTransformDomainSpacing
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:222
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::m_RequiredTransformDomainOrigin
OriginType m_RequiredTransformDomainOrigin
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:220
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::m_RequiredTransformDomainMeshSize
MeshSizeType m_RequiredTransformDomainMeshSize
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:219
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::DirectionType
typename TimeVaryingVelocityFieldControlPointLatticeType::DirectionType DirectionType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:104
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::GetRequiredControlPointLatticeDirection
const DirectionType GetRequiredControlPointLatticeDirection() const
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:195
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::FixedParametersType
typename TransformType::FixedParametersType FixedParametersType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:89
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::m_RequiredTransformDomainSize
SizeType m_RequiredTransformDomainSize
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:223
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::SizeValueType
typename TimeVaryingVelocityFieldControlPointLatticeType::SizeValueType SizeValueType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:102
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::OriginType
typename TimeVaryingVelocityFieldControlPointLatticeType::PointType OriginType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:99
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::GetRequiredControlPointLatticeOrigin
const OriginType GetRequiredControlPointLatticeOrigin() const
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:152
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::RegionType
typename TimeVaryingVelocityFieldControlPointLatticeType::RegionType RegionType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:96
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::MeshSizeType
typename TimeVaryingVelocityFieldControlPointLatticeType::SizeType MeshSizeType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:103
itk::TransformParametersAdaptor
Base helper class intended for multi-resolution image registration.
Definition: itkTransformParametersAdaptor.h:55
itkTransformParametersAdaptor.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::ParametersType
typename TransformType::ParametersType ParametersType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:87
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::GetRequiredControlPointLatticeSize
const SizeType GetRequiredControlPointLatticeSize() const
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:182
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::IndexType
typename TimeVaryingVelocityFieldControlPointLatticeType::IndexType IndexType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:97
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::m_RequiredTransformDomainDirection
DirectionType m_RequiredTransformDomainDirection
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:221
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::VectorType
typename TimeVaryingVelocityFieldControlPointLatticeType::PixelType VectorType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:98
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::FixedParametersValueType
typename TransformType::FixedParametersValueType FixedParametersValueType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:90
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::TimeVaryingVelocityFieldControlPointLatticeType
typename TransformType::TimeVaryingVelocityFieldControlPointLatticeType TimeVaryingVelocityFieldControlPointLatticeType
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:93
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor
TimeVaryingBSplineVelocityFieldTransformParametersAdaptor is a helper class intended to definition.
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:66
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::TimeVaryingBSplineVelocityFieldTransformParametersAdaptor::GetRequiredControlPointLatticeSpacing
const SpacingType GetRequiredControlPointLatticeSpacing() const
Definition: itkTimeVaryingBSplineVelocityFieldTransformParametersAdaptor.h:165