ITK  5.2.0
Insight Toolkit
SphinxExamples/src/Core/Common/StoreNonPixelDataInImage/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.
*
*=========================================================================*/
#include <itkImage.h>
using ImageType = itk::Image<unsigned char, 2>;
static void
CreateImage(ImageType::Pointer image);
int
main(int, char *[])
{
// Create an image
ImageType::Pointer image = ImageType::New();
CreateImage(image);
// Store some data in it
itk::EncapsulateMetaData<float>(dictionary, "ASimpleFloat", 1.2);
image->SetMetaDataDictionary(dictionary);
// View all of the data
dictionary.Print(std::cout);
// View the data individually
auto itr = dictionary.Begin();
while (itr != dictionary.End())
{
std::cout << "Key = " << itr->first << std::endl;
std::cout << "Value = ";
itr->second->Print(std::cout);
std::cout << std::endl;
++itr;
}
// Write the image (and the data) to a file
using WriterType = itk::ImageFileWriter<ImageType>;
WriterType::Pointer writer = WriterType::New();
writer->SetFileName("test.mhd");
writer->SetInput(image);
writer->Update();
// Read the image (and data) from the file
using ReaderType = itk::ImageFileReader<ImageType>;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName("test.mhd");
// Display the data
std::cout << "Data read from file:" << std::endl;
reader->GetMetaDataDictionary().Print(std::cout);
return EXIT_SUCCESS;
}
void
CreateImage(ImageType::Pointer image)
{
start.Fill(0);
size.Fill(10);
region.SetSize(size);
region.SetIndex(start);
image->SetRegions(region);
image->Allocate();
itk::ImageRegionIterator<ImageType> imageIterator(image, image->GetLargestPossibleRegion());
while (!imageIterator.IsAtEnd())
{
imageIterator.Set(20);
++imageIterator;
}
}
itkImageFileReader.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itkImageRegionIterator.h
itk::Index::Fill
void Fill(IndexValueType value)
Definition: itkIndex.h:270
itk::Size::Fill
void Fill(SizeValueType value)
Definition: itkSize.h:211
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::MetaDataDictionary::Print
virtual void Print(std::ostream &os) const
itk::MetaDataDictionary
Provides a mechanism for storing a collection of arbitrary data types.
Definition: itkMetaDataDictionary.h:53
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::MetaDataDictionary::Begin
Iterator Begin()
itk::ImageRegionIterator::Set
void Set(const PixelType &value) const
Definition: itkImageRegionIterator.h:116
itkImageFileWriter.h
itkMetaDataObject.h
itkMetaDataDictionary.h
itk::MetaDataDictionary::End
Iterator End()
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:174