ITK  6.0.0
Insight Toolkit
itkLevelSetDomainPartitionImage.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 itkLevelSetDomainPartitionImage_h
19#define itkLevelSetDomainPartitionImage_h
20
23
24namespace itk
25{
32template <typename TImage>
33class ITK_TEMPLATE_EXPORT LevelSetDomainPartitionImage : public LevelSetDomainPartitionBase<TImage>
34{
35public:
36 ITK_DISALLOW_COPY_AND_MOVE(LevelSetDomainPartitionImage);
37
42
43 static constexpr unsigned int ImageDimension = TImage::ImageDimension;
44
46 itkNewMacro(Self);
47
48 itkOverrideGetNameOfClassMacro(LevelSetDomainPartitionImage);
49
50 using ImageType = TImage;
53 using PixelType = typename ImageType::PixelType;
55 using SizeType = typename ImageType::SizeType;
57 using SpacingType = typename ImageType::SpacingType;
61
62 using typename Superclass::IdentifierListType;
63
75
76 using LevelSetDomainRegionVectorType = std::vector<RegionType>;
77
80 itkSetConstObjectMacro(Image, ImageType);
81 itkGetConstObjectMacro(Image, ImageType);
85 itkGetModifiableObjectMacro(ListDomain, ListImageType);
86
87 void
91
94 void
96
97protected:
99 ~LevelSetDomainPartitionImage() override = default;
100
103 void
105
107 ListImagePointer m_ListDomain{};
108 LevelSetDomainRegionVectorType m_LevelSetDomainRegionVector{};
109};
110} // end namespace itk
111
112#ifndef ITK_MANUAL_INSTANTIATION
113# include "itkLevelSetDomainPartitionImage.hxx"
114#endif
115
116#endif
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Templated n-dimensional image class.
Definition: itkImage.h:89
Helper class used to partition domain and efficiently compute overlap.
Helper class used to partition domain and efficiently compute overlap.
typename ListImageType::IndexType ListIndexType
typename ListIndexType::IndexValueType ListIndexValueType
typename ImageType::ConstPointer ImageConstPointer
const LevelSetDomainRegionVectorType & GetLevelSetDomainRegionVector() const
typename ListImageType::ConstPointer ListImageConstPointer
typename IndexType::IndexValueType IndexValueType
typename ImageType::SpacingType SpacingType
typename ListImageType::PointType ListPointType
void SetLevelSetDomainRegionVector(const LevelSetDomainRegionVectorType &domain)
typename ListImageType::SpacingType ListSpacingType
typename SizeType::SizeValueType SizeValueType
typename ListImageType::SizeType ListSizeType
typename ListImageType::RegionType ListRegionType
typename ListImageType::Pointer ListImagePointer
std::vector< RegionType > LevelSetDomainRegionVectorType
typename ListSizeType::SizeValueType ListSizeValueType
~LevelSetDomainPartitionImage() override=default
Light weight base class for most itk classes.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long IndexValueType
Definition: itkIntTypes.h:93
unsigned long SizeValueType
Definition: itkIntTypes.h:86