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 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273
|
/*=========================================================================
*
* 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.
*
*=========================================================================*/
#ifndef itkImageFileWriter_h
#define itkImageFileWriter_h
#include "ITKIOImageBaseExport.h"
#include "itkProcessObject.h"
#include "itkImageIOBase.h"
#include "itkMacro.h"
#include "itkMetaProgrammingLibrary.h"
namespace itk
{
/** \brief Base exception class for IO problems during writing.
*
* \class ImageFileWriterException
* \ingroup ITKIOImageBase
*/
class ITKIOImageBase_EXPORT ImageFileWriterException : public ExceptionObject
{
public:
ITK_DEFAULT_COPY_AND_MOVE(ImageFileWriterException);
/** \see LightObject::GetNameOfClass() */
itkOverrideGetNameOfClassMacro(ImageFileWriterException);
/** Constructor. */
ImageFileWriterException(const char * file,
unsigned int line,
const char * message = "Error in IO",
const char * loc = "Unknown")
: ExceptionObject(file, line, message, loc)
{}
/** Constructor. */
ImageFileWriterException(const std::string & file,
unsigned int line,
const char * message = "Error in IO",
const char * loc = "Unknown")
: ExceptionObject(file, line, message, loc)
{}
/** Has to have empty throw(). */
~ImageFileWriterException() noexcept override;
};
/** \class ImageFileWriter
* \brief Writes image data to a single file.
*
* ImageFileWriter writes its input data to a single output file.
* ImageFileWriter interfaces with an ImageIO class to write out the
* data. If you wish to write data into a series of files (e.g., a
* slice per file) use ImageSeriesWriter.
*
* A plugable factory pattern is used that allows different kinds of writers
* to be registered (even at run time) without having to modify the
* code in this class. You can either manually instantiate the ImageIO
* object and associate it with the ImageFileWriter, or let the class
* figure it out from the extension. Normally just setting the filename
* with a suitable suffix (".png", ".jpg", etc) and setting the input
* to the writer is enough to get the writer to work properly.
*
* \sa ImageSeriesReader
* \sa ImageIOBase
*
* \ingroup IOFilters
* \ingroup ITKIOImageBase
*
* \sphinx
* \sphinxexample{IO/ImageBase/WriteAnImage,Write An image}
* \endsphinx
*/
template <typename TInputImage>
class ITK_TEMPLATE_EXPORT ImageFileWriter : public ProcessObject
{
public:
ITK_DISALLOW_COPY_AND_MOVE(ImageFileWriter);
/** Standard class type aliases. */
using Self = ImageFileWriter;
using Superclass = ProcessObject;
using Pointer = SmartPointer<Self>;
using ConstPointer = SmartPointer<const Self>;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** \see LightObject::GetNameOfClass() */
itkOverrideGetNameOfClassMacro(ImageFileWriter);
/** Some convenient type alias. */
using InputImageType = TInputImage;
using InputImagePointer = typename InputImageType::Pointer;
using InputImageRegionType = typename InputImageType::RegionType;
using InputImagePixelType = typename InputImageType::PixelType;
/** Set/Get the image input of this writer. */
using Superclass::SetInput;
void
SetInput(const InputImageType * input);
const InputImageType *
GetInput();
const InputImageType *
GetInput(unsigned int idx);
/** Specify the name of the output file to write. */
itkSetStringMacro(FileName);
itkGetStringMacro(FileName);
/** Set/Get the ImageIO helper class. Usually this is created via the object
* factory mechanism that determines whether a particular ImageIO can
* write a certain file. This method provides a way to get the ImageIO
* instance that is created, or one can be manually set where the
* IO factory mechanism may not work (for example, raw image files or
* image files with non-standard filename suffix's.
* If the user specifies the ImageIO, we assume she makes the
* correct choice and will allow a file to be created regardless of
* the file extension. If the factory has set the ImageIO, the
* extension must be supported by the specified ImageIO. */
void
SetImageIO(ImageIOBase * io)
{
if (this->m_ImageIO != io)
{
this->Modified();
this->m_ImageIO = io;
}
m_FactorySpecifiedImageIO = false;
}
itkGetModifiableObjectMacro(ImageIO, ImageIOBase);
/** A special version of the Update() method for writers. It
* invokes start and end events and handles releasing data. It
* eventually calls GenerateData() which does the actual writing.
* Note: the write method will write data specified by the
* IORegion. If not set, then then the whole image is written. Note
* that the region will be cropped to fit the input image's
* LargestPossibleRegion. */
virtual void
Write();
/** Specify the region to write. If left nullptr, then the whole image
* is written. */
void
SetIORegion(const ImageIORegion & region);
const ImageIORegion &
GetIORegion() const
{
return m_PasteIORegion;
}
/** Set/Get the number of pieces to divide the input. The upstream pipeline
* will try to be executed this many times. */
itkSetMacro(NumberOfStreamDivisions, unsigned int);
itkGetConstReferenceMacro(NumberOfStreamDivisions, unsigned int);
/** Aliased to the Write() method to be consistent with the rest of the
* pipeline. */
void
Update() override
{
this->Write();
}
/** \brief Writes the entire image to file.
*
* Updates the pipeline, streaming it the NumberOfStreamDivisions times.
* Existing PasteIORegion is reset.
*/
void
UpdateLargestPossibleRegion() override
{
m_PasteIORegion = ImageIORegion(TInputImage::ImageDimension);
m_UserSpecifiedIORegion = false;
this->Write();
}
/** Set the compression On or Off */
itkSetMacro(UseCompression, bool);
itkGetConstReferenceMacro(UseCompression, bool);
itkBooleanMacro(UseCompression);
/** Set the compression level. \sa ImageIOBase for details.
* Set to a negative number to use ImageIO's default compression level. */
itkSetMacro(CompressionLevel, int);
itkGetConstReferenceMacro(CompressionLevel, int);
/** By default the MetaDataDictionary is taken from the input image and
* passed to the ImageIO. In some cases, however, a user may prefer to
* introduce her/his own MetaDataDictionary. This is often the case of
* the ImageSeriesWriter. This flag defined whether the MetaDataDictionary
* to use will be the one from the input image or the one already set in
* the ImageIO object. */
itkSetMacro(UseInputMetaDataDictionary, bool);
itkGetConstReferenceMacro(UseInputMetaDataDictionary, bool);
itkBooleanMacro(UseInputMetaDataDictionary);
protected:
ImageFileWriter() = default;
~ImageFileWriter() override = default;
void
PrintSelf(std::ostream & os, Indent indent) const override;
/** Does the real work. */
void
GenerateData() override;
private:
std::string m_FileName{};
ImageIOBase::Pointer m_ImageIO{};
bool m_UserSpecifiedImageIO{ false };
ImageIORegion m_PasteIORegion{ TInputImage::ImageDimension };
unsigned int m_NumberOfStreamDivisions{ 1 };
bool m_UserSpecifiedIORegion{ false };
bool m_FactorySpecifiedImageIO{ false }; // did factory mechanism set the ImageIO?
bool m_UseCompression{ false };
int m_CompressionLevel{ -1 };
bool m_UseInputMetaDataDictionary{ true };
};
/** Convenience function for writing an image.
*
* The image parameter may be a either SmartPointer or a raw pointer and const or non-const.
* */
template <typename TImagePointer>
ITK_TEMPLATE_EXPORT void
WriteImage(TImagePointer && image, const std::string & filename, bool compress = false)
{
using NonReferenceImagePointer = std::remove_reference_t<TImagePointer>;
static_assert(std::is_pointer_v<NonReferenceImagePointer> || mpl::IsSmartPointer<NonReferenceImagePointer>::Value,
"WriteImage requires a raw pointer or SmartPointer.");
using ImageType = std::remove_const_t<std::remove_reference_t<decltype(*image)>>;
auto writer = ImageFileWriter<ImageType>::New();
writer->SetInput(image);
writer->SetFileName(filename);
writer->SetUseCompression(compress);
writer->Update();
}
} // end namespace itk
#ifndef ITK_MANUAL_INSTANTIATION
# include "itkImageFileWriter.hxx"
#endif
#if defined ITK_IMAGEIO_FACTORY_REGISTER_MANAGER || defined ITK_IO_FACTORY_REGISTER_MANAGER
# include "itkImageIOFactoryRegisterManager.h"
#endif
#endif // itkImageFileWriter_h
|