File: ConfidenceConnected3D.cxx

package info (click to toggle)
insighttoolkit5 5.4.5-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 704,592 kB
  • sloc: cpp: 784,579; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,934; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 460; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (137 lines) | stat: -rw-r--r-- 3,945 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
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/

//  Software Guide : BeginCommandLineArgs
//  INPUTS:  {brainweb165a10f17.mha}
//  ARGUMENTS: {WhiteMatterSegmentation.mhd}
//  Software Guide : EndCommandLineArgs


#include "itkConfidenceConnectedImageFilter.h"
#include "itkCastImageFilter.h"
#include "itkCurvatureFlowImageFilter.h"
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"

// Software Guide : BeginLatex
//
// This example is a 3D version of the previous ConfidenceConnected example.
// In this particular case, we are extracting the white matter from an input
// Brain MRI dataset.
//
// Software Guide : EndLatex


int
main(int argc, char * argv[])
{
  if (argc < 3)
  {
    std::cerr << "Missing Parameters " << std::endl;
    std::cerr << "Usage: " << argv[0];
    std::cerr << " inputImage  outputImage " << std::endl;
    return EXIT_FAILURE;
  }


  using InternalPixelType = float;
  constexpr unsigned int Dimension = 3;
  using InternalImageType = itk::Image<InternalPixelType, Dimension>;

  using OutputPixelType = unsigned char;
  using OutputImageType = itk::Image<OutputPixelType, Dimension>;

  using CastingFilterType =
    itk::CastImageFilter<InternalImageType, OutputImageType>;
  auto caster = CastingFilterType::New();


  using ReaderType = itk::ImageFileReader<InternalImageType>;
  using WriterType = itk::ImageFileWriter<OutputImageType>;

  auto reader = ReaderType::New();
  auto writer = WriterType::New();

  reader->SetFileName(argv[1]);
  writer->SetFileName(argv[2]);

  using CurvatureFlowImageFilterType =
    itk::CurvatureFlowImageFilter<InternalImageType, InternalImageType>;
  auto smoothing = CurvatureFlowImageFilterType::New();

  using ConnectedFilterType =
    itk::ConfidenceConnectedImageFilter<InternalImageType, InternalImageType>;
  auto confidenceConnected = ConnectedFilterType::New();

  smoothing->SetInput(reader->GetOutput());
  confidenceConnected->SetInput(smoothing->GetOutput());
  caster->SetInput(confidenceConnected->GetOutput());
  writer->SetInput(caster->GetOutput());

  smoothing->SetNumberOfIterations(2);
  smoothing->SetTimeStep(0.05);

  confidenceConnected->SetMultiplier(2.5);
  confidenceConnected->SetNumberOfIterations(5);
  confidenceConnected->SetInitialNeighborhoodRadius(2);
  confidenceConnected->SetReplaceValue(255);

  InternalImageType::IndexType index1;
  index1[0] = 118;
  index1[1] = 133;
  index1[2] = 92;
  confidenceConnected->AddSeed(index1);

  InternalImageType::IndexType index2;
  index2[0] = 63;
  index2[1] = 135;
  index2[2] = 94;
  confidenceConnected->AddSeed(index2);

  InternalImageType::IndexType index3;
  index3[0] = 63;
  index3[1] = 157;
  index3[2] = 90;
  confidenceConnected->AddSeed(index3);

  InternalImageType::IndexType index4;
  index4[0] = 111;
  index4[1] = 150;
  index4[2] = 90;
  confidenceConnected->AddSeed(index4);

  InternalImageType::IndexType index5;
  index5[0] = 111;
  index5[1] = 50;
  index5[2] = 88;
  confidenceConnected->AddSeed(index5);

  try
  {
    writer->Update();
  }
  catch (const itk::ExceptionObject & excep)
  {
    std::cerr << "Exception caught !" << std::endl;
    std::cerr << excep << std::endl;
    return EXIT_FAILURE;
  }


  return EXIT_SUCCESS;
}