ITK  5.4.0
Insight Toolkit
itkMultiphaseFiniteDifferenceImageFilter.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 itkMultiphaseFiniteDifferenceImageFilter_h
19#define itkMultiphaseFiniteDifferenceImageFilter_h
20
23#include "vnl/vnl_vector.h"
25
26#include "itkListSample.h"
27#include "itkKdTreeGenerator.h"
28
29namespace itk
30{
157template <typename TInputImage,
158 typename TFeatureImage,
159 typename TOutputImage,
160 typename TFiniteDifferenceFunction = FiniteDifferenceFunction<TOutputImage>,
161 typename TIdCell = unsigned int>
162class ITK_TEMPLATE_EXPORT MultiphaseFiniteDifferenceImageFilter : public InPlaceImageFilter<TFeatureImage, TOutputImage>
163{
164public:
165 ITK_DISALLOW_COPY_AND_MOVE(MultiphaseFiniteDifferenceImageFilter);
166
172
174 itkOverrideGetNameOfClassMacro(MultiphaseFiniteDifferenceImageFilter);
175
177 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
178
180 using InputImageType = TInputImage;
183 using InputCoordRepType = typename InputPointType::CoordRepType;
189 using InputPixelType = typename InputImageType::PixelType;
190 using InputSpacingType = typename InputImageType::SpacingType;
192
193 using FeatureImageType = TFeatureImage;
197 using FeatureSpacingType = typename FeatureImageType::SpacingType;
199 using FeaturePixelType = typename FeatureImageType::PixelType;
200
201 using OutputImageType = TOutputImage;
203 using OutputPixelType = typename OutputImageType::PixelType;
209
210 using IdCellType = TIdCell;
211 using VectorIdCellType = std::vector<IdCellType>;
212
216 using FiniteDifferenceFunctionType = TFiniteDifferenceFunction;
218 using TimeStepType = typename FiniteDifferenceFunctionType::TimeStepType;
219 using TimeStepVectorType = typename std::vector<TimeStepType>;
220 using RadiusType = typename FiniteDifferenceFunctionType::RadiusType;
221
228
234 GetDifferenceFunction(const IdCellType & functionIndex) const
235 {
236 if (functionIndex < m_FunctionCount)
237 {
238 return (this->m_DifferenceFunctions[functionIndex]);
239 }
240 else
241 {
242 return nullptr;
243 }
244 }
251 virtual void
253 {
254 if (functionIndex < m_FunctionCount)
255 {
256 this->m_DifferenceFunctions[functionIndex] = function;
257 }
258 }
262 itkSetMacro(NumberOfIterations, unsigned int);
263 itkGetConstReferenceMacro(NumberOfIterations, unsigned int);
268 itkSetMacro(UseImageSpacing, bool);
269 itkBooleanMacro(UseImageSpacing);
270 itkGetConstReferenceMacro(UseImageSpacing, bool);
275 itkSetMacro(MaximumRMSError, double);
276 itkGetConstReferenceMacro(MaximumRMSError, double);
281 itkSetMacro(RMSChange, double);
282 itkGetConstReferenceMacro(RMSChange, double);
286 itkSetMacro(InitializedState, bool);
287 itkGetConstReferenceMacro(InitializedState, bool);
288 itkBooleanMacro(InitializedState);
293 itkSetMacro(ManualReinitialization, bool);
294 itkGetConstReferenceMacro(ManualReinitialization, bool);
295 itkBooleanMacro(ManualReinitialization);
299 itkSetMacro(ElapsedIterations, unsigned int);
300
302 itkGetConstReferenceMacro(ElapsedIterations, unsigned int);
303
304 void
305 SetLevelSet(const IdCellType & i, const InputImageType * levelSet)
306 {
307 m_LevelSet[i] = InputImageType::New();
308 m_LevelSet[i]->SetRequestedRegion(levelSet->GetRequestedRegion());
309 m_LevelSet[i]->SetBufferedRegion(levelSet->GetBufferedRegion());
310 m_LevelSet[i]->SetLargestPossibleRegion(levelSet->GetLargestPossibleRegion());
311 m_LevelSet[i]->Allocate();
312 m_LevelSet[i]->CopyInformation(levelSet);
313
314 ImageRegionConstIterator<InputImageType> in(levelSet, levelSet->GetBufferedRegion());
315 ImageRegionIterator<InputImageType> cp(m_LevelSet[i], levelSet->GetBufferedRegion());
316
317 in.GoToBegin();
318 cp.GoToBegin();
319
320 while (!in.IsAtEnd())
321 {
322 cp.Set(in.Get());
323 ++in;
324 ++cp;
325 }
326 }
327
328 InputImagePointer
330 {
331 if (i >= m_FunctionCount)
332 {
333 itkExceptionMacro("Request for level set #" << i << " but there are only " << m_FunctionCount);
334 }
335 else
336 {
337 return m_LevelSet[i];
338 }
339 }
340
341 void
343 {
344 this->m_Lookup = lookup;
345 }
346
347 void
349 {
350 this->m_KdTree = kdtree;
351 }
352
353 void
355 {
356 m_FunctionCount = n;
357
358 m_DifferenceFunctions.resize(m_FunctionCount, nullptr);
359
360 RadiusType radius;
361 radius.Fill(1);
362
363 for (unsigned int i = 0; i < this->m_FunctionCount; ++i)
364 {
365 this->m_DifferenceFunctions[i] = FiniteDifferenceFunctionType::New();
366 this->m_DifferenceFunctions[i]->Initialize(radius);
367 }
368
369 // Initialize the images
370 m_LevelSet.resize(m_FunctionCount, nullptr);
371
372 // Initialize the lookup table
373 this->m_Lookup.resize(m_FunctionCount);
374
375 IdCellType k = 1;
376
377 {
378 auto it = this->m_Lookup.begin();
379 auto _end = this->m_Lookup.end();
380 while (it != _end)
381 {
382 *it = k;
383 ++it;
384 ++k;
385 }
386 }
387 }
388
389protected:
391 {
392 this->m_KdTree = nullptr;
393 this->m_ElapsedIterations = 0;
394 this->m_MaximumRMSError = itk::Math::eps;
395 this->m_RMSChange = NumericTraits<double>::max();
396 this->m_UseImageSpacing = true;
397 this->m_ManualReinitialization = false;
398 this->m_InitializedState = false;
399 this->m_NumberOfIterations = NumericTraits<unsigned int>::max();
400 this->m_FunctionCount = 0;
401 this->InPlaceOff();
402 }
403
405
406 IdCellType m_FunctionCount{};
407 std::vector<InputImagePointer> m_LevelSet{};
409 KdTreePointer m_KdTree{};
410
411 unsigned int m_ElapsedIterations{};
412 double m_MaximumRMSError{};
413 double m_RMSChange{};
414 unsigned int m_NumberOfIterations{};
415
417 std::vector<FiniteDifferenceFunctionPointer> m_DifferenceFunctions{};
418
421 bool m_UseImageSpacing{};
422
423 void
424 PrintSelf(std::ostream & os, Indent indent) const override;
425
427 virtual void
429
433 virtual void
435
441 virtual TimeStepType
443
447 virtual void
449
453 void
454 GenerateData() override;
455
467 void
469
472 virtual bool
474
484 virtual bool
485 ThreadedHalt(void * itkNotUsed(threadInfo))
486 {
487 return this->Halt();
488 }
489
495 virtual void
497 {}
498
505 virtual void
507 {
508 for (IdCellType i = 0; i < this->m_FunctionCount; ++i)
509 {
510 this->m_DifferenceFunctions[i]->InitializeIteration();
511 }
512 }
528 inline TimeStepType
529 ResolveTimeStep(const TimeStepVectorType & timeStepList, const std::vector<uint8_t> & valid);
530
533 virtual void
535 {}
536
537private:
540 bool m_ManualReinitialization{};
541
543 bool m_InitializedState{};
544};
545} // end namespace itk
546
547#ifndef ITK_MANUAL_INSTANTIATION
548# include "itkMultiphaseFiniteDifferenceImageFilter.hxx"
549#endif
550
551#endif
A multi-dimensional iterator templated over image type that walks a region of pixels.
A multi-dimensional iterator templated over image type that walks a region of pixels.
void Set(const PixelType &value) const
Base class for all process objects that output image data.
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
Base class for filters that take an image as input and overwrite that image as the output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
virtual void ApplyUpdate(TimeStepType dt)=0
virtual void SetDifferenceFunction(const IdCellType &functionIndex, FiniteDifferenceFunctionPointer function)
void SetLevelSet(const IdCellType &i, const InputImageType *levelSet)
typename FiniteDifferenceFunctionType::Pointer FiniteDifferenceFunctionPointer
typename FiniteDifferenceFunctionType::TimeStepType TimeStepType
~MultiphaseFiniteDifferenceImageFilter() override=default
virtual const FiniteDifferenceFunctionPointer GetDifferenceFunction(const IdCellType &functionIndex) const
typename FiniteDifferenceFunctionType::RadiusType RadiusType
TimeStepType ResolveTimeStep(const TimeStepVectorType &timeStepList, const std::vector< uint8_t > &valid)
void PrintSelf(std::ostream &os, Indent indent) const override
virtual TimeStepType CalculateChange()=0
static constexpr T max(const T &)
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
This class generates a KdTree object without centroid information.
This class provides methods for k-nearest neighbor search and related data structures for a k-d tree.
Definition: itkKdTree.h:528
This class is the native implementation of the a Sample with an STL container.
Definition: itkListSample.h:52
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:63
static Pointer New()
SmartPointer< Self > Pointer
static constexpr double eps
Definition: itkMath.h:97
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long IndexValueType
Definition: itkIntTypes.h:90
unsigned long SizeValueType
Definition: itkIntTypes.h:83
long OffsetValueType
Definition: itkIntTypes.h:94