|
| 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 | + |
| 19 | +// Software Guide : BeginCommandLineArgs |
| 20 | +// INPUTS: {BrainProtonDensitySlice.png} |
| 21 | +// ARGUMENTS: {DerivativeImageFilterFloatOutput.mhd} |
| 22 | +// OUTPUTS: {DerivativeImageFilterOutput.png} |
| 23 | +// ARGUMENTS: 1 0 |
| 24 | +// Software Guide : EndCommandLineArgs |
| 25 | + |
| 26 | +// Software Guide : BeginLatex |
| 27 | +// |
| 28 | +// The \doxygen{DerivativeImageFilter} is used for computing the partial |
| 29 | +// derivative of an image, the derivative of an image along a particular |
| 30 | +// axial direction. |
| 31 | +// |
| 32 | +// \index{itk::DerivativeImageFilter} |
| 33 | +// |
| 34 | +// Software Guide : EndLatex |
| 35 | + |
| 36 | + |
| 37 | +#include "itkImage.h" |
| 38 | +#include "itkImageFileReader.h" |
| 39 | +#include "itkImageFileWriter.h" |
| 40 | +#include "itkRescaleIntensityImageFilter.h" |
| 41 | + |
| 42 | +// Software Guide : BeginLatex |
| 43 | +// |
| 44 | +// The header file corresponding to this filter should be included first. |
| 45 | +// |
| 46 | +// \index{itk::DerivativeImageFilter!header} |
| 47 | +// |
| 48 | +// Software Guide : EndLatex |
| 49 | + |
| 50 | + |
| 51 | +// Software Guide : BeginCodeSnippet |
| 52 | +#include "itkDerivativeImageFilter.h" |
| 53 | +// Software Guide : EndCodeSnippet |
| 54 | + |
| 55 | + |
| 56 | +int |
| 57 | +main(int argc, char * argv[]) |
| 58 | +{ |
| 59 | + if (argc < 6) |
| 60 | + { |
| 61 | + std::cerr << "Usage: " << std::endl; |
| 62 | + std::cerr << argv[0] << " inputImageFile outputImageFile normalizedOutputImageFile "; |
| 63 | + std::cerr << " derivativeOrder direction" << std::endl; |
| 64 | + return EXIT_FAILURE; |
| 65 | + } |
| 66 | + |
| 67 | + |
| 68 | + // Software Guide : BeginLatex |
| 69 | + // |
| 70 | + // Next, the pixel types for the input and output images must be defined |
| 71 | + // and, with them, the image types can be instantiated. Note that it is |
| 72 | + // important to select a signed type for the image, since the values of the |
| 73 | + // derivatives will be positive as well as negative. |
| 74 | + // |
| 75 | + // Software Guide : EndLatex |
| 76 | + |
| 77 | + // Software Guide : BeginCodeSnippet |
| 78 | + using InputPixelType = float; |
| 79 | + using OutputPixelType = float; |
| 80 | + |
| 81 | + constexpr unsigned int Dimension = 2; |
| 82 | + |
| 83 | + using InputImageType = itk::Image<InputPixelType, Dimension>; |
| 84 | + using OutputImageType = itk::Image<OutputPixelType, Dimension>; |
| 85 | + // Software Guide : EndCodeSnippet |
| 86 | + |
| 87 | + |
| 88 | + using ReaderType = itk::ImageFileReader<InputImageType>; |
| 89 | + using WriterType = itk::ImageFileWriter<OutputImageType>; |
| 90 | + |
| 91 | + auto reader = ReaderType::New(); |
| 92 | + auto writer = WriterType::New(); |
| 93 | + |
| 94 | + reader->SetFileName(argv[1]); |
| 95 | + writer->SetFileName(argv[2]); |
| 96 | + |
| 97 | + // Software Guide : BeginLatex |
| 98 | + // |
| 99 | + // Using the image types, it is now possible to define the filter type |
| 100 | + // and create the filter object. |
| 101 | + // |
| 102 | + // \index{itk::DerivativeImageFilter!instantiation} |
| 103 | + // \index{itk::DerivativeImageFilter!New()} |
| 104 | + // \index{itk::DerivativeImageFilter!Pointer} |
| 105 | + // |
| 106 | + // Software Guide : EndLatex |
| 107 | + |
| 108 | + // Software Guide : BeginCodeSnippet |
| 109 | + using FilterType = itk::DerivativeImageFilter<InputImageType, OutputImageType>; |
| 110 | + |
| 111 | + auto filter = FilterType::New(); |
| 112 | + // Software Guide : EndCodeSnippet |
| 113 | + |
| 114 | + |
| 115 | + // Software Guide : BeginLatex |
| 116 | + // |
| 117 | + // The order of the derivative is selected with the \code{SetOrder()} |
| 118 | + // method. The direction along which the derivative will be computed is |
| 119 | + // selected with the \code{SetDirection()} method. |
| 120 | + // |
| 121 | + // \index{itk::DerivativeImageFilter!SetOrder()} |
| 122 | + // \index{itk::DerivativeImageFilter!SetDirection()} |
| 123 | + // |
| 124 | + // Software Guide : EndLatex |
| 125 | + |
| 126 | + // Software Guide : BeginCodeSnippet |
| 127 | + filter->SetOrder(std::stoi(argv[4])); |
| 128 | + filter->SetDirection(std::stoi(argv[5])); |
| 129 | + // Software Guide : EndCodeSnippet |
| 130 | + |
| 131 | + |
| 132 | + // Software Guide : BeginLatex |
| 133 | + // |
| 134 | + // The input to the filter can be taken from any other filter, for example |
| 135 | + // a reader. The output can be passed down the pipeline to other filters, |
| 136 | + // for example, a writer. An \code{Update()} call on any downstream filter |
| 137 | + // will trigger the execution of the derivative filter. |
| 138 | + // |
| 139 | + // \index{itk::DerivativeImageFilter!SetInput()} |
| 140 | + // \index{itk::DerivativeImageFilter!GetOutput()} |
| 141 | + // |
| 142 | + // Software Guide : EndLatex |
| 143 | + |
| 144 | + |
| 145 | + // Software Guide : BeginCodeSnippet |
| 146 | + filter->SetInput(reader->GetOutput()); |
| 147 | + writer->SetInput(filter->GetOutput()); |
| 148 | + writer->Update(); |
| 149 | + // Software Guide : EndCodeSnippet |
| 150 | + |
| 151 | + |
| 152 | + // Software Guide : BeginLatex |
| 153 | + // |
| 154 | + // \begin{figure} |
| 155 | + // \center |
| 156 | + // \includegraphics[width=0.44\textwidth]{BrainProtonDensitySlice} |
| 157 | + // \includegraphics[width=0.44\textwidth]{DerivativeImageFilterOutput} |
| 158 | + // \itkcaption[Effect of the Derivative filter.]{Effect of the Derivative |
| 159 | + // filter on a slice from a MRI proton density brain image.} |
| 160 | + // \label{fig:DerivativeImageFilterOutput} |
| 161 | + // \end{figure} |
| 162 | + // |
| 163 | + // Figure \ref{fig:DerivativeImageFilterOutput} illustrates the effect of |
| 164 | + // the DerivativeImageFilter on a slice of MRI brain image. The derivative |
| 165 | + // is taken along the $x$ direction. The sensitivity to noise in the image |
| 166 | + // is evident from this result. |
| 167 | + // |
| 168 | + // Software Guide : EndLatex |
| 169 | + |
| 170 | + |
| 171 | + using WriteImageType = itk::Image<unsigned char, Dimension>; |
| 172 | + |
| 173 | + using NormalizeFilterType = itk::RescaleIntensityImageFilter<OutputImageType, WriteImageType>; |
| 174 | + |
| 175 | + using NormalizedWriterType = itk::ImageFileWriter<WriteImageType>; |
| 176 | + |
| 177 | + auto normalizer = NormalizeFilterType::New(); |
| 178 | + auto normalizedWriter = NormalizedWriterType::New(); |
| 179 | + |
| 180 | + normalizer->SetInput(filter->GetOutput()); |
| 181 | + normalizedWriter->SetInput(normalizer->GetOutput()); |
| 182 | + |
| 183 | + normalizer->SetOutputMinimum(0); |
| 184 | + normalizer->SetOutputMaximum(255); |
| 185 | + |
| 186 | + normalizedWriter->SetFileName(argv[3]); |
| 187 | + normalizedWriter->Update(); |
| 188 | + |
| 189 | + return EXIT_SUCCESS; |
| 190 | +} |
0 commit comments