ITK  6.0.0
Insight Toolkit
itkLevelSetDomainMapImageFilter.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 
19 #ifndef itkLevelSetDomainMapImageFilter_h
20 #define itkLevelSetDomainMapImageFilter_h
21 
22 #include "itkImageToImageFilter.h"
23 #include "itkImageRegionIterator.h"
25 #include <list>
26 #include <utility>
27 #include <vector>
28 
29 namespace itk
30 {
39 template <typename TInputImage, typename TOutputImage>
40 class ITK_TEMPLATE_EXPORT LevelSetDomainMapImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
41 {
42 public:
43  ITK_DISALLOW_COPY_AND_MOVE(LevelSetDomainMapImageFilter);
44 
49 
50  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
51 
53  itkNewMacro(Self);
54 
56  itkOverrideGetNameOfClassMacro(LevelSetDomainMapImageFilter);
57 
58  using InputImageType = TInputImage;
60  using InputImagePixelType = typename InputImageType::PixelType;
65 
66  using OutputImageType = TOutputImage;
69  using OutputImagePixelType = typename OutputImageType::PixelType;
70 
74 
78 
84  {
85  public:
86  LevelSetDomain() = default;
87 
89  : m_Region(reg)
90  , m_IdList(std::move(iList))
91  {}
92 
93  const InputImageRegionType *
94  GetRegion() const
95  {
96  return &(this->m_Region);
97  }
98 
99  const InputImagePixelType *
100  GetIdList() const
101  {
102  return &(this->m_IdList);
103  }
104 
105  private:
108  };
109 
111  using DomainMapType = std::map<IdentifierType, LevelSetDomain>;
112 
116  const DomainMapType &
117  GetDomainMap() const;
118 
119 protected:
121  ~LevelSetDomainMapImageFilter() override = default;
122 
126  ComputeConsistentRegion(const InputImageRegionType & inputRegion) const;
127 
130  void
131  GenerateData() override;
132 
134  void
135  PrintSelf(std::ostream & os, Indent indent) const override;
136 
137 private:
138  DomainMapType m_DomainMap{};
139 
140  const InputImageType * m_InputImage{};
141  OutputImageType * m_OutputImage{};
142 };
143 
144 } /* namespace itk */
145 
146 #ifndef ITK_MANUAL_INSTANTIATION
147 # include "itkLevelSetDomainMapImageFilter.hxx"
148 #endif
149 
150 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::LevelSetDomainMapImageFilter::OutputImageType
TOutputImage OutputImageType
Definition: itkLevelSetDomainMapImageFilter.h:66
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::ImageRegionConstIteratorWithIndex
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
Definition: itkImageRegionConstIteratorWithIndex.h:130
itk::LevelSetDomainMapImageFilter::OutputImageIndexType
typename OutputImageType::IndexType OutputImageIndexType
Definition: itkLevelSetDomainMapImageFilter.h:68
itk::LevelSetDomainMapImageFilter::LevelSetDomain::GetIdList
const InputImagePixelType * GetIdList() const
Definition: itkLevelSetDomainMapImageFilter.h:100
itk::LevelSetDomainMapImageFilter::LevelSetDomain
Specifies an image region where an unique std::list of level sets Id's are defined.
Definition: itkLevelSetDomainMapImageFilter.h:83
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::LevelSetDomainMapImageFilter::InputImageSizeType
typename InputImageType::SizeType InputImageSizeType
Definition: itkLevelSetDomainMapImageFilter.h:62
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::LevelSetDomainMapImageFilter::InputImageType
TInputImage InputImageType
Definition: itkLevelSetDomainMapImageFilter.h:58
itkImageRegionIteratorWithIndex.h
itk::LevelSetDomainMapImageFilter::LevelSetDomain::GetRegion
const InputImageRegionType * GetRegion() const
Definition: itkLevelSetDomainMapImageFilter.h:94
itk::ImageToImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkImageToImageFilter.h:133
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::LevelSetDomainMapImageFilter::InputImageSizeValueType
typename InputImageSizeType::SizeValueType InputImageSizeValueType
Definition: itkLevelSetDomainMapImageFilter.h:63
itkImageToImageFilter.h
itk::ImageRegionIteratorWithIndex
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Definition: itkImageRegionIteratorWithIndex.h:73
itk::LevelSetDomainMapImageFilter::LevelSetDomain::LevelSetDomain
LevelSetDomain(const InputImageRegionType &reg, InputImagePixelType iList)
Definition: itkLevelSetDomainMapImageFilter.h:88
itk::LevelSetDomainMapImageFilter::LevelSetDomain::m_IdList
InputImagePixelType m_IdList
Definition: itkLevelSetDomainMapImageFilter.h:107
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::LevelSetDomainMapImageFilter::LevelSetDomain::m_Region
InputImageRegionType m_Region
Definition: itkLevelSetDomainMapImageFilter.h:106
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::LevelSetDomainMapImageFilter::InputImageIndexType
typename InputImageType::IndexType InputImageIndexType
Definition: itkLevelSetDomainMapImageFilter.h:64
itk::LevelSetDomainMapImageFilter
Definition: itkLevelSetDomainMapImageFilter.h:40
itk::LevelSetDomainMapImageFilter::DomainMapType
std::map< IdentifierType, LevelSetDomain > DomainMapType
Definition: itkLevelSetDomainMapImageFilter.h:111
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86