go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkGridScheduleComputer.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 itkGridScheduleComputer_h
19#define itkGridScheduleComputer_h
20
21#include "itkObject.h"
22#include "itkImageBase.h"
23#include "itkTransform.h"
24
25namespace itk
26{
27
39template <typename TTransformScalarType, unsigned int VImageDimension>
40class ITK_TEMPLATE_EXPORT GridScheduleComputer : public Object
41{
42public:
44
47 using Superclass = Object;
48 using Pointer = SmartPointer<Self>;
49 using ConstPointer = SmartPointer<const Self>;
50
52 itkNewMacro(Self);
53
55 itkTypeMacro(GridScheduleComputer, Object);
56
58 itkStaticConstMacro(Dimension, unsigned int, VImageDimension);
59
61 using TransformScalarType = TTransformScalarType;
62 using ImageBaseType = ImageBase<Self::Dimension>;
63 using PointType = typename ImageBaseType::PointType;
64 using OriginType = typename ImageBaseType::PointType;
65 using SpacingType = typename ImageBaseType::SpacingType;
66 using DirectionType = typename ImageBaseType::DirectionType;
67 using SizeType = typename ImageBaseType::SizeType;
68 using RegionType = typename ImageBaseType::RegionType;
70 using VectorOriginType = std::vector<OriginType>;
71 using VectorSpacingType = std::vector<SpacingType>;
72 using VectorDirectionType = std::vector<DirectionType>;
73 using VectorRegionType = std::vector<RegionType>;
74 using VectorGridSpacingFactorType = std::vector<GridSpacingFactorType>;
75
77 using TransformType = Transform<TransformScalarType, Self::Dimension, Self::Dimension>;
78 using TransformPointer = typename TransformType::Pointer;
79 using TransformConstPointer = typename TransformType::ConstPointer;
80
82 itkSetMacro(ImageOrigin, OriginType);
83
85 itkGetConstMacro(ImageOrigin, OriginType);
86
88 itkSetMacro(ImageSpacing, SpacingType);
89
91 itkGetConstMacro(ImageSpacing, SpacingType);
92
94 itkSetMacro(ImageDirection, DirectionType);
95
97 itkGetConstMacro(ImageDirection, DirectionType);
98
100 itkSetMacro(ImageRegion, RegionType);
101
103 itkGetConstMacro(ImageRegion, RegionType);
104
106 itkSetClampMacro(BSplineOrder, unsigned int, 0, 5);
107
109 itkGetConstMacro(BSplineOrder, unsigned int);
110
112 itkSetMacro(FinalGridSpacing, SpacingType);
113
115 itkGetConstMacro(FinalGridSpacing, SpacingType);
116
118 virtual void
119 SetDefaultSchedule(unsigned int levels, double upsamplingFactor);
120
122 virtual void
124
126 virtual void
128
130 itkSetConstObjectMacro(InitialTransform, TransformType);
131
133 virtual void
135
137 virtual void
138 GetBSplineGrid(unsigned int level,
139 RegionType & gridRegion,
140 SpacingType & gridSpacing,
141 OriginType & gridOrigin,
142 DirectionType & gridDirection);
143
144protected:
147
149 ~GridScheduleComputer() override = default;
150
152 VectorSpacingType m_GridSpacings{};
153 VectorOriginType m_GridOrigins{};
154 VectorDirectionType m_GridDirections{};
155 VectorRegionType m_GridRegions{};
156 TransformConstPointer m_InitialTransform{};
157 VectorGridSpacingFactorType m_GridSpacingFactors{};
158
160 void
161 PrintSelf(std::ostream & os, Indent indent) const override;
162
164 itkGetConstMacro(NumberOfLevels, unsigned int);
165
167 virtual void
169 SpacingType & imageSpacing,
170 DirectionType & imageDirection,
171 SpacingType & finalGridSpacing) const;
172
173private:
175 OriginType m_ImageOrigin{};
176 SpacingType m_ImageSpacing{};
177 RegionType m_ImageRegion{};
178 DirectionType m_ImageDirection{};
179 unsigned int m_BSplineOrder{};
180 unsigned int m_NumberOfLevels{};
181 SpacingType m_FinalGridSpacing{};
182
184 itkSetClampMacro(UpsamplingFactor, float, 1.0, NumericTraits<float>::max());
185
187 float m_UpsamplingFactor{};
188};
189
190} // end namespace itk
191
192#ifndef ITK_MANUAL_INSTANTIATION
193# include "itkGridScheduleComputer.hxx"
194#endif
195
196#endif // end #ifndef itkGridScheduleComputer_h
This class computes all information about the B-spline grid, given the image information and the desi...
SmartPointer< const Self > ConstPointer
std::vector< RegionType > VectorRegionType
typename ImageBaseType::RegionType RegionType
Transform< TransformScalarType, Self::Dimension, Self::Dimension > TransformType
typename ImageBaseType::DirectionType DirectionType
typename TransformType::Pointer TransformPointer
std::vector< DirectionType > VectorDirectionType
typename ImageBaseType::SpacingType SpacingType
typename ImageBaseType::PointType PointType
std::vector< OriginType > VectorOriginType
virtual void SetSchedule(const VectorGridSpacingFactorType &schedule)
typename ImageBaseType::SizeType SizeType
~GridScheduleComputer() override=default
virtual void GetBSplineGrid(unsigned int level, RegionType &gridRegion, SpacingType &gridSpacing, OriginType &gridOrigin, DirectionType &gridDirection)
std::vector< GridSpacingFactorType > VectorGridSpacingFactorType
typename TransformType::ConstPointer TransformConstPointer
virtual void ApplyInitialTransform(OriginType &imageOrigin, SpacingType &imageSpacing, DirectionType &imageDirection, SpacingType &finalGridSpacing) const
TTransformScalarType TransformScalarType
virtual void SetDefaultSchedule(unsigned int levels, double upsamplingFactor)
itkStaticConstMacro(Dimension, unsigned int, VImageDimension)
virtual void ComputeBSplineGrid()
ITK_DISALLOW_COPY_AND_MOVE(GridScheduleComputer)
ImageBase< Self::Dimension > ImageBaseType
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void GetSchedule(VectorGridSpacingFactorType &schedule) const
typename ImageBaseType::PointType OriginType
std::vector< SpacingType > VectorSpacingType


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