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 
25 namespace itk
26 {
72 template <typename TInputImage,
73  typename TFeatureImage,
74  typename TSharedData = ConstrainedRegionBasedLevelSetFunctionSharedData<
75  TInputImage,
76  TFeatureImage,
77  ScalarChanAndVeseLevelSetFunctionData<TInputImage, TFeatureImage>>>
78 class ITK_TEMPLATE_EXPORT ScalarChanAndVeseLevelSetFunction
79  : public ScalarRegionBasedLevelSetFunction<TInputImage, TFeatureImage, TSharedData>
80 {
81 public:
82  ITK_DISALLOW_COPY_AND_MOVE(ScalarChanAndVeseLevelSetFunction);
83 
88 
90  itkNewMacro(Self);
91 
93  itkOverrideGetNameOfClassMacro(ScalarChanAndVeseLevelSetFunction);
94 
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;
105  using typename Superclass::InputRegionType;
106  using typename Superclass::InputPointType;
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;
127  using typename Superclass::ConstImageIteratorType;
130 
131  using typename Superclass::ListPixelType;
132  using typename Superclass::ListPixelConstIterator;
133  using typename Superclass::ListPixelIterator;
134  using typename Superclass::ListImageType;
135 
136 protected:
138  : Superclass()
139  {}
140  ~ScalarChanAndVeseLevelSetFunction() override = default;
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
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputIndexType
typename InputImageType::IndexType InputIndexType
Definition: itkRegionBasedLevelSetFunction.h:127
itk::FiniteDifferenceFunction< TInputImage >::RadiusType
typename ConstNeighborhoodIterator< TInputImage >::RadiusType RadiusType
Definition: itkFiniteDifferenceFunction.h:97
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::ImageRegionConstIteratorWithIndex
A multi-dimensional iterator templated over image type that walks an image region and is specialized ...
Definition: itkImageRegionConstIteratorWithIndex.h:130
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputSizeValueType
typename InputImageType::SizeValueType InputSizeValueType
Definition: itkRegionBasedLevelSetFunction.h:130
itkConstrainedRegionBasedLevelSetFunctionSharedData.h
itk::GTest::TypedefsAndConstructors::Dimension2::VectorType
ImageBaseType::SpacingType VectorType
Definition: itkGTestTypedefsAndConstructors.h:53
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputPointType
typename InputImageType::PointType InputPointType
Definition: itkRegionBasedLevelSetFunction.h:132
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::ScalarRegionBasedLevelSetFunction::ListPixelIterator
typename ListPixelType::iterator ListPixelIterator
Definition: itkScalarRegionBasedLevelSetFunction.h:117
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputSizeType
typename InputImageType::SizeType InputSizeType
Definition: itkRegionBasedLevelSetFunction.h:129
itk::SmartPointer< Self >
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::SharedDataType
TSharedData SharedDataType
Definition: itkRegionBasedLevelSetFunction.h:141
itk::ScalarChanAndVeseLevelSetFunction::FeatureImageType
TFeatureImage FeatureImageType
Definition: itkScalarChanAndVeseLevelSetFunction.h:108
itk::ScalarChanAndVeseLevelSetFunction::FeatureImageConstPointer
typename FeatureImageType::ConstPointer FeatureImageConstPointer
Definition: itkScalarChanAndVeseLevelSetFunction.h:109
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::FeatureOffsetType
typename FeatureImageType::OffsetType FeatureOffsetType
Definition: itkRegionBasedLevelSetFunction.h:139
itk::ScalarRegionBasedLevelSetFunction
LevelSet function that computes a speed image based on regional integrals.
Definition: itkScalarRegionBasedLevelSetFunction.h:63
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::FeaturePixelType
typename FeatureImageType::PixelType FeaturePixelType
Definition: itkRegionBasedLevelSetFunction.h:136
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputPixelType
typename InputImageType::PixelType InputPixelType
Definition: itkRegionBasedLevelSetFunction.h:126
itk::FiniteDifferenceFunction< TInputImage >::PixelType
typename ImageType::PixelType PixelType
Definition: itkFiniteDifferenceFunction.h:83
itk::ScalarChanAndVeseLevelSetFunction::ScalarChanAndVeseLevelSetFunction
ScalarChanAndVeseLevelSetFunction()
Definition: itkScalarChanAndVeseLevelSetFunction.h:137
itk::ScalarRegionBasedLevelSetFunction::ListPixelType
std::list< unsigned int > ListPixelType
Definition: itkScalarRegionBasedLevelSetFunction.h:115
itk::ImageRegionIteratorWithIndex
A multi-dimensional iterator templated over image type that walks pixels within a region and is speci...
Definition: itkImageRegionIteratorWithIndex.h:73
itk::ScalarChanAndVeseLevelSetFunction
LevelSet function that computes a speed image based on regional integrals of probabilities.
Definition: itkScalarChanAndVeseLevelSetFunction.h:78
itkScalarChanAndVeseLevelSetFunctionData.h
itk::ScalarChanAndVeseLevelSetFunction::InputImageType
TInputImage InputImageType
Definition: itkScalarChanAndVeseLevelSetFunction.h:97
itk::ScalarRegionBasedLevelSetFunction::ListPixelConstIterator
typename ListPixelType::const_iterator ListPixelConstIterator
Definition: itkScalarRegionBasedLevelSetFunction.h:116
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkRegionBasedLevelSetFunction.h:124
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::FeatureIndexType
typename FeatureImageType::IndexType FeatureIndexType
Definition: itkRegionBasedLevelSetFunction.h:137
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputIndexValueType
typename InputImageType::IndexValueType InputIndexValueType
Definition: itkRegionBasedLevelSetFunction.h:128
itk::ImageRegionConstIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionConstIterator.h:109
itkScalarRegionBasedLevelSetFunction.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::ScalarValueType
PixelType ScalarValueType
Definition: itkRegionBasedLevelSetFunction.h:86
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::SharedDataPointer
typename SharedDataType::Pointer SharedDataPointer
Definition: itkRegionBasedLevelSetFunction.h:142
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkRegionBasedLevelSetFunction.h:125
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputRegionType
typename InputImageType::RegionType InputRegionType
Definition: itkRegionBasedLevelSetFunction.h:131
itk::FiniteDifferenceFunction< TInputImage >::TimeStepType
double TimeStepType
Definition: itkFiniteDifferenceFunction.h:90