ITK  6.0.0
Insight Toolkit
SphinxExamples/src/Core/Common/MakePartOfImageTransparent/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
*
* 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.
*
*=========================================================================*/
#include "itkImage.h"
#include "itkTIFFImageIO.h"
#include "itkRGBAPixel.h"
int
main(int argc, char * argv[])
{
std::string outputFilename;
if (argc > 1)
{
outputFilename = argv[1];
}
else
{
outputFilename = "test.tif";
}
using PixelType = itk::RGBAPixel<unsigned char>;
using ImageType = itk::Image<PixelType, 2>;
start[0] = 0;
start[1] = 0;
size[0] = 200;
size[1] = 300;
region.SetSize(size);
region.SetIndex(start);
auto image = ImageType::New();
image->SetRegions(region);
image->Allocate();
itk::ImageRegionIterator<ImageType> imageIterator(image, region);
while (!imageIterator.IsAtEnd())
{
ImageType::PixelType pixel = imageIterator.Get();
if (imageIterator.GetIndex()[0] > 100)
{
pixel.SetRed(0);
pixel.SetGreen(255);
pixel.SetBlue(0);
// pixel.SetAlpha(255); // invisible
pixel.SetAlpha(122);
}
else
{
pixel.SetRed(255);
pixel.SetGreen(0);
pixel.SetBlue(0);
pixel.SetAlpha(static_cast<unsigned char>(0.5 * 255));
}
imageIterator.Set(pixel);
++imageIterator;
}
using WriterType = itk::ImageFileWriter<ImageType>;
using TIFFIOType = itk::TIFFImageIO;
auto writer = WriterType::New();
auto tiffIO = TIFFIOType::New();
tiffIO->SetPixelType(itk::IOPixelEnum::RGBA);
writer->SetFileName(outputFilename);
writer->SetInput(image);
writer->SetImageIO(tiffIO);
writer->Update();
return EXIT_SUCCESS;
}
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkRGBAPixel.h
itkImageRegionIterator.h
itk::ImageConstIterator::IsAtEnd
bool IsAtEnd() const
Definition: itkImageConstIterator.h:377
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:306
itk::ImageFileWriter
Writes image data to a single file.
Definition: itkImageFileWriter.h:90
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::TIFFImageIO
ImageIO object for reading and writing TIFF images.
Definition: itkTIFFImageIO.h:48
itk::ImageRegionIterator::Set
void Set(const PixelType &value) const
Definition: itkImageRegionIterator.h:116
itkImageFileWriter.h
itk::Size::SetSize
void SetSize(const SizeValueType val[VDimension])
Definition: itkSize.h:179
itk::RGBAPixel
Represent Red, Green, Blue and Alpha components for color images.
Definition: itkRGBAPixel.h:59
itkTIFFImageIO.h
itk::ImageConstIterator::Get
PixelType Get() const
Definition: itkImageConstIterator.h:336
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
New
static Pointer New()
itk::CommonEnums::IOPixel::RGBA