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 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230
|
/*
* Copyright (C) 2005-2017 Centre National d'Etudes Spatiales (CNES)
*
* This file is part of Orfeo Toolbox
*
* https://www.orfeo-toolbox.org/
*
* 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
*
* 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 : BeginLatex
//
// This example illustrates the use of the
// \doxygen{otb}{MorphologicalSegmentationFilter}. This filter
// performs a segmentation of the details \code{supFilter} and
// \code{infFilter} extracted with the morphological pyramid. The
// segmentation algorithm used is based on seeds extraction using the
// \doxygen{otb}{ImageToPointSetFilter}, followed by a connected
// threshold segmentation using the
// \doxygen{itk}{ConnectedThresholdImageFilter}. The threshold for seeds
// extraction and segmentation are computed using quantiles. A pre
// processing step is applied by multiplying the full resolution
// brighter details (resp. darker details) with the original image
// (resp. the inverted original image). This performs an enhancement of
// the regions contour precision. The details from the pyramid are set
// via the \code{SetBrighterDetails()} and \code{SetDarkerDetails()}
// methods. The brighter and darker details depend on the filter used
// in the pyramid analysis. If the
// \doxygen{otb}{OpeningClosingMorphologicalFilter} filter is used,
// then the brighter details are those from the \code{supFilter} image
// list, whereas if the
// \doxygen{otb}{ClosingOpeningMorphologicalFilter} filter is used,
// the brighter details are those from the \code{infFilter} list. The
// output of the segmentation filter is a single segmentation images
// list, containing first the brighter details segmentation from
// higher scale to lower, and then the darker details in the same
// order. The attention of the user is drawn to the fact that since
// the label filter used internally will deal with a large number of
// labels, the \code{OutputPixelType} is required to be sufficiently
// precise. Unsigned short or Unsigned long would be a good choice,
// unless the user has a very good reason to think that a less precise
// type will be sufficient. The first step to use this filter is to
// include its header file.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "otbMorphologicalPyramidSegmentationFilter.h"
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The mathematical morphology filters to be used have also to be
// included here, as well as the morphological pyramid analysis filter.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
#include "otbOpeningClosingMorphologicalFilter.h"
#include "itkBinaryBallStructuringElement.h"
#include "otbMorphologicalPyramidAnalysisFilter.h"
// Software Guide : EndCodeSnippet
#include "otbImageFileReader.h"
#include "otbImageFileWriter.h"
#include "otbImage.h"
#include "itkMacro.h"
int main(int itkNotUsed(argc), char * argv[])
{
const char* inputFilename = argv[1];
const char* outputFilenamePrefix = argv[2];
const char * outputFilenameSuffix = argv[3];
const unsigned int numberOfLevels = atoi(argv[4]);
const double decimationRatio = atof(argv[5]);
const float seedsQuantile = atof(argv[6]);
const float segmentationQuantile = atof(argv[7]);
const unsigned int minObjectSize = atoi(argv[8]);
// Software Guide : BeginLatex
//
// As usual, we start by defining the types for the pixels, the
// images, the reader and the writer. We also define the types needed
// for the morphological pyramid analysis.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
const unsigned int Dimension = 2;
typedef unsigned char InputPixelType;
typedef unsigned short OutputPixelType;
typedef otb::Image<InputPixelType, Dimension> InputImageType;
typedef otb::Image<OutputPixelType, Dimension> OutputImageType;
typedef otb::ImageFileReader<InputImageType> ReaderType;
typedef otb::ImageFileWriter<OutputImageType> WriterType;
typedef itk::BinaryBallStructuringElement<InputPixelType, Dimension>
StructuringElementType;
typedef otb::OpeningClosingMorphologicalFilter<InputImageType,
InputImageType,
StructuringElementType>
OpeningClosingFilterType;
typedef otb::MorphologicalPyramidAnalysisFilter<InputImageType,
InputImageType,
OpeningClosingFilterType>
PyramidFilterType;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We can now define the type for the
// \doxygen{otb}{MorphologicalPyramidSegmentationFilter} which is
// templated over the input and output image types.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef otb::MorphologicalPyramidSegmentationFilter<InputImageType,
OutputImageType>
SegmentationFilterType;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Since the output of the segmentation filter is a list of images,
// we define an iterator type which will be used to access the
// segmented images.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
typedef SegmentationFilterType::OutputImageListIteratorType
OutputListIteratorType;
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The following code snippet shows how to read the input image and
// perform the morphological pyramid analysis.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(inputFilename);
PyramidFilterType::Pointer pyramid = PyramidFilterType::New();
pyramid->SetNumberOfLevels(numberOfLevels);
pyramid->SetDecimationRatio(decimationRatio);
pyramid->SetInput(reader->GetOutput());
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// We can now instantiate the segmentation filter and set its
// parameters. As one can see, the \code{SetReferenceImage()} is used to
// pass the original image in order to obtain sharp region
// boundaries. Using the \code{SetBrighterDetails()} and
// \code{SetDarkerDetails()} the output of the analysis is passed to the
// filter. Finally, the parameters for the segmentation are set by
// using the \code{SetSeedsQuantile()},
// \code{SetConnectedThresholdQuantile()} and
// \code{SetMinimumObjectSize()} methods.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
SegmentationFilterType::Pointer segmentation = SegmentationFilterType::New();
segmentation->SetReferenceImage(reader->GetOutput());
segmentation->SetBrighterDetails(pyramid->GetSupFilter());
segmentation->SetDarkerDetails(pyramid->GetInfFilter());
segmentation->SetSeedsQuantile(seedsQuantile);
segmentation->SetConnectedThresholdQuantile(segmentationQuantile);
segmentation->SetMinimumObjectSize(minObjectSize);
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The pipeline is executed bu calling the \code{Update()} method.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
segmentation->Update();
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// Finally, we get an iterator to the list generated as output for the
// segmentation and we use it to iterate through the list and write
// the images to files.
//
// Software Guide : EndLatex
// Software Guide : BeginCodeSnippet
OutputListIteratorType it = segmentation->GetOutput()->Begin();
WriterType::Pointer writer;
int index = 1;
std::stringstream oss;
while (it != segmentation->GetOutput()->End())
{
oss << outputFilenamePrefix << index << "." << outputFilenameSuffix;
writer = WriterType::New();
writer->SetInput(it.Get());
writer->SetFileName(oss.str().c_str());
writer->Update();
std::cout << oss.str() << " file written." << std::endl;
oss.str("");
++index;
++it;
}
// Software Guide : EndCodeSnippet
// Software Guide : BeginLatex
//
// The user will pay attention to the fact that the list contains
// first the brighter details segmentation from
// higher scale to lower, and then the darker details in the same
// order.
//
// Software Guide : EndLatex
return EXIT_SUCCESS;
}
|