ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkMRFImageFilter.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 itkMRFImageFilter_h
19#define itkMRFImageFilter_h
20
21#include "vnl/vnl_vector.h"
22#include "vnl/vnl_matrix.h"
23
25
27
30#include "itkSize.h"
31#include "ITKMarkovRandomFieldsClassifiersExport.h"
32
33namespace itk
34{
40{
41public:
47 enum class MRFStop : uint8_t
48 {
51 };
52};
53
54// Define how to print enumeration
55extern ITKMarkovRandomFieldsClassifiers_EXPORT std::ostream &
56 operator<<(std::ostream & out, const MRFImageFilterEnums::MRFStop value);
147template <typename TInputImage, typename TClassifiedImage>
148class ITK_TEMPLATE_EXPORT MRFImageFilter : public ImageToImageFilter<TInputImage, TClassifiedImage>
149{
150public:
151 ITK_DISALLOW_COPY_AND_MOVE(MRFImageFilter);
153
159 using typename Superclass::OutputImagePointer;
160
162 itkNewMacro(Self);
163
165 itkOverrideGetNameOfClassMacro(MRFImageFilter);
166
168 using InputImageType = TInputImage;
169 using InputImagePointer = typename TInputImage::Pointer;
170 using InputImageConstPointer = typename TInputImage::ConstPointer;
171
173 using InputImagePixelType = typename TInputImage::PixelType;
174
176 using InputImageRegionType = typename TInputImage::RegionType;
177
181
183 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
184
186 using TrainingImagePointer = typename TClassifiedImage::Pointer;
187
189 using TrainingImagePixelType = typename TClassifiedImage::PixelType;
190
193 using LabelledImagePointer = typename TClassifiedImage::Pointer;
194
197 using LabelledImagePixelType = typename TClassifiedImage::PixelType;
198
201 using LabelledImageRegionType = typename TClassifiedImage::RegionType;
202
204 using LabelledImageIndexType = typename TClassifiedImage::IndexType;
205 using IndexValueType = typename LabelledImageIndexType::IndexValueType;
206
208 using LabelledImageOffsetType = typename TClassifiedImage::OffsetType;
209
212
214 static constexpr unsigned int ClassifiedImageDimension = TClassifiedImage::ImageDimension;
215
218
220 using SizeType = typename TInputImage::SizeType;
221
223 using NeighborhoodRadiusType = typename TInputImage::SizeType;
224
227
229
231
233
234 using InputImageFaceListIterator = typename InputImageFaceListType::iterator;
235
238
240
242
244
245 using LabelledImageFaceListIterator = typename LabelledImageFaceListType::iterator;
246
248 void
249 SetClassifier(typename ClassifierType::Pointer ptrToClassifier);
250
252 itkSetMacro(NumberOfClasses, unsigned int);
253 itkGetConstMacro(NumberOfClasses, unsigned int);
255
258 itkSetMacro(MaximumNumberOfIterations, unsigned int);
259 itkGetConstMacro(MaximumNumberOfIterations, unsigned int);
261
264 itkSetMacro(ErrorTolerance, double);
265 itkGetConstMacro(ErrorTolerance, double);
267
270 itkSetMacro(SmoothingFactor, double);
271 itkGetConstMacro(SmoothingFactor, double);
273
275 void
277
282 void
284
285 void
287
291 {
293
294 for (int i = 0; i < InputImageDimension; ++i)
295 {
296 radius[i] = m_InputImageNeighborhoodRadius[i];
297 }
298 return radius;
299 }
300
306 virtual void
307 SetMRFNeighborhoodWeight(std::vector<double> betaMatrix);
308
309 virtual std::vector<double>
314
316#if !defined(ITK_LEGACY_REMOVE)
318 // We need to expose the enum values at the class level
319 // for backwards compatibility
320 static constexpr MRFStopEnum MaximumNumberOfIterations = MRFStopEnum::MaximumNumberOfIterations;
321 static constexpr MRFStopEnum ErrorTolerance = MRFStopEnum::ErrorTolerance;
322#endif
323
326 itkGetConstReferenceMacro(StopCondition, MRFStopEnum);
327
328 /* Get macro for number of iterations */
329 itkGetConstReferenceMacro(NumberOfIterations, unsigned int);
330
331 itkConceptMacro(UnsignedIntConvertibleToClassifiedCheck,
333 itkConceptMacro(ClassifiedConvertibleToUnsignedIntCheck,
335 itkConceptMacro(ClassifiedConvertibleToIntCheck, (Concept::Convertible<LabelledImagePixelType, int>));
336 itkConceptMacro(IntConvertibleToClassifiedCheck, (Concept::Convertible<int, LabelledImagePixelType>));
338
339protected:
341 ~MRFImageFilter() override = default;
342 void
343 PrintSelf(std::ostream & os, Indent indent) const override;
344
346 void
348
351 virtual void
353
355 virtual void
357
363
366
368 virtual void
371 LabelStatusImageNeighborhoodIterator & labelStatusIter);
372
373 void
374 GenerateData() override;
375
376 void
378
379 void
381
382 void
384
385private:
386 using InputImageSizeType = typename TInputImage::SizeType;
387
389
391
393
394 using LabelStatusImageFaceListIterator = typename LabelStatusImageFaceListType::iterator;
395
399
400 unsigned int m_NumberOfClasses{ 0 };
401 unsigned int m_MaximumNumberOfIterations{ 50 };
402 unsigned int m_KernelSize{};
403
408 double m_ErrorTolerance{ 0.2 };
409 double m_SmoothingFactor{ 1 };
410 double * m_ClassProbability{ nullptr }; // Class likelihood
411 unsigned int m_NumberOfIterations{ 0 };
412 MRFStopEnum m_StopCondition{ MRFStopEnum::MaximumNumberOfIterations };
413
415
416 std::vector<double> m_MRFNeighborhoodWeight{};
417 std::vector<double> m_NeighborInfluence{};
418 std::vector<double> m_MahalanobisDistance{};
419 std::vector<double> m_DummyVector{};
420
423
427 virtual void
429
431 void
433}; // class MRFImageFilter
434} // namespace itk
435
436#ifndef ITK_MANUAL_INSTANTIATION
437# include "itkMRFImageFilter.hxx"
438#endif
439
440#endif
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Base class for all data objects in ITK.
Base class for the ImageClassifierBase object.
A multi-dimensional iterator templated over image type that walks a region of pixels.
A multi-dimensional iterator templated over image type that walks a region of pixels.
typename OutputImageType::Pointer OutputImagePointer
Templated n-dimensional image class.
Definition itkImage.h:89
ImageRegion< VImageDimension > RegionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Contains all enum classes in MRFImageFilter class;.
ImageRegionIterator< TClassifiedImage > LabelledImageRegionIterator
std::vector< double > m_MRFNeighborhoodWeight
ImageRegionConstIterator< TInputImage > InputImageRegionConstIterator
virtual void SetMRFNeighborhoodWeight(std::vector< double > betaMatrix)
NeighborhoodIterator< LabelStatusImageType > LabelStatusImageNeighborhoodIterator
LabelledImageNeighborhoodRadiusType m_LabelledImageNeighborhoodRadius
ImageToImageFilter< TInputImage, TClassifiedImage > Superclass
typename TClassifiedImage::OffsetType LabelledImageOffsetType
SmartPointer< Self > Pointer
typename LabelledImageIndexType::IndexValueType IndexValueType
void GenerateData() override
virtual void DoNeighborhoodOperation(const InputImageNeighborhoodIterator &imageIter, LabelledImageNeighborhoodIterator &labelledIter, LabelStatusImageNeighborhoodIterator &labelStatusIter)
ClassifierType::Pointer m_ClassifierPtr
unsigned int m_MaximumNumberOfIterations
void SetNeighborhoodRadius(const SizeValueType)
ImageRegionIterator< TInputImage > InputImageRegionIterator
typename InputImageNeighborhoodIterator::RadiusType InputImageNeighborhoodRadiusType
void SetClassifier(typename ClassifierType::Pointer ptrToClassifier)
std::vector< double > m_MahalanobisDistance
typename TInputImage::RegionType InputImageRegionType
typename TClassifiedImage::Pointer TrainingImagePointer
void SetNeighborhoodRadius(const NeighborhoodRadiusType &)
LabelStatusImagePointer m_LabelStatusImage
LabelStatusImageNeighborhoodRadiusType m_LabelStatusImageNeighborhoodRadius
void GenerateOutputInformation() override
typename InputImageFacesCalculator::FaceListType InputImageFaceListType
typename TInputImage::SizeType SizeType
ImageRegionIterator< LabelStatusImageType > LabelStatusImageIterator
InputImageNeighborhoodRadiusType m_InputImageNeighborhoodRadius
unsigned int m_NumberOfIterations
virtual void MinimizeFunctional()
const NeighborhoodRadiusType GetNeighborhoodRadius() const
typename TInputImage::ConstPointer InputImageConstPointer
NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TClassifiedImage > LabelledImageFacesCalculator
static constexpr unsigned int ClassifiedImageDimension
typename InputImageFaceListType::iterator InputImageFaceListIterator
typename LabelStatusImageType::RegionType LabelStatusRegionType
virtual void SetDefaultMRFNeighborhoodWeight()
typename TInputImage::SizeType InputImageSizeType
SmartPointer< const Self > ConstPointer
void GenerateInputRequestedRegion() override
Image< int, Self::InputImageDimension > LabelStatusImageType
typename TClassifiedImage::PixelType LabelledImagePixelType
typename TInputImage::SizeType NeighborhoodRadiusType
NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TInputImage > InputImageFacesCalculator
void PrintSelf(std::ostream &os, Indent indent) const override
static constexpr unsigned int InputImageDimension
NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< LabelStatusImageType > LabelStatusImageFacesCalculator
typename LabelStatusImageFacesCalculator::FaceListType LabelStatusImageFaceListType
NeighborhoodIterator< TClassifiedImage > LabelledImageNeighborhoodIterator
void SetNeighborhoodRadius(const SizeValueType *radiusArray)
typename TClassifiedImage::IndexType LabelledImageIndexType
~MRFImageFilter() override=default
typename TInputImage::Pointer InputImagePointer
typename TClassifiedImage::Pointer LabelledImagePointer
unsigned int m_NumberOfClasses
typename TClassifiedImage::RegionType LabelledImageRegionType
MRFImageFilterEnums::MRFStop MRFStopEnum
std::vector< double > m_DummyVector
typename TClassifiedImage::PixelType TrainingImagePixelType
typename LabelledImageFaceListType::iterator LabelledImageFaceListIterator
void EnlargeOutputRequestedRegion(DataObject *) override
typename TInputImage::PixelType InputImagePixelType
virtual void ApplyMRFImageFilter()
typename LabelStatusImageFaceListType::iterator LabelStatusImageFaceListIterator
typename LabelStatusImageNeighborhoodIterator::RadiusType LabelStatusImageNeighborhoodRadiusType
typename LabelStatusImageType::IndexType LabelStatusIndexType
ImageClassifierBase< TInputImage, TClassifiedImage > ClassifierType
typename LabelledImageFacesCalculator::FaceListType LabelledImageFaceListType
std::vector< double > m_NeighborInfluence
typename LabelStatusImageType::Pointer LabelStatusImagePointer
typename LabelledImageNeighborhoodIterator::RadiusType LabelledImageNeighborhoodRadiusType
ConstNeighborhoodIterator< TInputImage > InputImageNeighborhoodIterator
virtual std::vector< double > GetMRFNeighborhoodWeight()
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
Splits an image into a main region and several "face" regions which are used to handle computations o...