ITK  6.0.0
Insight Toolkit
itkKLMRegionGrowImageFilter.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 itkKLMRegionGrowImageFilter_h
19#define itkKLMRegionGrowImageFilter_h
20
21#include "itkImage.h"
25#include "itkConceptChecking.h"
26#include <algorithm>
27#include <functional>
28
29namespace itk
30{
166template <typename TInputImage, typename TOutputImage>
167class ITK_TEMPLATE_EXPORT KLMRegionGrowImageFilter : public RegionGrowImageFilter<TInputImage, TOutputImage>
168{
169public:
170 ITK_DISALLOW_COPY_AND_MOVE(KLMRegionGrowImageFilter);
171
177
179 itkNewMacro(Self);
180
182 itkOverrideGetNameOfClassMacro(KLMRegionGrowImageFilter);
183
185 using InputImageType = TInputImage;
188
190 using InputImagePixelType = typename TInputImage::PixelType;
191
194
196 static constexpr unsigned int InputImageVectorDimension = InputImagePixelType::Dimension;
197
200
204
207
210 using typename Superclass::GridSizeType;
211
213 using OutputImageType = TOutputImage;
215
217 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
218
220 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
221
223 using OutputImagePixelType = typename TOutputImage::PixelType;
224
227
229 static constexpr unsigned int OutputImageVectorDimension = OutputImagePixelType::Dimension;
230
233
236
239
241 static constexpr RegionLabelType LabelImageDimension = InputImageDimension;
242
245
248
251
254
257
259 using MeanRegionIntensityType = vnl_vector<double>;
260
263
266
270 itkSetMacro(MaximumLambda, double);
271 itkGetConstReferenceMacro(MaximumLambda, double);
275 itkSetMacro(NumberOfRegions, unsigned int);
276 itkGetConstReferenceMacro(NumberOfRegions, unsigned int);
282
284 void
286
288 void
290
291#ifdef ITK_USE_CONCEPT_CHECKING
292 // Begin concept checking
295# if defined(THIS_CONCEPT_FAILS_ON_GCC)
297 itkConceptMacro(SameVectorDimension,
299# endif
300 // End concept checking
301#endif
302
303protected:
305 ~KLMRegionGrowImageFilter() override = default;
306 void
307 PrintSelf(std::ostream & os, Indent indent) const override;
308
312 void
313 GenerateData() override;
314
318 void
320
325 void
327
330 void
332
336 void
337 MergeRegions() override;
338
340 virtual void
342
344 void
346
348 void
350
354
360 virtual void
362
366 virtual void
368
369private:
373
374 double m_MaximumLambda{ 1000 };
375 unsigned int m_NumberOfRegions{ 0 };
376
379 double m_InternalLambda{ 0 };
380 unsigned int m_InitialNumberOfRegions{ 0 };
381 double m_TotalBorderLength{ 0.0 };
382
383 std::vector<KLMSegmentationRegionPtr> m_RegionsPointer{};
384 std::vector<KLMSegmentationBorderPtr> m_BordersPointer{};
385 std::vector<KLMSegmentationBorderArrayPtr> m_BordersDynamicPointer{};
386 KLMSegmentationBorderArrayPtr * m_BorderCandidate{ nullptr };
387
388 MeanRegionIntensityType m_InitialRegionMean{};
389 double m_InitialRegionArea{ 0 };
390}; // class KLMRegionGrowImageFilter
391} // namespace itk
392
393#ifndef ITK_MANUAL_INSTANTIATION
394# include "itkKLMRegionGrowImageFilter.hxx"
395#endif
396
397#endif
Base class for all data objects in ITK.
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.
Base class for all process objects that output image data.
typename OutputImageType::PixelType OutputImagePixelType
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::Pointer InputImagePointer
typename InputImageType::PixelType InputImagePixelType
Templated n-dimensional image class.
Definition: itkImage.h:89
TPixel PixelType
Definition: itkImage.h:108
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Object maintaining a reference to a list of borders associated with a region.
Base class for a region growing object that performs energy-based region growing for multiband images...
void ApplyRegionGrowImageFilter() override
typename TOutputImage::IndexType OutputImageIndexType
typename KLMSegmentationRegion::RegionLabelType RegionLabelType
typename TInputImage::IndexType InputImageIndexType
virtual void GenerateOutputImage()
void GenerateInputRequestedRegion() override
LabelImagePointer GenerateLabelledImage(LabelImageType *labelImagePtr)
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void InitializeRegionParameters(InputRegionType region)
~KLMRegionGrowImageFilter() override=default
void EnlargeOutputRequestedRegion(DataObject *) override
typename KLMSegmentationBorder::Pointer KLMSegmentationBorderPtr
typename TInputImage::PixelType::VectorType InputImageVectorType
typename TInputImage::SizeType InputImageSizeType
typename LabelImageType::Pointer LabelImagePointer
typename TInputImage::RegionType InputRegionType
typename LabelImageType::IndexType LabelImageIndexType
typename LabelImageType::PixelType LabelImagePixelType
typename TOutputImage::PixelType::VectorType OutputImageVectorType
LabelImagePointer GetLabelledImage()
typename KLMSegmentationRegion::Pointer KLMSegmentationRegionPtr
Base class for KLMSegmentationBorder object.
Superclass::RegionLabelType RegionLabelType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for RegionGrowImageFilter object.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
class ITK_FORWARD_EXPORT KLMSegmentationBorder