ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkLevelSetEvolution.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 itkLevelSetEvolution_h
19#define itkLevelSetEvolution_h
20
23
26
29
32
35
36namespace itk
37{
47template <typename TEquationContainer, typename TLevelSet>
48class ITK_TEMPLATE_EXPORT LevelSetEvolution{};
49
50template <typename TEquationContainer, typename TImage>
51class ITK_TEMPLATE_EXPORT LevelSetEvolution<TEquationContainer, LevelSetDenseImage<TImage>>
52 : public LevelSetEvolutionBase<TEquationContainer, LevelSetDenseImage<TImage>>
53{
54public:
55 ITK_DISALLOW_COPY_AND_MOVE(LevelSetEvolution);
56
58
63
65 itkNewMacro(Self);
66
68 itkOverrideGetNameOfClassMacro(LevelSetEvolution);
69
72 using typename Superclass::TermContainerType;
74
75 using typename Superclass::TermType;
76 using typename Superclass::TermPointer;
77
78 using typename Superclass::InputImageType;
82 using typename Superclass::InputPixelRealType;
83
84 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
85
87
89
91
92 using typename Superclass::LevelSetOutputType;
94 using typename Superclass::LevelSetDataType;
95
96 using typename Superclass::IdListType;
97 using typename Superclass::IdListIterator;
99 using typename Superclass::IdListImageType;
100 using typename Superclass::CacheImageType;
102
105
108
111
113
115
117
119 void
120 SetNumberOfWorkUnits(const ThreadIdType numberOfThreads);
121
125
126 ~LevelSetEvolution() override = default;
127
128protected:
130
133 void
135
137 void
139
141 void
143
145 void
146 UpdateLevelSets() override;
147
149 void
150 UpdateEquations() override;
151
153 void
155
156 typename LevelSetContainerType::Pointer m_UpdateBuffer{};
157
158 friend class LevelSetEvolutionComputeIterationThreader<LevelSetType,
159 ThreadedImageRegionPartitioner<TImage::ImageDimension>,
160 Self>;
163 ThreadedImageRegionPartitioner<TImage::ImageDimension>,
164 Self>;
165 typename SplitLevelSetComputeIterationThreaderType::Pointer m_SplitLevelSetComputeIterationThreader{};
166
167 using DomainMapConstIteratorType = typename DomainMapImageFilterType::DomainMapType::const_iterator;
172 typename SplitDomainMapComputeIterationThreaderType::Pointer m_SplitDomainMapComputeIterationThreader{};
173
174 friend class LevelSetEvolutionUpdateLevelSetsThreader<LevelSetType,
175 ThreadedImageRegionPartitioner<TImage::ImageDimension>,
176 Self>;
179 ThreadedImageRegionPartitioner<TImage::ImageDimension>,
180 Self>;
181 typename SplitLevelSetUpdateLevelSetsThreaderType::Pointer m_SplitLevelSetUpdateLevelSetsThreader{};
182
185};
186
187
188template <typename TEquationContainer, typename TOutput, unsigned int VDimension>
189class ITK_TEMPLATE_EXPORT LevelSetEvolution<TEquationContainer, WhitakerSparseLevelSetImage<TOutput, VDimension>>
190 : public LevelSetEvolutionBase<TEquationContainer, WhitakerSparseLevelSetImage<TOutput, VDimension>>
191{
192public:
193 ITK_DISALLOW_COPY_AND_MOVE(LevelSetEvolution);
194
196
201
203 itkNewMacro(Self);
204
206 itkOverrideGetNameOfClassMacro(LevelSetEvolution);
207
210 using typename Superclass::TermContainerType;
212
213 using typename Superclass::TermType;
214 using typename Superclass::TermPointer;
215
216 using typename Superclass::InputImageType;
217 using typename Superclass::InputImagePixelType;
220 using typename Superclass::InputPixelRealType;
221
222 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
223
226
227 using typename Superclass::LevelSetInputType;
228 using typename Superclass::LevelSetOutputType;
230 using typename Superclass::LevelSetDataType;
231
233
236
237
238 using typename Superclass::IdListType;
239 using typename Superclass::IdListIterator;
240 using typename Superclass::IdListImageType;
241 using typename Superclass::CacheImageType;
243
246
248
252
254 void
255 SetNumberOfWorkUnits(const ThreadIdType numberOfWorkUnits);
256
260
261protected:
264
265 using NodePairType = std::pair<LevelSetInputType, LevelSetOutputType>;
266
267 // For sparse case, the update buffer needs to be the size of the active layer
268 std::map<IdentifierType, LevelSetLayerType *> m_UpdateBuffer{};
269
272 void
274
276 void
278
280 void
282
284 void
285 UpdateLevelSets() override;
286
288 void
289 UpdateEquations() override;
290
295 typename SplitLevelSetComputeIterationThreaderType::Pointer m_SplitLevelSetComputeIterationThreader{};
296};
297
298
299// Shi
300template <typename TEquationContainer, unsigned int VDimension>
301class ITK_TEMPLATE_EXPORT LevelSetEvolution<TEquationContainer, ShiSparseLevelSetImage<VDimension>>
302 : public LevelSetEvolutionBase<TEquationContainer, ShiSparseLevelSetImage<VDimension>>
303{
304public:
305 ITK_DISALLOW_COPY_AND_MOVE(LevelSetEvolution);
306
308
313
315 itkNewMacro(Self);
316
318 itkOverrideGetNameOfClassMacro(LevelSetEvolution);
319
322 using typename Superclass::TermContainerType;
324
325 using typename Superclass::TermType;
326 using typename Superclass::TermPointer;
327
328 using typename Superclass::InputImageType;
329 using typename Superclass::InputImagePixelType;
332 using typename Superclass::InputPixelRealType;
333
334 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
335
338
339 using typename Superclass::LevelSetInputType;
340 using typename Superclass::LevelSetOutputType;
342 using typename Superclass::LevelSetDataType;
343
345
348
349
350 using typename Superclass::IdListType;
351 using typename Superclass::IdListIterator;
352 using typename Superclass::IdListImageType;
353 using typename Superclass::CacheImageType;
355
358
360
363
364 LevelSetEvolution() = default;
365 ~LevelSetEvolution() override = default;
366
367protected:
369 void
370 UpdateLevelSets() override;
371
373 void
374 UpdateEquations() override;
375};
376
377// Malcolm
378template <typename TEquationContainer, unsigned int VDimension>
379class ITK_TEMPLATE_EXPORT LevelSetEvolution<TEquationContainer, MalcolmSparseLevelSetImage<VDimension>>
380 : public LevelSetEvolutionBase<TEquationContainer, MalcolmSparseLevelSetImage<VDimension>>
381{
382public:
383 ITK_DISALLOW_COPY_AND_MOVE(LevelSetEvolution);
384
386
391
393 itkNewMacro(Self);
394
396 itkOverrideGetNameOfClassMacro(LevelSetEvolution);
397
400 using typename Superclass::TermContainerType;
402
403 using typename Superclass::TermType;
404 using typename Superclass::TermPointer;
405
406 using typename Superclass::InputImageType;
407 using typename Superclass::InputImagePixelType;
410 using typename Superclass::InputPixelRealType;
411
412 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
413
416
417 using typename Superclass::LevelSetInputType;
418 using typename Superclass::LevelSetOutputType;
420 using typename Superclass::LevelSetDataType;
421
424
427
428
429 using typename Superclass::IdListType;
430 using typename Superclass::IdListIterator;
431 using typename Superclass::IdListImageType;
432 using typename Superclass::CacheImageType;
434
437
439
442
443 LevelSetEvolution() = default;
444 ~LevelSetEvolution() override = default;
445
446protected:
447 void
448 UpdateLevelSets() override;
449 void
450 UpdateEquations() override;
451};
452} // namespace itk
453
454#ifndef ITK_MANUAL_INSTANTIATION
455# include "itkLevelSetEvolution.hxx"
456#endif
457
458#endif // itkLevelSetEvolution_h
Binarize an input image by thresholding.
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 pixels within a region and is speci...
Base class for the "dense" representation of a level-set function on one image.
typename TermContainerType::LevelSetContainerType LevelSetContainerType
typename LevelSetContainerType::IdListImageType IdListImageType
typename NumericTraits< InputImagePixelType >::RealType InputPixelRealType
typename LevelSetContainerType::CacheImageType CacheImageType
typename LevelSetContainerType::IdListConstIterator IdListConstIterator
typename LevelSetContainerType::DomainMapImageFilterType DomainMapImageFilterType
typename LevelSetContainerType::IdListIterator IdListIterator
typename LevelSetContainerType::LevelSetIdentifierType LevelSetIdentifierType
typename EquationContainerType::TermContainerType TermContainerType
LevelSetEvolutionStoppingCriterion< LevelSetContainerType > StoppingCriterionType
LevelSetEvolutionUpdateLevelSetsThreader< LevelSetType, ThreadedImageRegionPartitioner< TImage::ImageDimension >, Self > SplitLevelSetUpdateLevelSetsThreaderType
LevelSetEvolutionComputeIterationThreader< LevelSetType, ThreadedDomainMapPartitionerType, Self > SplitDomainMapComputeIterationThreaderType
LevelSetEvolutionComputeIterationThreader< LevelSetType, ThreadedImageRegionPartitioner< TImage::ImageDimension >, Self > SplitLevelSetComputeIterationThreaderType
SplitLevelSetComputeIterationThreaderType::Pointer m_SplitLevelSetComputeIterationThreader
ImageRegionConstIteratorWithIndex< LevelSetImageType > LevelSetImageConstIteratorType
SignedMaurerDistanceMapImageFilter< LevelSetImageType, LevelSetImageType > MaurerType
ImageRegionConstIteratorWithIndex< InputImageType > InputImageConstIteratorType
SplitDomainMapComputeIterationThreaderType::Pointer m_SplitDomainMapComputeIterationThreader
ThreadedIteratorRangePartitioner< DomainMapConstIteratorType > ThreadedDomainMapPartitionerType
ImageRegionIteratorWithIndex< LevelSetImageType > LevelSetImageIteratorType
typename DomainMapImageFilterType::DomainMapType::const_iterator DomainMapConstIteratorType
SplitLevelSetUpdateLevelSetsThreaderType::Pointer m_SplitLevelSetUpdateLevelSetsThreader
BinaryThresholdImageFilter< LevelSetImageType, LevelSetImageType > ThresholdFilterType
LevelSetEvolutionBase< TEquationContainer, LevelSetType > Superclass
UpdateMalcolmSparseLevelSet< ImageDimension, EquationContainerType > UpdateLevelSetFilterType
UpdateShiSparseLevelSet< ImageDimension, EquationContainerType > UpdateLevelSetFilterType
ThreadedIteratorRangePartitioner< typename LevelSetType::LayerConstIterator > SplitLevelSetPartitionerType
UpdateWhitakerSparseLevelSet< ImageDimension, LevelSetOutputType, EquationContainerType > UpdateLevelSetFilterType
LevelSetEvolutionComputeIterationThreader< LevelSetType, SplitLevelSetPartitionerType, Self > SplitLevelSetComputeIterationThreaderType
Class for iterating and evolving the level-set function.
std::map< InputType, OutputType, Functor::LexicographicCompare > LayerType
Derived class for the Malcolm representation of level-set function.
Derived class for the shi representation of level-set function.
This filter calculates the Euclidean distance transform of a binary image in linear time for arbitrar...
Implements transparent reference counting.
Class for partitioning of an ImageRegion.
Partitions an iterator range for threading.
Base class for updating the Malcolm representation of level-set function.
Base class for updating the Shi representation of level-set function.
Base class for updating the level-set function.
Derived class for the sparse-field representation of level-set function.
std::map< InputType, OutputType, Functor::LexicographicCompare > LayerType
LabelMap< LabelObjectType > LabelMapType
typename LabelMapType::Pointer LabelMapPointer
AddImageFilter Self
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType