ITK  6.0.0
Insight Toolkit
itkRegionBasedLevelSetFunctionData.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 itkRegionBasedLevelSetFunctionData_h
19#define itkRegionBasedLevelSetFunctionData_h
20
21#include "itkLightObject.h"
22
23#include "itkVector.h"
24#include "itkListSample.h"
25#include "itkKdTreeGenerator.h"
26
29
30namespace itk
31{
64template <typename TInputImage, typename TFeatureImage>
65class ITK_TEMPLATE_EXPORT RegionBasedLevelSetFunctionData : public LightObject
66{
67public:
68 ITK_DISALLOW_COPY_AND_MOVE(RegionBasedLevelSetFunctionData);
69
74
75 static constexpr unsigned int ImageDimension = TFeatureImage::ImageDimension;
76
78 itkNewMacro(Self);
79
80 itkOverrideGetNameOfClassMacro(RegionBasedLevelSetFunctionData);
81
82 using InputImageType = TInputImage;
85 using InputPixelType = typename InputImageType::PixelType;
89 using InputSpacingType = typename InputImageType::SpacingType;
93
94 using FeatureImageType = TFeatureImage;
97 using FeaturePixelType = typename FeatureImageType::PixelType;
101 using FeatureSpacingType = typename FeatureImageType::SpacingType;
104
105 // Allocates m_HeavisideFunctionOfLevelSetImage to have same origin,
106 // spacing and size as image. Also sets the m_Start and m_End indices.
107 void
109
110 // Checks if the given index lies in the domain of the current
111 // level-set function. The domain is defined by the start and end indices.
112 template <typename TIndex>
113 bool
114 VerifyInsideRegion(const TIndex & featureIndex)
115 {
116 for (unsigned int j = 0; j < ImageDimension; ++j)
117 {
118 if ((featureIndex[j] < static_cast<InputIndexValueType>(this->m_Start[j])) ||
119 (featureIndex[j] > static_cast<InputIndexValueType>(this->m_End[j])))
120 {
121 return false;
122 }
123 }
124 return true;
125 }
126
127 // Get the index into the domain of the current level-set function
128 InputIndexType
129 GetIndex(const FeatureIndexType & featureIndex);
130
131 // Get the index in the domain of the feature image
133 GetFeatureIndex(const InputIndexType & inputIndex);
134
135 double m_WeightedNumberOfPixelsInsideLevelSet{};
136 double m_WeightedNumberOfPixelsOutsideLevelSet{};
137
138 InputImagePointer m_HeavisideFunctionOfLevelSetImage{};
139 InputIndexType m_Start{};
141
142protected:
145};
146} // end namespace itk
147
148#ifndef ITK_MANUAL_INSTANTIATION
149# include "itkRegionBasedLevelSetFunctionData.hxx"
150#endif
151#endif
Light weight base class for most itk classes.
Helper class used to share data in the ScalarChanAndVeseLevelSetFunction.
typename InputImageType::PixelType InputPixelType
typename FeatureSizeType::SizeValueType FeatureSizeValueType
typename InputImageType::PointType InputPointType
typename InputImageType::RegionType InputRegionType
typename InputImageType::ConstPointer InputImageConstPointer
InputIndexType GetIndex(const FeatureIndexType &featureIndex)
void CreateHeavisideFunctionOfLevelSetImage(const InputImageType *image)
typename FeatureImageType::SpacingType FeatureSpacingType
typename InputImageType::SpacingType InputSpacingType
typename InputImageType::IndexType InputIndexType
typename FeatureImageType::RegionType FeatureRegionType
typename FeatureImageType::IndexType FeatureIndexType
FeatureIndexType GetFeatureIndex(const InputIndexType &inputIndex)
typename InputIndexType::IndexValueType InputIndexValueType
typename FeatureImageType::ConstPointer FeatureImageConstPointer
typename FeatureImageType::PointType FeaturePointType
typename FeatureImageType::PixelType FeaturePixelType
typename InputSizeType::SizeValueType InputSizeValueType
typename FeatureImageType::SizeType FeatureSizeType
typename FeatureImageType::Pointer FeatureImagePointer
~RegionBasedLevelSetFunctionData() override=default
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