ITK  6.0.0
Insight Toolkit
itkImportImageContainer.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 itkImportImageContainer_h
19 #define itkImportImageContainer_h
20 
21 #include "itkObject.h"
22 #include "itkObjectFactory.h"
23 #include <utility>
24 
25 namespace itk
26 {
44 template <typename TElementIdentifier, typename TElement>
45 class ITK_TEMPLATE_EXPORT ImportImageContainer : public Object
46 {
47 public:
48  ITK_DISALLOW_COPY_AND_MOVE(ImportImageContainer);
49 
52  using Superclass = Object;
55 
57  using ElementIdentifier = TElementIdentifier;
58  using Element = TElement;
59 
61  itkNewMacro(Self);
62 
64  itkOverrideGetNameOfClassMacro(ImportImageContainer);
65 
67  TElement *
69  {
70  return m_ImportPointer;
71  }
72 
79  void
80  SetImportPointer(TElement * ptr, TElementIdentifier num, bool LetContainerManageMemory = false);
81 
83  TElement &
85  {
86  return m_ImportPointer[id];
87  }
88 
90  const TElement &
91  operator[](const ElementIdentifier id) const
92  {
93  return m_ImportPointer[id];
94  }
95 
98  TElement *
100  {
101  return m_ImportPointer;
102  }
103 
105  ElementIdentifier
106  Capacity() const
107  {
108  return m_Capacity;
109  }
110 
112  ElementIdentifier
113  Size() const
114  {
115  return m_Size;
116  }
117 
133  void
134  Reserve(ElementIdentifier size, const bool UseValueInitialization = false);
135 
142  void
143  Squeeze();
144 
146  void
147  Initialize();
148 
158  itkSetMacro(ContainerManageMemory, bool);
159  itkGetConstMacro(ContainerManageMemory, bool);
160  itkBooleanMacro(ContainerManageMemory);
163 protected:
164  ImportImageContainer() = default;
165  ~ImportImageContainer() override;
166 
170  void
171  PrintSelf(std::ostream & os, Indent indent) const override;
172 
177  virtual TElement *
178  AllocateElements(ElementIdentifier size, bool UseValueInitialization = false) const;
179 
180  virtual void
181  DeallocateManagedMemory();
182 
183  /* Set the m_Size member that represents the number of elements
184  * currently stored in the container. Use this function with great
185  * care since it only changes the m_Size member and not the actual size
186  * of the import pointer m_ImportPointer. It should typically
187  * be used only to override AllocateElements and
188  * DeallocateManagedMemory. */
189  itkSetMacro(Size, TElementIdentifier);
190 
191  /* Set the m_Capacity member that represents the capacity of
192  * the current container. Use this function with great care
193  * since it only changes the m_Capacity member and not the actual
194  * capacity of the import pointer m_ImportPointer. It should typically
195  * be used only to override AllocateElements and
196  * DeallocateManagedMemory. */
197  itkSetMacro(Capacity, TElementIdentifier);
198 
199  /* Set the m_ImportPointer member. Use this function with great care
200  * since it only changes the m_ImportPointer member but not the m_Size
201  * and m_Capacity members. It should typically be used only to override
202  * AllocateElements and DeallocateManagedMemory. */
203  void
204  SetImportPointer(TElement * ptr)
205  {
206  m_ImportPointer = ptr;
207  }
208 
209 private:
210  TElement * m_ImportPointer{};
211  TElementIdentifier m_Size{};
212  TElementIdentifier m_Capacity{};
213  bool m_ContainerManageMemory{ true };
214 };
215 } // end namespace itk
216 
217 #ifndef ITK_MANUAL_INSTANTIATION
218 # include "itkImportImageContainer.hxx"
219 #endif
220 
221 #endif
itkObjectFactory.h
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itk::ImportImageContainer::Capacity
ElementIdentifier Capacity() const
Definition: itkImportImageContainer.h:106
itk::ImportImageContainer::SetImportPointer
void SetImportPointer(TElement *ptr)
Definition: itkImportImageContainer.h:204
itk::ImportImageContainer::ElementIdentifier
TElementIdentifier ElementIdentifier
Definition: itkImportImageContainer.h:57
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImportImageContainer::operator[]
const TElement & operator[](const ElementIdentifier id) const
Definition: itkImportImageContainer.h:91
itk::ImportImageContainer::GetBufferPointer
TElement * GetBufferPointer()
Definition: itkImportImageContainer.h:99
itk::ImportImageContainer
Defines an itk::Image front-end to a standard C-array.
Definition: itkImportImageContainer.h:45
itkObject.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ImportImageContainer::GetImportPointer
TElement * GetImportPointer()
Definition: itkImportImageContainer.h:68
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::ImportImageContainer::Size
ElementIdentifier Size() const
Definition: itkImportImageContainer.h:113
itk::ImportImageContainer::operator[]
TElement & operator[](const ElementIdentifier id)
Definition: itkImportImageContainer.h:84
itk::ImportImageContainer::Element
TElement Element
Definition: itkImportImageContainer.h:58