ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkImageSeriesWriter.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * https://www.apache.org/licenses/LICENSE-2.0.txt
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 *
17 *=========================================================================*/
18#ifndef itkImageSeriesWriter_h
19#define itkImageSeriesWriter_h
20#include "ITKIOImageBaseExport.h"
21
22#include "itkImageRegion.h"
23#include "itkImageFileWriter.h"
24#include <vector>
25#include <string>
26
27namespace itk
28{
33class ITKIOImageBase_EXPORT ImageSeriesWriterException : public ExceptionObject
34{
35public:
37 ~ImageSeriesWriterException() noexcept override;
38
40 itkOverrideGetNameOfClassMacro(ImageSeriesWriterException);
41
43 ImageSeriesWriterException(std::string file, unsigned int line, std::string message = "Error in IO")
44 : ExceptionObject(std::move(file), line, std::move(message))
45 {}
46};
47
74template <typename TInputImage, typename TOutputImage>
75class ITK_TEMPLATE_EXPORT ImageSeriesWriter : public ProcessObject
76{
77public:
78 ITK_DISALLOW_COPY_AND_MOVE(ImageSeriesWriter);
79
85
87 itkNewMacro(Self);
88
90 itkOverrideGetNameOfClassMacro(ImageSeriesWriter);
91
93 using InputImageType = TInputImage;
94 using InputImageRegionType = typename InputImageType::RegionType;
95 using OutputImageType = TOutputImage;
96 using OutputImageRegionType = typename OutputImageType::RegionType;
98 using FileNamesContainer = std::vector<std::string>;
99
103 using DictionaryArrayType = std::vector<DictionaryRawPointer>;
105
108 void
109 SetInput(const InputImageType * input);
110
111 const InputImageType *
113
114 const InputImageType *
115 GetInput(unsigned int idx);
116
124 itkSetObjectMacro(ImageIO, ImageIOBase);
125 itkGetModifiableObjectMacro(ImageIO, ImageIOBase);
127
132 virtual void
134
137 void
138 Update() override
139 {
140 this->Write();
141 }
142
146 itkSetMacro(StartIndex, SizeValueType);
147 itkGetConstMacro(StartIndex, SizeValueType);
149
153 itkSetMacro(IncrementIndex, SizeValueType);
154 itkGetConstMacro(IncrementIndex, SizeValueType);
156
162 itkSetStringMacro(SeriesFormat);
163 itkGetStringMacro(SeriesFormat);
165
169 void
171 {
172 if (m_FileNames != name)
173 {
174 m_FileNames = name;
175 this->Modified();
176 }
177 }
178
179 const FileNamesContainer &
181 {
182 return m_FileNames;
183 }
184
185
188 void
189 SetFileName(const std::string & name)
190 {
191 m_FileNames.clear();
192 m_FileNames.push_back(name);
193 this->Modified();
194 }
195
198 void
199 AddFileName(const std::string & name)
200 {
201 m_FileNames.push_back(name);
202 this->Modified();
203 }
204
207 itkSetMacro(MetaDataDictionaryArray, DictionaryArrayRawPointer);
208
211 itkSetMacro(UseCompression, bool);
212 itkGetConstReferenceMacro(UseCompression, bool);
213 itkBooleanMacro(UseCompression);
215protected:
217 ~ImageSeriesWriter() override = default;
218 void
219 PrintSelf(std::ostream & os, Indent indent) const override;
220
222 void
223 GenerateData() override;
224
227 void
229
231
232 // track whether the ImageIO is user specified
234
235private:
238
244 std::string m_SeriesFormat{};
247
249
252
253 // These two methods provide now a common implementation for the
254 // GenerateNumericFileNamesAndWrite() and avoid the duplication of code that
255 // was leaving one of the code branches out of date.
256 void
258
259 void
261};
262} // end namespace itk
263
264#ifndef ITK_MANUAL_INSTANTIATION
265# include "itkImageSeriesWriter.hxx"
266#endif
267
268#endif // itkImageSeriesWriter_h
ExceptionObject() noexcept=default
Writes image data to a single file.
Abstract superclass defines image IO interface.
SmartPointer< Self > Pointer
ImageSeriesWriterException(std::string file, unsigned int line, std::string message="Error in IO")
~ImageSeriesWriterException() noexcept override
void GenerateNumericFileNamesAndWrite()
std::vector< DictionaryRawPointer > DictionaryArrayType
FileNamesContainer m_FileNames
void SetFileName(const std::string &name)
const InputImageType * GetInput(unsigned int idx)
~ImageSeriesWriter() override=default
void GenerateData() override
const DictionaryArrayType * DictionaryArrayRawPointer
void SetInput(const InputImageType *input)
const InputImageType * GetInput()
SmartPointer< Self > Pointer
void AddFileName(const std::string &name)
std::vector< std::string > FileNamesContainer
const FileNamesContainer & GetFileNames() const
virtual void Write()
SmartPointer< const Self > ConstPointer
typename InputImageType::RegionType InputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
void SetFileNames(const FileNamesContainer &name)
void PrintSelf(std::ostream &os, Indent indent) const override
MetaDataDictionary * DictionaryRawPointer
ImageIOBase::Pointer m_ImageIO
DictionaryArrayRawPointer m_MetaDataDictionaryArray
MetaDataDictionary DictionaryType
ImageFileWriter< TOutputImage > WriterType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Provides a mechanism for storing a collection of arbitrary data types.
virtual void Modified() const
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86
STL namespace.