ITK  5.4.0
Insight Toolkit
itkImageRegistrationMethodImageSource.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 itkImageRegistrationMethodImageSource_h
19#define itkImageRegistrationMethodImageSource_h
20#include "itkImage.h"
24
34namespace itk
35{
36
37namespace testhelper
38{
39
40template <typename TFixedPixelType, typename TMovingPixelType, unsigned int VDimension>
42{
43public:
49
51 itkNewMacro(Self);
52
54 itkOverrideGetNameOfClassMacro(ImageRegistrationMethodImageSource);
55
56
59
60 const MovingImageType *
62 {
64 }
65
66 const FixedImageType *
68 {
69 return m_FixedImage.GetPointer();
70 }
71
72 const ParametersType &
74 {
75 return m_Parameters;
76 }
77
78
79 void
81 {
82 const typename MovingImageType::RegionType region(size);
83
84 m_MovingImage->SetRegions(region);
85 m_MovingImage->Allocate();
86
87 m_FixedImage->SetRegions(region);
88 m_FixedImage->Allocate();
89
90 /* Fill images with a 2D gaussian*/
91 using MovingImageIteratorType = itk::ImageRegionIteratorWithIndex<MovingImageType>;
92
93 using FixedImageIteratorType = itk::ImageRegionIteratorWithIndex<FixedImageType>;
94
95
97 center[0] = static_cast<double>(region.GetSize()[0]) / 2.0;
98 center[1] = static_cast<double>(region.GetSize()[1]) / 2.0;
99
100 const double s = static_cast<double>(region.GetSize()[0]) / 2.0;
101
104
105 /* Set the displacement */
106 itk::Vector<double, 2> displacement;
107 displacement[0] = m_Parameters[0];
108 displacement[1] = m_Parameters[1];
109
110
111 MovingImageIteratorType ri(m_MovingImage, region);
112 FixedImageIteratorType ti(m_FixedImage, region);
113 while (!ri.IsAtEnd())
114 {
115 p[0] = ri.GetIndex()[0];
116 p[1] = ri.GetIndex()[1];
117 d = p - center;
118 d += displacement;
119 const double x = d[0];
120 const double y = d[1];
121 const double value = 200.0 * std::exp(-(x * x + y * y) / (s * s));
122 ri.Set(static_cast<typename MovingImageType::PixelType>(value));
123 ++ri;
124 }
125
126
127 while (!ti.IsAtEnd())
128 {
129 p[0] = ti.GetIndex()[0];
130 p[1] = ti.GetIndex()[1];
131 d = p - center;
132 const double x = d[0];
133 const double y = d[1];
134 const double value = 200.0 * std::exp(-(x * x + y * y) / (s * s));
135 ti.Set(static_cast<typename FixedImageType::PixelType>(value));
136 ++ti;
137 }
138 }
139
140protected:
142 {
146 m_Parameters[0] = 7.0;
147 m_Parameters[1] = 3.0;
148 }
149
150private:
153
155};
156
157} // end namespace testhelper
158
159} // end namespace itk
160#endif
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
An image region represents a structured region of data.
const SizeType & GetSize() const
Templated n-dimensional image class.
Definition: itkImage.h:89
static Pointer New()
TPixel PixelType
Definition: itkImage.h:108
Light weight base class for most itk classes.
Base class for most ITK classes.
Definition: itkObject.h:62
ObjectType * GetPointer() const noexcept
void GenerateImages(const typename MovingImageType::SizeType &size)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....