ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMetaImageIO.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 itkMetaImageIO_h
19 #define itkMetaImageIO_h
20 #include "ITKIOMetaExport.h"
21 
22 
23 #include <fstream>
24 #include "itkImageIOBase.h"
25 #include "itkSingletonMacro.h"
26 #include "metaObject.h"
27 #include "metaImage.h"
28 
29 namespace itk
30 {
41 class ITKIOMeta_EXPORT MetaImageIO:public ImageIOBase
42 {
43 public:
44  ITK_DISALLOW_COPY_AND_ASSIGN(MetaImageIO);
45 
47  using Self = MetaImageIO;
50 
52  itkNewMacro(Self);
53 
55  itkTypeMacro(MetaImageIO, Superclass);
56 
62  bool SupportsDimension(unsigned long) override
63  {
64  return true;
65  }
66 
67  /*-------- This part of the interfaces deals with reading data. ----- */
68 
71  bool CanReadFile(const char *) override;
72 
74  void ReadImageInformation() override;
75 
77  void Read(void *buffer) override;
78 
79  MetaImage * GetMetaImagePointer();
80 
81  /*-------- This part of the interfaces deals with writing data. ----- */
82 
85  bool CanWriteFile(const char *) override;
86 
88  void WriteImageInformation() override;
89 
92  void Write(const void *buffer) override;
93 
97  virtual void SetDataFileName(const char *filename);
98 
101  virtual void SetDoublePrecision(unsigned int precision)
102  {
103  m_MetaImage.SetDoublePrecision(precision);
104  }
105 
111  GenerateStreamableReadRegionFromRequestedRegion(const ImageIORegion & requested) const override;
112 
113  unsigned int
114  GetActualNumberOfSplitsForWriting(unsigned int numberOfRequestedSplits,
115  const ImageIORegion & pasteRegion,
116  const ImageIORegion & largestPossibleRegion) override;
117 
119  GetSplitRegionForWriting(unsigned int ithPiece,
120  unsigned int numberOfActualSplits,
121  const ImageIORegion & pasteRegion,
122  const ImageIORegion & largestPossibleRegion) override;
123 
127  bool CanStreamRead() override
128  {
129  if ( m_MetaImage.CompressedData() )
130  {
131  return false;
132  }
133  return true;
134  }
136 
142  bool CanStreamWrite() override
143  {
144  if ( this->GetUseCompression() )
145  {
146  return false;
147  }
148  return true;
149  }
151 
155  itkSetMacro(SubSamplingFactor, unsigned int);
156  itkGetConstMacro(SubSamplingFactor, unsigned int);
158 
169  static void SetDefaultDoublePrecision(unsigned int precision);
170  static unsigned int GetDefaultDoublePrecision();
172 
173 protected:
174  MetaImageIO();
175  ~MetaImageIO() override;
176  void PrintSelf(std::ostream & os, Indent indent) const override;
177 
178 private:
180  itkGetGlobalDeclarationMacro(unsigned int, DefaultDoublePrecision);
181 
182  MetaImage m_MetaImage;
183 
184  unsigned int m_SubSamplingFactor;
185 
186  static unsigned int * m_DefaultDoublePrecision;
187 };
188 } // end namespace itk
189 
190 #endif // itkMetaImageIO_h
Light weight base class for most itk classes.
virtual void SetDoublePrecision(unsigned int precision)
MetaImage m_MetaImage
An ImageIORegion represents a structured region of data.
bool CanStreamWrite() override
Abstract superclass defines image IO interface.
bool SupportsDimension(unsigned long) override
static unsigned int * m_DefaultDoublePrecision
#define itkGetGlobalDeclarationMacro(Type, VarName)
bool CanStreamRead() override
unsigned int m_SubSamplingFactor
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Base class for most ITK classes.
Definition: itkObject.h:60
Read MetaImage file format.