ITK  6.0.0
Insight Toolkit
itkScalarRegionBasedLevelSetFunction.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 itkScalarRegionBasedLevelSetFunction_h
19#define itkScalarRegionBasedLevelSetFunction_h
20
25
26namespace itk
27{
62template <typename TInputImage, typename TFeatureImage, typename TSharedData>
63class ITK_TEMPLATE_EXPORT ScalarRegionBasedLevelSetFunction
64 : public RegionBasedLevelSetFunction<TInputImage, TFeatureImage, TSharedData>
65{
66public:
67 ITK_DISALLOW_COPY_AND_MOVE(ScalarRegionBasedLevelSetFunction);
68
73
74 // itkNewMacro() is not provided since this is an abstract class.
75
77 itkOverrideGetNameOfClassMacro(ScalarRegionBasedLevelSetFunction);
78
79 static constexpr unsigned int ImageDimension = TFeatureImage::ImageDimension;
80
81 using typename Superclass::InputImageType;
83 using typename Superclass::InputImagePointer;
84 using typename Superclass::InputPixelType;
85 using typename Superclass::InputIndexType;
87 using typename Superclass::InputSizeType;
88 using typename Superclass::InputSizeValueType;
89 using typename Superclass::InputRegionType;
90 using typename Superclass::InputPointType;
91
92 using typename Superclass::FeatureImageType;
94 using typename Superclass::FeaturePixelType;
95 using typename Superclass::FeatureIndexType;
96 using typename Superclass::FeatureOffsetType;
97
98 using typename Superclass::ScalarValueType;
99 using typename Superclass::NeighborhoodType;
100 using typename Superclass::FloatOffsetType;
101 using typename Superclass::RadiusType;
102 using typename Superclass::TimeStepType;
103 using typename Superclass::GlobalDataStruct;
104 using typename Superclass::PixelType;
105 using typename Superclass::VectorType;
106
107 using typename Superclass::SharedDataType;
108 using typename Superclass::SharedDataPointer;
109
114
115 using ListPixelType = std::list<unsigned int>;
116 using ListPixelConstIterator = typename ListPixelType::const_iterator;
117 using ListPixelIterator = typename ListPixelType::iterator;
119
124 void
125 UpdatePixel(const unsigned int idx,
127 InputPixelType & newValue,
128 bool & status);
129
130protected:
132 : Superclass()
133 {}
135
139 ComputeOverlapParameters(const FeatureIndexType & featIndex, ScalarValueType & product) override;
140
144 virtual void
145 UpdateSharedDataInsideParameters(const unsigned int & iId,
146 const FeaturePixelType & iVal,
147 const ScalarValueType & iChange) = 0;
148
149 virtual void
150 UpdateSharedDataOutsideParameters(const unsigned int & iId,
151 const FeaturePixelType & iVal,
152 const ScalarValueType & iChange) = 0;
153};
154} // namespace itk
155
156#ifndef ITK_MANUAL_INSTANTIATION
157# include "itkScalarRegionBasedLevelSetFunction.hxx"
158#endif
159
160#endif
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
typename ConstNeighborhoodIterator< TInputImage >::RadiusType RadiusType
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
A multi-dimensional iterator templated over image type that walks a region of pixels.
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
Light weight base class for most itk classes.
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
LevelSet function that computes a speed image based on regional integrals.
LevelSet function that computes a speed image based on regional integrals.
void UpdatePixel(const unsigned int idx, NeighborhoodIterator< TInputImage > &iterator, InputPixelType &newValue, bool &status)
virtual void UpdateSharedDataInsideParameters(const unsigned int &iId, const FeaturePixelType &iVal, const ScalarValueType &iChange)=0
virtual void UpdateSharedDataOutsideParameters(const unsigned int &iId, const FeaturePixelType &iVal, const ScalarValueType &iChange)=0
typename FeatureImageType::ConstPointer FeatureImageConstPointer
ScalarValueType ComputeOverlapParameters(const FeatureIndexType &featIndex, ScalarValueType &product) override
~ScalarRegionBasedLevelSetFunction() override=default
typename ListPixelType::const_iterator ListPixelConstIterator
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:63
SmartPointer< const Self > ConstPointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....