ITK  5.2.0
Insight Toolkit
SphinxExamples/src/Core/Common/UseParallelizeImageRegion/Code.cxx
/*=========================================================================
*
* 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
*
* 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.
*
*=========================================================================*/
constexpr unsigned int Dimension = 2;
using InputPixelType = unsigned char;
using OutputPixelType = itk::RGBAPixel<unsigned char>;
OutputImageType::Pointer
segmentationAndCustomColorization(InputImageType::Pointer inImage)
{
using WatershedFilterType = itk::WatershedImageFilter<InputImageType>;
WatershedFilterType::Pointer watershed = WatershedFilterType::New();
watershed->SetThreshold(0.05);
watershed->SetLevel(0.3);
watershed->SetInput(inImage);
watershed->Update();
LabeledImageType::Pointer image = watershed->GetOutput();
image->DisconnectPipeline();
OutputImageType::Pointer outImage = OutputImageType::New();
outImage->CopyInformation(image);
outImage->SetRegions(image->GetBufferedRegion());
outImage->Allocate(true);
// ParallelizeImageRegion invokes the provided lambda function in parallel
// each invocation will contain a piece of the overall region
mt->ParallelizeImageRegion<Dimension>(
image->GetBufferedRegion(),
// the lambda will have access to outer variables 'image' and 'outImage'
// it will have parameter 'region', which needs to be processed
[image, outImage](const LabeledImageType::RegionType & region) {
for (; !iIt.IsAtEnd(); ++iIt, ++oIt)
{
OutputPixelType p;
p.SetAlpha(iIt.Get());
static_assert(Dimension <= 3, "Dimension has to be 2 or 3");
for (unsigned d = 0; d < Dimension; d++)
{
p.SetElement(d, ind[d]);
}
oIt.Set(p);
}
},
nullptr); // we don't have a filter whose progress needs to be updated
return outImage;
}
int
main(int argc, char * argv[])
{
if (argc != 3)
{
std::cerr << "Usage: " << std::endl;
std::cerr << argv[0];
std::cerr << " <InputFileName>";
std::cerr << " <OutputFileName>" << std::endl;
return EXIT_FAILURE;
}
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(argv[1]);
WriterType::Pointer writer = WriterType::New();
writer->SetFileName(argv[2]);
writer->UseCompressionOn();
try
{
OutputImageType::Pointer outImage = segmentationAndCustomColorization(reader->GetOutput());
writer->SetInput(outImage);
writer->Update();
}
catch (itk::ExceptionObject & error)
{
std::cerr << "Error: " << error << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
itk::MultiThreaderBase::New
static Pointer New()
itkImageFileReader.h
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::WatershedImageFilter
A low-level image analysis algorithm that automatically produces a hierarchy of segmented,...
Definition: itkWatershedImageFilter.h:149
itk::ImageConstIterator::IsAtEnd
bool IsAtEnd() const
Definition: itkImageConstIterator.h:384
itk::ImageFileReader
Data source that reads image data from a single file.
Definition: itkImageFileReader.h:75
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageConstIterator::GetIndex
const IndexType GetIndex() const
Definition: itkImageConstIterator.h:313
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:88
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageRegionIterator::Set
void Set(const PixelType &value) const
Definition: itkImageRegionIterator.h:116
itkImageFileWriter.h
itkWatershedImageFilter.h
itk::RGBAPixel
Represent Red, Green, Blue and Alpha components for color images.
Definition: itkRGBAPixel.h:59
itk::ImageConstIterator::Get
PixelType Get() const
Definition: itkImageConstIterator.h:343
itk::ImageRegionConstIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionConstIterator.h:109
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44