ITK  6.0.0
Insight Toolkit
itkTimeVaryingBSplineVelocityFieldImageRegistrationMethod.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 itkTimeVaryingBSplineVelocityFieldImageRegistrationMethod_h
19#define itkTimeVaryingBSplineVelocityFieldImageRegistrationMethod_h
20
22
26
27namespace itk
28{
29
86template <typename TFixedImage,
87 typename TMovingImage,
88 typename TOutputTransform = TimeVaryingBSplineVelocityFieldTransform<double, TFixedImage::ImageDimension>,
89 typename TVirtualImage = TFixedImage,
90 typename TPointSet = PointSet<unsigned int, TFixedImage::ImageDimension>>
92 : public ImageRegistrationMethodv4<TFixedImage, TMovingImage, TOutputTransform, TVirtualImage, TPointSet>
93{
94public:
96
102
104 itkNewMacro(Self);
105
107 static constexpr unsigned int ImageDimension = TFixedImage::ImageDimension;
108
111
113 using FixedImageType = TFixedImage;
115 using typename Superclass::FixedImagesContainerType;
116 using MovingImageType = TMovingImage;
118 using typename Superclass::MovingImagesContainerType;
119
120 using InputPointSetType = TPointSet;
122 using typename Superclass::PointSetsContainerType;
123
125 using typename Superclass::ImageMetricType;
129
130 using VirtualImageType = typename Superclass::VirtualImageType;
131 using typename Superclass::VirtualImageBaseType;
132 using typename Superclass::VirtualImageBaseConstPointer;
133
134 using typename Superclass::MetricType;
135 using typename Superclass::MultiMetricType;
137 using typename Superclass::PointSetMetricType;
138
141 using typename Superclass::FixedImageMaskType;
143 using typename Superclass::FixedImageMasksContainerType;
144 using typename Superclass::MovingImageMaskType;
146 using typename Superclass::MovingImageMasksContainerType;
147
148 using typename Superclass::InitialTransformType;
149 using OutputTransformType = TOutputTransform;
151 using RealType = typename OutputTransformType::ScalarType;
152 using DerivativeType = typename OutputTransformType::DerivativeType;
153 using DerivativeValueType = typename DerivativeType::ValueType;
154 using DisplacementFieldType = typename OutputTransformType::DisplacementFieldType;
156
158
160 typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticeType;
162 typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticePointer;
163 using TimeVaryingVelocityFieldType = typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticeType;
165 typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticePointer;
166 using DisplacementVectorType = typename TimeVaryingVelocityFieldControlPointLatticeType::PixelType;
167
168 using typename Superclass::CompositeTransformType;
170
171 using typename Superclass::DecoratedOutputTransformType;
173
175
178
184 using WeightsElementType = typename WeightsContainerType::Element;
187
189 itkSetMacro(LearningRate, RealType);
190 itkGetConstMacro(LearningRate, RealType);
194 itkSetMacro(NumberOfIterationsPerLevel, NumberOfIterationsArrayType);
195 itkGetConstMacro(NumberOfIterationsPerLevel, NumberOfIterationsArrayType);
199 itkSetMacro(ConvergenceThreshold, RealType);
200 itkGetConstMacro(ConvergenceThreshold, RealType);
204 itkSetMacro(ConvergenceWindowSize, unsigned int);
205 itkGetConstMacro(ConvergenceWindowSize, unsigned int);
209 itkSetMacro(NumberOfTimePointSamples, SizeValueType);
210 itkGetConstMacro(NumberOfTimePointSamples, SizeValueType);
213protected:
216 void
217 PrintSelf(std::ostream & os, Indent indent) const override;
218
220 void
221 GenerateData() override;
222
227 virtual void
229
231 void
233
234 void
240 const TransformBaseType *,
243 const TransformBaseType *,
245
246private:
247 DisplacementFieldTransformPointer m_IdentityDisplacementFieldTransform{};
248
249 RealType m_LearningRate{};
250
251 RealType m_ConvergenceThreshold{};
252 unsigned int m_ConvergenceWindowSize{ 10 };
253
254 NumberOfIterationsArrayType m_NumberOfIterationsPerLevel{};
255
256 SizeValueType m_NumberOfTimePointSamples{ 4 };
257
258 WeightsElementType m_BoundaryWeight{};
259};
260} // end namespace itk
261
262#ifndef ITK_MANUAL_INSTANTIATION
263# include "itkTimeVaryingBSplineVelocityFieldImageRegistrationMethod.hxx"
264#endif
265
266#endif
Image filter which provides a B-spline output approximation.
typename Superclass::TransformType TransformType
A templated class holding a point in n-Dimensional image space.
Provides local/dense/high-dimensionality transformation via a a displacement field.
Implementation of an image mask as spatial object.
Interface method for the current registration framework.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
Interface method for the current registration framework using the time varying velocity field transfo...
typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticeType TimeVaryingVelocityFieldType
typename TimeVaryingVelocityFieldControlPointLatticeType::PixelType DisplacementVectorType
void PrintSelf(std::ostream &os, Indent indent) const override
void AttachMetricGradientPointSetAtSpecificTimePoint(const RealType, VelocityFieldPointSetType *, WeightsContainerType *, const FixedImagesContainerType, const PointSetsContainerType, const TransformBaseType *, const MovingImagesContainerType, const PointSetsContainerType, const TransformBaseType *, const FixedImageMasksContainerType)
void GetMetricDerivativePointSetForAllTimePoints(VelocityFieldPointSetType *, WeightsContainerType *)
typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticePointer TimeVaryingVelocityFieldPointer
typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticeType TimeVaryingVelocityFieldControlPointLatticeType
typename OutputTransformType::TimeVaryingVelocityFieldControlPointLatticePointer TimeVaryingVelocityFieldControlPointLatticePointer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86