ITK  6.0.0
Insight Toolkit
itkScalarChanAndVeseLevelSetFunction.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 itkScalarChanAndVeseLevelSetFunction_h
19#define itkScalarChanAndVeseLevelSetFunction_h
20
24
25namespace itk
26{
72template <typename TInputImage,
73 typename TFeatureImage,
74 typename TSharedData = ConstrainedRegionBasedLevelSetFunctionSharedData<
75 TInputImage,
76 TFeatureImage,
78class ITK_TEMPLATE_EXPORT ScalarChanAndVeseLevelSetFunction
79 : public ScalarRegionBasedLevelSetFunction<TInputImage, TFeatureImage, TSharedData>
80{
81public:
82 ITK_DISALLOW_COPY_AND_MOVE(ScalarChanAndVeseLevelSetFunction);
88
90 itkNewMacro(Self);
91
93 itkOverrideGetNameOfClassMacro(ScalarChanAndVeseLevelSetFunction);
95 static constexpr unsigned int ImageDimension = TFeatureImage::ImageDimension;
96
97 using InputImageType = TInputImage;
99 using typename Superclass::InputImagePointer;
100 using typename Superclass::InputPixelType;
101 using typename Superclass::InputIndexType;
102 using typename Superclass::InputIndexValueType;
103 using typename Superclass::InputSizeType;
104 using typename Superclass::InputSizeValueType;
107
108 using FeatureImageType = TFeatureImage;
110 using typename Superclass::FeaturePixelType;
111 using typename Superclass::FeatureIndexType;
112 using typename Superclass::FeatureOffsetType;
113
114 using typename Superclass::ScalarValueType;
115 using typename Superclass::NeighborhoodType;
116 using typename Superclass::FloatOffsetType;
117 using typename Superclass::RadiusType;
118 using typename Superclass::TimeStepType;
119 using typename Superclass::GlobalDataStruct;
120 using typename Superclass::PixelType;
121 using typename Superclass::VectorType;
122
123 using typename Superclass::SharedDataType;
124 using typename Superclass::SharedDataPointer;
125
126 using typename Superclass::ImageIteratorType;
130
131 using typename Superclass::ListPixelType;
133 using typename Superclass::ListPixelIterator;
135
136protected:
138 : Superclass()
139 {}
141
142 void
143 ComputeParameters() override;
144
145 void
146 UpdateSharedDataParameters() override;
147
148 ScalarValueType
149 ComputeInternalTerm(const FeaturePixelType & iValue, const FeatureIndexType & iIdx) override;
150
151 ScalarValueType
152 ComputeExternalTerm(const FeaturePixelType & iValue, const FeatureIndexType & iIdx) override;
153
154 void
155 UpdateSharedDataInsideParameters(const unsigned int & iId,
156 const FeaturePixelType & iVal,
157 const ScalarValueType & iChange) override;
158
159 void
160 UpdateSharedDataOutsideParameters(const unsigned int & iId,
161 const FeaturePixelType & iVal,
162 const ScalarValueType & iChange) override;
163};
164} // namespace itk
165
166#ifndef ITK_MANUAL_INSTANTIATION
167# include "itkScalarChanAndVeseLevelSetFunction.hxx"
168#endif
169
170#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
Helper class used to share data in the ScalarChanAndVeseLevelSetFunction.
LevelSet function that computes a speed image based on regional integrals of probabilities.
typename FeatureImageType::ConstPointer FeatureImageConstPointer
LevelSet function that computes a speed image based on regional integrals.
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....