ITK  5.4.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 
26 namespace itk
27 {
62 template <typename TInputImage, typename TFeatureImage, typename TSharedData>
63 class ITK_TEMPLATE_EXPORT ScalarRegionBasedLevelSetFunction
64  : public RegionBasedLevelSetFunction<TInputImage, TFeatureImage, TSharedData>
65 {
66 public:
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;
86  using typename Superclass::InputIndexValueType;
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 
130 protected:
132  : Superclass()
133  {}
134  ~ScalarRegionBasedLevelSetFunction() override = default;
135 
138  ScalarValueType
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
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
itkNeighborhoodIterator.h
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputSizeValueType
typename InputImageType::SizeValueType InputSizeValueType
Definition: itkRegionBasedLevelSetFunction.h:130
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 >
itkImageRegionIteratorWithIndex.h
itkImageRegionConstIterator.h
itk::RegionBasedLevelSetFunction
LevelSet function that computes a speed image based on regional integrals.
Definition: itkRegionBasedLevelSetFunction.h:64
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::SharedDataType
TSharedData SharedDataType
Definition: itkRegionBasedLevelSetFunction.h:141
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::ScalarRegionBasedLevelSetFunction::ScalarRegionBasedLevelSetFunction
ScalarRegionBasedLevelSetFunction()
Definition: itkScalarRegionBasedLevelSetFunction.h:131
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
itkRegionBasedLevelSetFunction.h
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::NeighborhoodIterator
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Definition: itkNeighborhoodIterator.h:212
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: itkAnnulusOperator.h:24
itk::ScalarRegionBasedLevelSetFunction::FeatureImageConstPointer
typename FeatureImageType::ConstPointer FeatureImageConstPointer
Definition: itkScalarRegionBasedLevelSetFunction.h:93
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::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::FeatureImageType
TFeatureImage FeatureImageType
Definition: itkRegionBasedLevelSetFunction.h:134
itk::ImageRegionConstIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionConstIterator.h:109
itk::RegionBasedLevelSetFunction< TInputImage, TFeatureImage, TSharedData >::InputImageType
TInputImage InputImageType
Definition: itkRegionBasedLevelSetFunction.h:123
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 >::TimeStepType
double TimeStepType
Definition: itkRegionBasedLevelSetFunction.h:83
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