File: DerivativeImageFilter.cxx

package info (click to toggle)
insighttoolkit4 4.13.3withdata-dfsg1-4
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 489,260 kB
  • sloc: cpp: 557,342; ansic: 146,850; fortran: 34,788; python: 16,572; sh: 2,187; lisp: 2,070; tcl: 993; java: 362; perl: 200; makefile: 129; csh: 81; pascal: 69; xml: 19; ruby: 10
file content (192 lines) | stat: -rw-r--r-- 6,302 bytes parent folder | download | duplicates (5)
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
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
/*=========================================================================
 *
 *  Copyright Insight Software Consortium
 *
 *  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
 *
 *         http://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:  {BrainProtonDensitySlice.png}
//    ARGUMENTS: {DerivativeImageFilterFloatOutput.mhd}
//    OUTPUTS: {DerivativeImageFilterOutput.png}
//    ARGUMENTS:    1 0
//  Software Guide : EndCommandLineArgs

//  Software Guide : BeginLatex
//
//  The \doxygen{DerivativeImageFilter} is used for computing the partial
//  derivative of an image, the derivative of an image along a particular axial
//  direction.
//
//  \index{itk::DerivativeImageFilter}
//
//  Software Guide : EndLatex


#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"
#include "itkRescaleIntensityImageFilter.h"

//  Software Guide : BeginLatex
//
//  The header file corresponding to this filter should be included first.
//
//  \index{itk::DerivativeImageFilter!header}
//
//  Software Guide : EndLatex


// Software Guide : BeginCodeSnippet
#include "itkDerivativeImageFilter.h"
// Software Guide : EndCodeSnippet


int main( int argc, char * argv[] )
{
  if( argc < 6 )
    {
    std::cerr << "Usage: " << std::endl;
    std::cerr << argv[0] << "  inputImageFile   outputImageFile  normalizedOutputImageFile ";
    std::cerr << " derivativeOrder direction" << std::endl;
    return EXIT_FAILURE;
    }


  //  Software Guide : BeginLatex
  //
  //  Next, the pixel types for the input and output images must be defined and, with
  //  them, the image types can be instantiated. Note that it is important to
  //  select a signed type for the image, since the values of the derivatives
  //  will be positive as well as negative.
  //
  //  Software Guide : EndLatex

  // Software Guide : BeginCodeSnippet
  typedef   float  InputPixelType;
  typedef   float  OutputPixelType;

  const unsigned int Dimension = 2;

  typedef itk::Image< InputPixelType,  Dimension >   InputImageType;
  typedef itk::Image< OutputPixelType, Dimension >   OutputImageType;
  // Software Guide : EndCodeSnippet


  typedef itk::ImageFileReader< InputImageType  >  ReaderType;
  typedef itk::ImageFileWriter< OutputImageType >  WriterType;

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

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

  //  Software Guide : BeginLatex
  //
  //  Using the image types, it is now possible to define the filter type
  //  and create the filter object.
  //
  //  \index{itk::DerivativeImageFilter!instantiation}
  //  \index{itk::DerivativeImageFilter!New()}
  //  \index{itk::DerivativeImageFilter!Pointer}
  //
  //  Software Guide : EndLatex

  // Software Guide : BeginCodeSnippet
  typedef itk::DerivativeImageFilter<
               InputImageType, OutputImageType >  FilterType;

  FilterType::Pointer filter = FilterType::New();
  // Software Guide : EndCodeSnippet


  //  Software Guide : BeginLatex
  //
  //  The order of the derivative is selected with the \code{SetOrder()}
  //  method.  The direction along which the derivative will be computed is
  //  selected with the \code{SetDirection()} method.
  //
  //  \index{itk::DerivativeImageFilter!SetOrder()}
  //  \index{itk::DerivativeImageFilter!SetDirection()}
  //
  //  Software Guide : EndLatex

  // Software Guide : BeginCodeSnippet
  filter->SetOrder(     atoi( argv[4] ) );
  filter->SetDirection( atoi( argv[5] ) );
  // Software Guide : EndCodeSnippet


  //  Software Guide : BeginLatex
  //
  //  The input to the filter can be taken from any other filter, for example
  //  a reader. The output can be passed down the pipeline to other filters,
  //  for example, a writer. An \code{Update()} call on any downstream filter will
  //  trigger the execution of the derivative filter.
  //
  //  \index{itk::DerivativeImageFilter!SetInput()}
  //  \index{itk::DerivativeImageFilter!GetOutput()}
  //
  //  Software Guide : EndLatex


  // Software Guide : BeginCodeSnippet
  filter->SetInput( reader->GetOutput() );
  writer->SetInput( filter->GetOutput() );
  writer->Update();
  // Software Guide : EndCodeSnippet


  //  Software Guide : BeginLatex
  //
  // \begin{figure}
  // \center
  // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySlice}
  // \includegraphics[width=0.44\textwidth]{DerivativeImageFilterOutput}
  // \itkcaption[Effect of the Derivative filter.]{Effect of the Derivative filter
  // on a slice from a MRI proton density brain image.}
  // \label{fig:DerivativeImageFilterOutput}
  // \end{figure}
  //
  //  Figure \ref{fig:DerivativeImageFilterOutput} illustrates the effect of
  //  the DerivativeImageFilter on a slice of MRI brain image. The derivative
  //  is taken along the $x$ direction.  The sensitivity to noise in the image
  //  is evident from this result.
  //
  //  Software Guide : EndLatex


  typedef itk::Image< unsigned char, Dimension >  WriteImageType;

  typedef itk::RescaleIntensityImageFilter<
                                  OutputImageType,
                                  WriteImageType >    NormalizeFilterType;

  typedef itk::ImageFileWriter< WriteImageType >       NormalizedWriterType;

  NormalizeFilterType::Pointer normalizer = NormalizeFilterType::New();
  NormalizedWriterType::Pointer normalizedWriter = NormalizedWriterType::New();

  normalizer->SetInput( filter->GetOutput() );
  normalizedWriter->SetInput( normalizer->GetOutput() );

  normalizer->SetOutputMinimum(   0 );
  normalizer->SetOutputMaximum( 255 );

  normalizedWriter->SetFileName( argv[3] );
  normalizedWriter->Update();

  return EXIT_SUCCESS;
}