ITK  5.2.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  * 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 itkImageRegistrationMethodImageSource_h
19 #define itkImageRegistrationMethodImageSource_h
20 #include "itkImage.h"
23 #include "itkOptimizerParameters.h"
24 
34 namespace itk
35 {
36 
37 namespace testhelper
38 {
39 
40 template <typename TFixedPixelType, typename TMovingPixelType, unsigned int NDimension>
42 {
43 public:
45  using Superclass = Object;
49 
51  itkNewMacro(Self);
52 
54  itkTypeMacro(Image, Object);
55 
56 
59 
60  const MovingImageType *
62  {
63  return m_MovingImage.GetPointer();
64  }
65 
66  const FixedImageType *
67  GetFixedImage() const
68  {
69  return m_FixedImage.GetPointer();
70  }
71 
72  const ParametersType &
74  {
75  return m_Parameters;
76  }
77 
78 
79  void
81  {
82  typename MovingImageType::IndexType index;
83  index.Fill(0);
84  typename MovingImageType::RegionType region;
85  region.SetSize(size);
86  region.SetIndex(index);
87 
88  m_MovingImage->SetLargestPossibleRegion(region);
89  m_MovingImage->SetBufferedRegion(region);
90  m_MovingImage->SetRequestedRegion(region);
91  m_MovingImage->Allocate();
92 
93  m_FixedImage->SetLargestPossibleRegion(region);
94  m_FixedImage->SetBufferedRegion(region);
95  m_FixedImage->SetRequestedRegion(region);
96  m_FixedImage->Allocate();
97 
98  /* Fill images with a 2D gaussian*/
99  using MovingImageIteratorType = itk::ImageRegionIteratorWithIndex<MovingImageType>;
100 
101  using FixedImageIteratorType = itk::ImageRegionIteratorWithIndex<FixedImageType>;
102 
103 
104  itk::Point<double, 2> center;
105  center[0] = (double)region.GetSize()[0] / 2.0;
106  center[1] = (double)region.GetSize()[1] / 2.0;
107 
108  const double s = (double)region.GetSize()[0] / 2.0;
109 
112 
113  /* Set the displacement */
114  itk::Vector<double, 2> displacement;
115  displacement[0] = m_Parameters[0];
116  displacement[1] = m_Parameters[1];
117 
118 
119  MovingImageIteratorType ri(m_MovingImage, region);
120  FixedImageIteratorType ti(m_FixedImage, region);
121  while (!ri.IsAtEnd())
122  {
123  p[0] = ri.GetIndex()[0];
124  p[1] = ri.GetIndex()[1];
125  d = p - center;
126  d += displacement;
127  const double x = d[0];
128  const double y = d[1];
129  const double value = 200.0 * std::exp(-(x * x + y * y) / (s * s));
130  ri.Set(static_cast<typename MovingImageType::PixelType>(value));
131  ++ri;
132  }
133 
134 
135  while (!ti.IsAtEnd())
136  {
137  p[0] = ti.GetIndex()[0];
138  p[1] = ti.GetIndex()[1];
139  d = p - center;
140  const double x = d[0];
141  const double y = d[1];
142  const double value = 200.0 * std::exp(-(x * x + y * y) / (s * s));
143  ti.Set(static_cast<typename FixedImageType::PixelType>(value));
144  ++ti;
145  }
146  }
147 
148 protected:
150  {
154  m_Parameters[0] = 7.0;
155  m_Parameters[1] = 3.0;
156  }
157 
158 private:
161 
163 };
164 
165 } // end namespace testhelper
166 
167 } // end namespace itk
168 #endif
itk::OptimizerParameters< double >
itk::testhelper::ImageRegistrationMethodImageSource::FixedImageType
itk::Image< TFixedPixelType, NDimension > FixedImageType
Definition: itkImageRegistrationMethodImageSource.h:58
itkCommandIterationUpdate.h
itk::testhelper::ImageRegistrationMethodImageSource::GenerateImages
void GenerateImages(const typename MovingImageType::SizeType &size)
Definition: itkImageRegistrationMethodImageSource.h:80
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::testhelper::ImageRegistrationMethodImageSource::m_MovingImage
MovingImageType::Pointer m_MovingImage
Definition: itkImageRegistrationMethodImageSource.h:160
itkImage.h
itk::SmartPointer< Self >
itkImageRegionIteratorWithIndex.h
itk::testhelper::ImageRegistrationMethodImageSource::m_FixedImage
FixedImageType::Pointer m_FixedImage
Definition: itkImageRegistrationMethodImageSource.h:159
itk::testhelper::ImageRegistrationMethodImageSource::GetActualParameters
const ParametersType & GetActualParameters() const
Definition: itkImageRegistrationMethodImageSource.h:73
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itkOptimizerParameters.h
itk::testhelper::ImageRegistrationMethodImageSource::ImageRegistrationMethodImageSource
ImageRegistrationMethodImageSource()
Definition: itkImageRegistrationMethodImageSource.h:149
itk::Image::RegionType
typename Superclass::RegionType RegionType
Definition: itkImage.h:150
itk::SmartPointer::GetPointer
ObjectType * GetPointer() const noexcept
Definition: itkSmartPointer.h:132
itk::testhelper::ImageRegistrationMethodImageSource::GetFixedImage
const FixedImageType * GetFixedImage() const
Definition: itkImageRegistrationMethodImageSource.h:67
itk::testhelper::ImageRegistrationMethodImageSource
Definition: itkImageRegistrationMethodImageSource.h:41
itk::ImageRegionIteratorWithIndex
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Definition: itkImageRegionIteratorWithIndex.h:73
itk::Image::IndexType
typename Superclass::IndexType IndexType
Definition: itkImage.h:132
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::Image::New
static Pointer New()
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
itk::testhelper::ImageRegistrationMethodImageSource::ParametersType
OptimizerParameters< double > ParametersType
Definition: itkImageRegistrationMethodImageSource.h:48
itk::Point< double, 2 >
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::Object::Object
Object()
itk::Image::SizeType
typename Superclass::SizeType SizeType
Definition: itkImage.h:139
itk::testhelper::ImageRegistrationMethodImageSource::m_Parameters
ParametersType m_Parameters
Definition: itkImageRegistrationMethodImageSource.h:162
itk::testhelper::ImageRegistrationMethodImageSource::GetMovingImage
const MovingImageType * GetMovingImage() const
Definition: itkImageRegistrationMethodImageSource.h:61