ITK  4.1.0
Insight Segmentation and Registration Toolkit
itkImageSeriesWriter.h
Go to the documentation of this file.
00001 /*=========================================================================
00002  *
00003  *  Copyright Insight Software Consortium
00004  *
00005  *  Licensed under the Apache License, Version 2.0 (the "License");
00006  *  you may not use this file except in compliance with the License.
00007  *  You may obtain a copy of the License at
00008  *
00009  *         http://www.apache.org/licenses/LICENSE-2.0.txt
00010  *
00011  *  Unless required by applicable law or agreed to in writing, software
00012  *  distributed under the License is distributed on an "AS IS" BASIS,
00013  *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
00014  *  See the License for the specific language governing permissions and
00015  *  limitations under the License.
00016  *
00017  *=========================================================================*/
00018 #ifndef __itkImageSeriesWriter_h
00019 #define __itkImageSeriesWriter_h
00020 
00021 #include "itkImageRegion.h"
00022 #include "itkImageFileWriter.h"
00023 #include <vector>
00024 #include <string>
00025 
00026 namespace itk
00027 {
00032 class ITK_ABI_EXPORT ImageSeriesWriterException:public ExceptionObject
00033 {
00034 public:
00036   itkTypeMacro(ImageSeriesWriterException, ExceptionObject);
00037 
00039   ImageSeriesWriterException(char *file, unsigned int line,
00040                              const char *message = "Error in IO"):
00041     ExceptionObject(file, line)
00042   {
00043     SetDescription(message);
00044   }
00045 
00047   ImageSeriesWriterException(const std::string & file, unsigned int line,
00048                              const char *message = "Error in IO"):
00049     ExceptionObject(file, line)
00050   {
00051     SetDescription(message);
00052   }
00053 };
00055 
00078 template< class TInputImage, class TOutputImage >
00079 class ITK_EXPORT ImageSeriesWriter:public ProcessObject
00080 {
00081 public:
00083   typedef ImageSeriesWriter          Self;
00084   typedef ProcessObject              Superclass;
00085   typedef SmartPointer< Self >       Pointer;
00086   typedef SmartPointer< const Self > ConstPointer;
00087 
00089   itkNewMacro(Self);
00090 
00092   itkTypeMacro(ImageSeriesWriter, ProcessObject);
00093 
00095   typedef TInputImage                          InputImageType;
00096   typedef typename InputImageType::RegionType  InputImageRegionType;
00097   typedef TOutputImage                         OutputImageType;
00098   typedef typename OutputImageType::RegionType OutputImageRegionType;
00099   typedef ImageFileWriter< TOutputImage >      WriterType;
00100   typedef std::vector< std::string >           FileNamesContainer;
00101 
00103   typedef MetaDataDictionary                  DictionaryType;
00104   typedef MetaDataDictionary *                DictionaryRawPointer;
00105   typedef std::vector< DictionaryRawPointer > DictionaryArrayType;
00106   typedef const DictionaryArrayType *         DictionaryArrayRawPointer;
00107 
00109   using Superclass::SetInput;
00110   void SetInput(const InputImageType *input);
00111 
00112   const InputImageType * GetInput(void);
00113 
00114   const InputImageType * GetInput(unsigned int idx);
00115 
00122   itkSetObjectMacro(ImageIO, ImageIOBase);
00123   itkGetObjectMacro(ImageIO, ImageIOBase);
00125 
00130   virtual void Write(void);
00131 
00134   virtual void Update()
00135   {
00136     this->Write();
00137   }
00138 
00141   itkSetMacro(StartIndex, SizeValueType);
00142   itkGetConstMacro(StartIndex, SizeValueType);
00144 
00147   itkSetMacro(IncrementIndex, SizeValueType);
00148   itkGetConstMacro(IncrementIndex, SizeValueType);
00150 
00155   itkSetStringMacro(SeriesFormat);
00156   itkGetStringMacro(SeriesFormat);
00158 
00161   void SetFileNames(const FileNamesContainer & name)
00162   {
00163     if ( m_FileNames != name )
00164       {
00165       m_FileNames = name;
00166       this->Modified();
00167       }
00168   }
00170 
00171   const FileNamesContainer & GetFileNames() const
00172   {
00173     return m_FileNames;
00174   }
00175 
00178   void SetFileName(std::string const & name)
00179   {
00180     m_FileNames.clear();
00181     m_FileNames.push_back(name);
00182     this->Modified();
00183   }
00185 
00188   void AddFileName(std::string const & name)
00189   {
00190     m_FileNames.push_back(name);
00191     this->Modified();
00192   }
00194 
00197   itkSetMacro(MetaDataDictionaryArray, DictionaryArrayRawPointer);
00198 
00200   itkSetMacro(UseCompression, bool);
00201   itkGetConstReferenceMacro(UseCompression, bool);
00202   itkBooleanMacro(UseCompression);
00203 protected:
00204   ImageSeriesWriter();
00205   ~ImageSeriesWriter();
00206   void PrintSelf(std::ostream & os, Indent indent) const;
00208 
00210   void GenerateData(void);
00211 
00214   void GenerateNumericFileNamesAndWrite(void);
00215 
00216   ImageIOBase::Pointer m_ImageIO;
00217 
00218   //track whether the ImageIO is user specified
00219   bool m_UserSpecifiedImageIO;
00220 private:
00221   ImageSeriesWriter(const Self &); //purposely not implemented
00222   void operator=(const Self &);    //purposely not implemented
00223 
00225   FileNamesContainer m_FileNames;
00226 
00232   std::string   m_SeriesFormat;
00233   SizeValueType m_StartIndex;
00234   SizeValueType m_IncrementIndex;
00235 
00236   bool m_UseCompression;
00237 
00239   DictionaryArrayRawPointer m_MetaDataDictionaryArray;
00240 
00241   // These two methods provide now a common implementation for the
00242   // GenerateNumericFileNamesAndWrite() and avoid the duplication of code that
00243   // was leaving one of the code branches out of date.
00244   void GenerateNumericFileNames(void);
00245 
00246   void WriteFiles();
00247 };
00248 } // end namespace itk
00249 
00250 #ifndef ITK_MANUAL_INSTANTIATION
00251 #include "itkImageSeriesWriter.hxx"
00252 #endif
00253 
00254 #endif // __itkImageSeriesWriter_h
00255