go home Home | Main Page | Topics | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
Loading...
Searching...
No Matches
elxResamplerBase.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 elxResamplerBase_h
19#define elxResamplerBase_h
20
22#include "elxMacro.h"
23
24#include "elxBaseComponentSE.h"
26#include "itkResampleImageFilter.h"
27#include "itkCastImageFilter.h"
28#include "elxProgressCommand.h"
29
30namespace elastix
31{
78
79template <typename TElastix>
80class ITK_TEMPLATE_EXPORT ResamplerBase : public BaseComponentSE<TElastix>
81{
82public:
84
88
91
93 using typename Superclass::ElastixType;
94 using typename Superclass::RegistrationType;
95
98 using InputImageType = typename ElastixType::MovingImageType;
99 using OutputImageType = typename ElastixType::MovingImageType;
100 // typedef typename ElastixType::FixedImageType OutputImageType;
102
104 using ITKBaseType = itk::ResampleImageFilter<InputImageType, OutputImageType, CoordinateType>;
105
107 using TransformType = typename ITKBaseType::TransformType;
108 using InterpolatorType = typename ITKBaseType::InterpolatorType;
109 using SizeType = typename ITKBaseType::SizeType;
110 using IndexType = typename ITKBaseType::IndexType;
111 using SpacingType = typename ITKBaseType::SpacingType;
112 using DirectionType = typename ITKBaseType::DirectionType;
113 using OriginPointType = typename ITKBaseType::OriginPointType;
114 using OutputPixelType = typename ITKBaseType::PixelType;
115
117 using ParameterMapType = typename ElastixType::ParameterMapType;
118
120 itkStaticConstMacro(ImageDimension, unsigned int, OutputImageType::ImageDimension);
121
125 {
126 return &(this->GetSelf());
127 }
128
129
131 const ITKBaseType *
133 {
134 return &(this->GetSelf());
135 }
136
137
141 virtual int
143 {
144 return 0;
145 }
146
153 void
155
159 void
161
165 void
167
171 void
173
175 virtual void
177
179 void
180 WriteToFile(std::ostream & transformationParameterInfo) const;
181
183 void
185
187 void
188 ResampleAndWriteResultImage(const std::string & filename, const bool showProgress);
189
191 virtual void
193
194protected:
196 ResamplerBase() = default;
198 ~ResamplerBase() override = default;
199
201 virtual void
203
204private:
206
207 virtual ParameterMapType
209 {
210 return {};
211 }
212
214 void
215 WriteResultImage(OutputImageType * imageimage, const std::string & filename, const bool showProgress);
216
218 void
220
222 template <typename TResultPixel>
223 itk::SmartPointer<itk::ImageBase<ImageDimension>>
224 CastImage(const InputImageType & inputImage) const
225 {
226 const auto castFilter =
227 itk::CastImageFilter<InputImageType, itk::Image<TResultPixel, InputImageType::ImageDimension>>::New();
228 castFilter->SetInput(&inputImage);
229 castFilter->Update();
230 return castFilter->GetOutput();
231 }
232
236 template <typename TResultPixel>
237 bool
238 CastImageIfPixelTypesMatch(const std::string_view resultImagePixelType,
239 const InputImageType & inputImage,
240 itk::SmartPointer<itk::ImageBase<ImageDimension>> & outputImage) const
241 {
242 if (resultImagePixelType == PixelTypeToFixedWidthString<TResultPixel>())
243 {
244 outputImage = CastImage<TResultPixel>(inputImage);
245 return true;
246 }
247 return false;
248 }
249
251 template <typename... TResultPixel>
252 itk::SmartPointer<itk::ImageBase<ImageDimension>>
253 CastImageAsSpecifiedByFixedWidthPixelType(const std::string_view resultImagePixelType,
254 const InputImageType & inputImage) const
255 {
256 if (itk::SmartPointer<itk::ImageBase<ImageDimension>> outputImage;
257 (CastImageIfPixelTypesMatch<TResultPixel>(resultImagePixelType, inputImage, outputImage) || ...))
258 {
259 return outputImage;
260 }
261 return nullptr;
262 }
263};
264
265} // end namespace elastix
266
267#ifndef ITK_MANUAL_INSTANTIATION
268# include "elxResamplerBase.hxx"
269#endif
270
271#endif // end #ifndef elxResamplerBase_h
typename ElastixType::RegistrationBaseType RegistrationType
virtual const itk::Object & GetSelf() const =0
void WriteResultImage(OutputImageType *imageimage, const std::string &filename, const bool showProgress)
typename ITKBaseType::IndexType IndexType
void CreateTransformParameterMap(ParameterMapType &parameterMap) const
typename ITKBaseType::SpacingType SpacingType
~ResamplerBase() override=default
void AfterRegistrationBase() override
typename ElastixType::MovingImageType OutputImageType
itk::ResampleImageFilter< InputImageType, OutputImageType, CoordinateType > ITKBaseType
itk::SmartPointer< itk::ImageBase< ImageDimension > > CastImageAsSpecifiedByFixedWidthPixelType(const std::string_view resultImagePixelType, const InputImageType &inputImage) const
itk::SmartPointer< itk::ImageBase< ImageDimension > > CastImage(const InputImageType &inputImage) const
void WriteToFile(std::ostream &transformationParameterInfo) const
ElastixBase::CoordinateType CoordinateType
BaseComponentSE< TElastix > Superclass
bool CastImageIfPixelTypesMatch(const std::string_view resultImagePixelType, const InputImageType &inputImage, itk::SmartPointer< itk::ImageBase< ImageDimension > > &outputImage) const
elxDeclarePureVirtualGetSelfMacro(ITKBaseType)
const ITKBaseType * GetAsITKBaseType() const
typename ElastixType::ParameterMapType ParameterMapType
typename ITKBaseType::TransformType TransformType
typename ITKBaseType::SizeType SizeType
ITKBaseType * GetAsITKBaseType()
void BeforeRegistrationBase() override
typename ITKBaseType::InterpolatorType InterpolatorType
virtual void CreateItkResultImage()
void AfterEachIterationBase() override
itkStaticConstMacro(ImageDimension, unsigned int, OutputImageType::ImageDimension)
virtual ParameterMapType CreateDerivedTransformParameterMap() const
typename ITKBaseType::OriginPointType OriginPointType
typename ITKBaseType::DirectionType DirectionType
ITK_DISALLOW_COPY_AND_MOVE(ResamplerBase)
virtual void SetComponents()
itkOverrideGetNameOfClassMacro(ResamplerBase)
virtual int BeforeAllTransformix()
typename ITKBaseType::PixelType OutputPixelType
void ResampleAndWriteResultImage(const std::string &filename, const bool showProgress)
typename ElastixType::MovingImageType InputImageType
void AfterEachResolutionBase() override
virtual void ReadFromFile()
std::string PixelTypeToFixedWidthString()


Generated on 1774142652 for elastix by doxygen 1.15.0 elastix logo