go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkGenericMultiResolutionPyramidImageFilter.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 itkGenericMultiResolutionPyramidImageFilter_h
19#define itkGenericMultiResolutionPyramidImageFilter_h
20
21#include "itkMultiResolutionPyramidImageFilter.h"
22#include "itkSmoothingRecursiveGaussianImageFilter.h"
23
24namespace itk
25{
116template <class TInputImage, class TOutputImage, class TPrecisionType = double>
118 : public MultiResolutionPyramidImageFilter<TInputImage, TOutputImage>
119{
120public:
122
125 using Superclass = MultiResolutionPyramidImageFilter<TInputImage, TOutputImage>;
126 using SuperSuperclass = typename Superclass::Superclass;
127 using Pointer = SmartPointer<Self>;
128 using ConstPointer = SmartPointer<const Self>;
129
131 itkNewMacro(Self);
132
134 itkTypeMacro(GenericMultiResolutionPyramidImageFilter, MultiResolutionPyramidImageFilter);
135
137 itkStaticConstMacro(ImageDimension, unsigned int, TInputImage::ImageDimension);
138 itkStaticConstMacro(OutputImageDimension, unsigned int, TOutputImage::ImageDimension);
139
141 using typename Superclass::ScheduleType;
142 using typename Superclass::InputImageType;
143 using typename Superclass::OutputImageType;
144 using typename Superclass::InputImagePointer;
145 using typename Superclass::OutputImagePointer;
146 using typename Superclass::InputImageConstPointer;
147 using SpacingType = typename Superclass::InputImageType::SpacingType;
148 using PixelType = typename InputImageType::PixelType;
149 using ScalarRealType = typename NumericTraits<PixelType>::ScalarRealType;
150
152 using SmoothingScheduleType = vnl_matrix<ScalarRealType>;
153 using RescaleScheduleType = ScheduleType;
154
156 using SigmaArrayType = FixedArray<ScalarRealType, Self::ImageDimension>;
158
166 void
167 SetSchedule(const ScheduleType & schedule) override;
168
176 virtual void
178
180 virtual void
182
184 const RescaleScheduleType &
186 {
187 return this->m_Schedule;
188 }
189
190
195 virtual void
197
199 virtual void
201
203 itkGetConstReferenceMacro(SmoothingSchedule, SmoothingScheduleType);
204
209 void
210 SetNumberOfLevels(unsigned int num) override;
211
215 virtual void
216 SetCurrentLevel(unsigned int level);
217
219 itkGetConstReferenceMacro(CurrentLevel, unsigned int);
220
222 virtual void
224
225 itkGetConstMacro(ComputeOnlyForCurrentLevel, bool);
226 itkBooleanMacro(ComputeOnlyForCurrentLevel);
227
228#ifdef ITK_USE_CONCEPT_CHECKING
230 itkConceptMacro(SameDimensionCheck, (Concept::SameDimension<ImageDimension, OutputImageDimension>));
231 itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<typename TOutputImage::PixelType>));
233#endif
234
235protected:
238
240 void
241 PrintSelf(std::ostream & os, Indent indent) const override;
242
250 void
252
258 void
259 GenerateOutputRequestedRegion(DataObject * output) override;
260
262 void
264
266 void
267 GenerateData() override;
268
270 void
272
273 SmoothingScheduleType m_SmoothingSchedule{};
274 unsigned int m_CurrentLevel{};
275 bool m_ComputeOnlyForCurrentLevel{};
276 bool m_SmoothingScheduleDefined{};
277
278private:
282 using SmootherType = SmoothingRecursiveGaussianImageFilter<InputImageType, OutputImageType>;
283
288 using ImageToImageFilterSameTypes = ImageToImageFilter<OutputImageType, OutputImageType>;
289 using ImageToImageFilterDifferentTypes = ImageToImageFilter<InputImageType, OutputImageType>;
290
294 bool
295 SetupSmoother(const unsigned int level,
296 typename SmootherType::Pointer & smoother,
297 const InputImageConstPointer & input);
298
302 int
303 SetupShrinkerOrResampler(const unsigned int level,
304 typename SmootherType::Pointer & smoother,
305 const bool sameType,
306 const InputImageConstPointer & input,
307 const OutputImagePointer & outputPtr,
308 typename ImageToImageFilterSameTypes::Pointer & rescaleSameTypes,
309 typename ImageToImageFilterDifferentTypes::Pointer & rescaleDifferentTypes);
310
312 void
313 DefineShrinkerOrResampler(const bool sameType,
314 const RescaleFactorArrayType & shrinkFactors,
315 const OutputImagePointer & outputPtr,
316 typename ImageToImageFilterSameTypes::Pointer & rescaleSameTypes,
317 typename ImageToImageFilterDifferentTypes::Pointer & rescaleDifferentTypes);
318
320 void
322
326 bool
327 ComputeForCurrentLevel(const unsigned int level) const;
328
330 double
331 GetDefaultSigma(const unsigned int level,
332 const unsigned int dim,
333 const unsigned int * factors,
334 const SpacingType & spacing) const;
335
337 void
338 GetSigma(const unsigned int level, SigmaArrayType & sigmaArray) const;
339
341 void
342 GetShrinkFactors(const unsigned int level, RescaleFactorArrayType & shrinkFactors) const;
343
347 bool
348 AreSigmasAllZeros(const SigmaArrayType & sigmaArray) const;
349
353 bool
354 AreRescaleFactorsAllOnes(const RescaleFactorArrayType & rescaleFactorArray) const;
355
357 bool
359
361 bool
363};
364
365} // namespace itk
366
367#ifndef ITK_MANUAL_INSTANTIATION
368# include "itkGenericMultiResolutionPyramidImageFilter.hxx"
369#endif
370
371#endif
Framework for creating images in a multi-resolution pyramid.
~GenericMultiResolutionPyramidImageFilter() override=default
ImageToImageFilter< OutputImageType, OutputImageType > ImageToImageFilterSameTypes
virtual void SetComputeOnlyForCurrentLevel(const bool _arg)
bool AreRescaleFactorsAllOnes(const RescaleFactorArrayType &rescaleFactorArray) const
FixedArray< ScalarRealType, Self::ImageDimension > SigmaArrayType
virtual void SetSmoothingSchedule(const SmoothingScheduleType &schedule)
MultiResolutionPyramidImageFilter< TInputImage, TOutputImage > Superclass
ITK_DISALLOW_COPY_AND_MOVE(GenericMultiResolutionPyramidImageFilter)
virtual void SetRescaleSchedule(const RescaleScheduleType &schedule)
bool ComputeForCurrentLevel(const unsigned int level) const
double GetDefaultSigma(const unsigned int level, const unsigned int dim, const unsigned int *factors, const SpacingType &spacing) const
void GetShrinkFactors(const unsigned int level, RescaleFactorArrayType &shrinkFactors) const
int SetupShrinkerOrResampler(const unsigned int level, typename SmootherType::Pointer &smoother, const bool sameType, const InputImageConstPointer &input, const OutputImagePointer &outputPtr, typename ImageToImageFilterSameTypes::Pointer &rescaleSameTypes, typename ImageToImageFilterDifferentTypes::Pointer &rescaleDifferentTypes)
itkStaticConstMacro(ImageDimension, unsigned int, TInputImage::ImageDimension)
virtual void SetCurrentLevel(unsigned int level)
SmoothingRecursiveGaussianImageFilter< InputImageType, OutputImageType > SmootherType
void GenerateOutputRequestedRegion(DataObject *output) override
itkStaticConstMacro(OutputImageDimension, unsigned int, TOutputImage::ImageDimension)
void GetSigma(const unsigned int level, SigmaArrayType &sigmaArray) const
ImageToImageFilter< InputImageType, OutputImageType > ImageToImageFilterDifferentTypes
void SetSchedule(const ScheduleType &schedule) override
bool AreSigmasAllZeros(const SigmaArrayType &sigmaArray) const
typename NumericTraits< PixelType >::ScalarRealType ScalarRealType
bool SetupSmoother(const unsigned int level, typename SmootherType::Pointer &smoother, const InputImageConstPointer &input)
void DefineShrinkerOrResampler(const bool sameType, const RescaleFactorArrayType &shrinkFactors, const OutputImagePointer &outputPtr, typename ImageToImageFilterSameTypes::Pointer &rescaleSameTypes, typename ImageToImageFilterDifferentTypes::Pointer &rescaleDifferentTypes)
void PrintSelf(std::ostream &os, Indent indent) const override
void SetNumberOfLevels(unsigned int num) override


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