go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
elxBSplineStackTransform.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright UMC Utrecht and contributors
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 elxBSplineStackTransform_h
19#define elxBSplineStackTransform_h
20
21#include "elxIncludes.h" // include first to avoid MSVS warning
22
27
31
32namespace elastix
33{
34
111template <class TElastix>
112class ITK_TEMPLATE_EXPORT BSplineStackTransform
113 : public itk::AdvancedCombinationTransform<typename elx::TransformBase<TElastix>::CoordRepType,
114 elx::TransformBase<TElastix>::FixedImageDimension>
115 , public TransformBase<TElastix>
116{
117public:
119
125 using Pointer = itk::SmartPointer<Self>;
126 using ConstPointer = itk::SmartPointer<const Self>;
127
129 itkNewMacro(Self);
130
133
138 elxClassNameMacro("BSplineStackTransform");
139
141 itkStaticConstMacro(SpaceDimension, unsigned int, Superclass2::FixedImageDimension);
142 itkStaticConstMacro(ReducedSpaceDimension, unsigned int, Superclass2::FixedImageDimension - 1);
143
149 Self::SpaceDimension>;
151
155 Self::ReducedSpaceDimension>;
157
159 using typename Superclass1::ParametersType;
160 using typename Superclass2::ParameterMapType;
161 using typename Superclass1::NumberOfParametersType;
162
164 using PixelType = typename BSplineTransformBaseType::PixelType;
173
175 using typename Superclass2::ElastixType;
176 using typename Superclass2::RegistrationType;
177 using typename Superclass2::CoordRepType;
178 using typename Superclass2::FixedImageType;
179 using typename Superclass2::MovingImageType;
182
184 using ReducedDimensionImageType = itk::Image<PixelType, Self::ReducedSpaceDimension>;
185 using ReducedDimensionRegionType = itk::ImageRegion<Self::ReducedSpaceDimension>;
186 using ReducedDimensionSizeType = typename ReducedDimensionRegionType::SizeType;
187 using ReducedDimensionIndexType = typename ReducedDimensionRegionType::IndexType;
188 using ReducedDimensionSpacingType = typename ReducedDimensionImageType::SpacingType;
189 using ReducedDimensionDirectionType = typename ReducedDimensionImageType::DirectionType;
190 using ReducedDimensionOriginType = typename ReducedDimensionImageType::PointType;
191
195 using GridScheduleType = typename GridScheduleComputerType ::VectorGridSpacingFactorType;
198
200 using ParameterValueType = std::vector<std::string>;
201
207 int
208 BeforeAll() override;
209
223 void
225
230 void
232
239 virtual void
241
243 void
244 ReadFromFile() override;
245
247 virtual void
248 SetOptimizerScales(const unsigned int edgeWidth);
249
250protected:
252 BSplineStackTransform() { this->Superclass1::SetCurrentTransform(m_StackTransform); }
253
255 ~BSplineStackTransform() override = default;
256
258 virtual void
260
261private:
263
269 void
271
275
279
281 const typename StackTransformType::Pointer m_StackTransform{ StackTransformType::New() };
282
285
289
291 unsigned int m_SplineOrder;
292
295 double m_StackOrigin, m_StackSpacing;
296
298 unsigned int
300};
301
302} // end namespace elastix
303
304#ifndef ITK_MANUAL_INSTANTIATION
305# include "elxBSplineStackTransform.hxx"
306#endif
307
308#endif // end #ifndef elxBSplineStackTransform_h
A B-spline transform based on the itkStackTransform.
itk::ImageRegion< Self::ReducedSpaceDimension > ReducedDimensionRegionType
typename ReducedDimensionRegionType::IndexType ReducedDimensionIndexType
typename BSplineTransformBaseType::Pointer BSplineTransformBasePointer
typename GridScheduleComputerType ::VectorGridSpacingFactorType GridScheduleType
typename ReducedDimensionImageType::SpacingType ReducedDimensionSpacingType
typename BSplineTransformBaseType::ImagePointer ImagePointer
typename Superclass2::CombinationTransformType CombinationTransformType
typename BSplineTransformBaseType::RegionType RegionType
ParameterMapType CreateDerivedTransformParameterMap() const override
typename ReducedDimensionImageType::PointType ReducedDimensionOriginType
GridScheduleComputerPointer m_GridScheduleComputer
typename Superclass2::ITKBaseType ITKBaseType
typename BSplineTransformBaseType::ImageType ImageType
typename ReducedDimensionImageType::DirectionType ReducedDimensionDirectionType
itk::SmartPointer< const Self > ConstPointer
typename GridScheduleComputerType::Pointer GridScheduleComputerPointer
ReducedDimensionBSplineTransformBasePointer m_DummySubTransform
itk::Image< PixelType, Self::ReducedSpaceDimension > ReducedDimensionImageType
ITK_DISALLOW_COPY_AND_MOVE(BSplineStackTransform)
typename BSplineTransformBaseType::IndexType IndexType
itkStaticConstMacro(SpaceDimension, unsigned int, Superclass2::FixedImageDimension)
typename TElastix::ParameterMapType ParameterMapType
typename BSplineTransformBaseType::OriginType OriginType
void BeforeRegistration() override
typename ReducedDimensionRegionType::SizeType ReducedDimensionSizeType
void BeforeEachResolution() override
typename BSplineTransformBaseType::SizeType SizeType
~BSplineStackTransform() override=default
std::vector< std::string > ParameterValueType
elxClassNameMacro("BSplineStackTransform")
virtual void PreComputeGridInformation()
typename BSplineTransformBaseType::DirectionType DirectionType
unsigned int InitializeBSplineTransform()
typename GridUpsamplerType::Pointer GridUpsamplerPointer
typename BSplineTransformBaseType::SpacingType SpacingType
virtual void SetOptimizerScales(const unsigned int edgeWidth)
typename ReducedDimensionBSplineTransformBaseType::Pointer ReducedDimensionBSplineTransformBasePointer
itkStaticConstMacro(ReducedSpaceDimension, unsigned int, Superclass2::FixedImageDimension - 1)
This class is the elastix base class for all Transforms.
typename TElastix::FixedImageType FixedImageType
typename TElastix::ParameterMapType ParameterMapType
typename TElastix::MovingImageType MovingImageType
typename ElastixType::RegistrationBaseType RegistrationType
Base class for deformable transform using a B-spline representation.
This class combines two transforms: an 'initial transform' with a 'current transform'.
itk::SmartPointer< BSplineStackTransform > Pointer
This class computes all information about the B-spline grid, given the image information and the desi...
Convenience class for upsampling a B-spline coefficient image.


Generated on 2024-07-17 for elastix by doxygen 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) elastix logo