File: vtkvmtkOtsuMultipleThresholdsImageFilter.cxx

package info (click to toggle)
vmtk 1.3%2Bdfsg-2.1%2Bdeb9u1
  • links: PTS, VCS
  • area: non-free
  • in suites: stretch
  • size: 8,932 kB
  • sloc: cpp: 82,947; ansic: 31,817; python: 21,462; perl: 381; makefile: 93; ruby: 41; sh: 19
file content (110 lines) | stat: -rw-r--r-- 3,476 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
/*=========================================================================

Program:   VMTK
Module:    $RCSfile: vtkvmtkOtsuMultipleThresholdsImageFilter.cxx,v $
Language:  C++
Date:      $Date: 2006/04/06 16:48:25 $
Version:   $Revision: 1.1 $

  Copyright (c) Luca Antiga, David Steinman. All rights reserved.
  See LICENCE file for details.

  Portions of this code are covered under the VTK copyright.
  See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm 
  for details.

  Portions of this code are covered under the ITK copyright.
  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm
  for details.

     This software is distributed WITHOUT ANY WARRANTY; without even 
     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
     PURPOSE.  See the above copyright notices for more information.

=========================================================================*/

#include "vtkvmtkOtsuMultipleThresholdsImageFilter.h"
#include "vtkImageData.h"
#include "vtkFloatArray.h"
#include "vtkInformationVector.h"
#include "vtkInformation.h"
#include "vtkObjectFactory.h"

#include "vtkvmtkITKFilterUtilities.h"

#include "itkOtsuMultipleThresholdsImageFilter.h"


vtkStandardNewMacro(vtkvmtkOtsuMultipleThresholdsImageFilter);

vtkvmtkOtsuMultipleThresholdsImageFilter::vtkvmtkOtsuMultipleThresholdsImageFilter()
{
  this->NumberOfHistogramBins = 128;
  this->NumberOfThresholds = 1;
  this->LabelOffset = 0;
  this->Thresholds = NULL;
}

vtkvmtkOtsuMultipleThresholdsImageFilter::~vtkvmtkOtsuMultipleThresholdsImageFilter()
{
  if (this->Thresholds)
    {
    this->Thresholds->Delete();
    this->Thresholds = NULL;
    }
}

int vtkvmtkOtsuMultipleThresholdsImageFilter::RequestInformation (
  vtkInformation * vtkNotUsed(request),
  vtkInformationVector **inputVector,
  vtkInformationVector *outputVector)
{
  vtkInformation *outInfo = outputVector->GetInformationObject(0);
  vtkDataObject::SetPointDataActiveScalarInfo(outInfo, VTK_UNSIGNED_SHORT, -1);

  return 1;
}

void vtkvmtkOtsuMultipleThresholdsImageFilter::SimpleExecute(vtkImageData *input, vtkImageData *output)
{
  typedef float PixelType;
  typedef unsigned short OutputPixelType;
  const int Dimension = 3;
  typedef itk::Image<PixelType, Dimension> ImageType;
  typedef itk::Image<OutputPixelType, Dimension> OutputImageType;

  ImageType::Pointer inImage = ImageType::New();

  vtkvmtkITKFilterUtilities::VTKToITKImage<ImageType>(input,inImage);

  typedef itk::OtsuMultipleThresholdsImageFilter<ImageType, OutputImageType> OtsuFilterType;
  typedef OtsuFilterType::ThresholdVectorType ThresholdVectorType;

  OtsuFilterType::Pointer imageFilter = OtsuFilterType::New();
  imageFilter->SetInput(inImage);
  imageFilter->SetNumberOfHistogramBins(this->NumberOfHistogramBins);
  imageFilter->SetNumberOfThresholds(this->NumberOfThresholds);
  imageFilter->SetLabelOffset(this->LabelOffset);
  imageFilter->Update();

  const ThresholdVectorType& thresholds = imageFilter->GetThresholds();

  if (this->Thresholds)
    {
    this->Thresholds->Delete();
    this->Thresholds = NULL;
    }

  this->Thresholds = vtkFloatArray::New();
  this->Thresholds->SetNumberOfTuples(thresholds.size());

  for (unsigned long j=0; j<thresholds.size(); j++)
    {
    this->Thresholds->SetValue(j,thresholds[j]);
    }

  OutputImageType::Pointer outputImage = imageFilter->GetOutput();

  vtkvmtkITKFilterUtilities::ITKToVTKImage<OutputImageType>(outputImage,output);
}