ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
52template <typename TInputImage,
53 typename TFeatureImage,
55 TInputImage,
56 TFeatureImage,
58class ITK_TEMPLATE_EXPORT ScalarChanAndVeseLevelSetFunction
59 : public ScalarRegionBasedLevelSetFunction<TInputImage, TFeatureImage, TSharedData>
60{
61public:
62 ITK_DISALLOW_COPY_AND_MOVE(ScalarChanAndVeseLevelSetFunction);
63
68
70 itkNewMacro(Self);
71
73 itkOverrideGetNameOfClassMacro(ScalarChanAndVeseLevelSetFunction);
74
75 static constexpr unsigned int ImageDimension = TFeatureImage::ImageDimension;
76
77 using InputImageType = TInputImage;
79 using typename Superclass::InputImagePointer;
80 using typename Superclass::InputPixelType;
81 using typename Superclass::InputIndexType;
83 using typename Superclass::InputSizeType;
84 using typename Superclass::InputSizeValueType;
85 using typename Superclass::InputRegionType;
86 using typename Superclass::InputPointType;
87
88 using FeatureImageType = TFeatureImage;
89 using FeatureImageConstPointer = typename FeatureImageType::ConstPointer;
90 using typename Superclass::FeaturePixelType;
91 using typename Superclass::FeatureIndexType;
92 using typename Superclass::FeatureOffsetType;
93
94 using typename Superclass::ScalarValueType;
95 using typename Superclass::NeighborhoodType;
96 using typename Superclass::FloatOffsetType;
97 using typename Superclass::RadiusType;
98 using typename Superclass::TimeStepType;
99 using typename Superclass::GlobalDataStruct;
100 using typename Superclass::PixelType;
101 using typename Superclass::VectorType;
102
103 using typename Superclass::SharedDataType;
104 using typename Superclass::SharedDataPointer;
105
106 using typename Superclass::ImageIteratorType;
110
111 using typename Superclass::ListPixelType;
113 using typename Superclass::ListPixelIterator;
114 using typename Superclass::ListImageType;
115
116protected:
121
122 void
124
125 void
127
128 ScalarValueType
129 ComputeInternalTerm(const FeaturePixelType & iValue, const FeatureIndexType & iIdx) override;
130
131 ScalarValueType
132 ComputeExternalTerm(const FeaturePixelType & iValue, const FeatureIndexType & iIdx) override;
133
134 void
135 UpdateSharedDataInsideParameters(const unsigned int & iId,
136 const FeaturePixelType & iVal,
137 const ScalarValueType & iChange) override;
138
139 void
140 UpdateSharedDataOutsideParameters(const unsigned int & iId,
141 const FeaturePixelType & iVal,
142 const ScalarValueType & iChange) override;
143};
144} // namespace itk
145
146#ifndef ITK_MANUAL_INSTANTIATION
147# include "itkScalarChanAndVeseLevelSetFunction.hxx"
148#endif
149
150#endif
Helper class used to share data in the ScalarChanAndVeseLevelSetFunction.
ConstNeighborhoodIterator< TImageType, DefaultBoundaryConditionType > NeighborhoodType
Vector< float, Self::ImageDimension > FloatOffsetType
typename ImageType::PixelType PixelType
typename ConstNeighborhoodIterator< TImageType >::RadiusType RadiusType
Helper class used to share data in the ScalarChanAndVeseLevelSetFunction.
void UpdateSharedDataInsideParameters(const unsigned int &iId, const FeaturePixelType &iVal, const ScalarValueType &iChange) override
typename FeatureImageType::ConstPointer FeatureImageConstPointer
ScalarRegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData > Superclass
~ScalarChanAndVeseLevelSetFunction() override=default
void UpdateSharedDataOutsideParameters(const unsigned int &iId, const FeaturePixelType &iVal, const ScalarValueType &iChange) override
ScalarValueType ComputeInternalTerm(const FeaturePixelType &iValue, const FeatureIndexType &iIdx) override
ScalarValueType ComputeExternalTerm(const FeaturePixelType &iValue, const FeatureIndexType &iIdx) override
LevelSet function that computes a speed image based on regional integrals.
ImageRegionConstIteratorWithIndex< InputImageType > ConstImageIteratorType
ImageRegionConstIterator< FeatureImageType > ConstFeatureIteratorType
ImageRegionIteratorWithIndex< FeatureImageType > FeatureImageIteratorType
ImageRegionIteratorWithIndex< InputImageType > ImageIteratorType
typename ListPixelType::const_iterator ListPixelConstIterator
Image< ListPixelType, Self::ImageDimension > ListImageType
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....