File: VesselnessMeasureImageFilter.cxx

package info (click to toggle)
insighttoolkit5 5.4.5-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 704,588 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: 461; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (78 lines) | stat: -rw-r--r-- 2,653 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
/*=========================================================================
 *
 *  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.
 *
 *=========================================================================*/

// The example takes an image (say MRA image), computes the vesselness measure
// of the image using the HessianRecursiveGaussianImageFilter and the
// Hessian3DToVesselnessMeasureImageFilter. The goal is to detect bright
// tubular structures in the image.

#include "itkHessian3DToVesselnessMeasureImageFilter.h"
#include "itkHessianRecursiveGaussianImageFilter.h"
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"

int
main(int argc, char * argv[])
{
  if (argc < 3)
  {
    std::cerr << "Usage: inputImage outputImage [sigma] [alpha_1] [alpha_2]"
              << std::endl;
  }

  constexpr unsigned int Dimension = 3;
  using InputPixelType = double;
  using OutputPixelType = float;

  using InputImageType = itk::Image<InputPixelType, Dimension>;
  using OutputImageType = itk::Image<OutputPixelType, Dimension>;

  using HessianFilterType =
    itk::HessianRecursiveGaussianImageFilter<InputImageType>;
  using VesselnessMeasureFilterType =
    itk::Hessian3DToVesselnessMeasureImageFilter<OutputPixelType>;
  using ReaderType = itk::ImageFileReader<InputImageType>;
  using WriterType = itk::ImageFileWriter<OutputImageType>;

  auto hessianFilter = HessianFilterType::New();
  auto vesselnessFilter = VesselnessMeasureFilterType::New();

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

  reader->SetFileName(argv[1]);
  hessianFilter->SetInput(reader->GetOutput());
  if (argc >= 4)
  {
    hessianFilter->SetSigma(static_cast<double>(std::stod(argv[3])));
  }
  vesselnessFilter->SetInput(hessianFilter->GetOutput());
  writer->SetInput(vesselnessFilter->GetOutput());
  writer->SetFileName(argv[2]);
  if (argc >= 5)
  {
    vesselnessFilter->SetAlpha1(static_cast<double>(std::stod(argv[4])));
  }
  if (argc >= 6)
  {
    vesselnessFilter->SetAlpha2(static_cast<double>(std::stod(argv[5])));
  }

  writer->Update();
  return EXIT_SUCCESS;
}