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 & operator[](const ElementIdentifier id) { return m_ImportPointer[id]; }
84 
86  const TElement & operator[](const ElementIdentifier id) const { return m_ImportPointer[id]; }
87 
90  TElement *
92  {
93  return m_ImportPointer;
94  }
95 
97  ElementIdentifier
98  Capacity() const
99  {
100  return m_Capacity;
101  }
102 
104  ElementIdentifier
105  Size() const
106  {
107  return m_Size;
108  }
109 
125  void
126  Reserve(ElementIdentifier size, const bool UseValueInitialization = false);
127 
134  void
135  Squeeze();
136 
138  void
139  Initialize();
140 
150  itkSetMacro(ContainerManageMemory, bool);
151  itkGetConstMacro(ContainerManageMemory, bool);
152  itkBooleanMacro(ContainerManageMemory);
155 protected:
156  ImportImageContainer() = default;
157  ~ImportImageContainer() override;
158 
162  void
163  PrintSelf(std::ostream & os, Indent indent) const override;
164 
169  virtual TElement *
170  AllocateElements(ElementIdentifier size, bool UseValueInitialization = false) const;
171 
172  virtual void
173  DeallocateManagedMemory();
174 
175  /* Set the m_Size member that represents the number of elements
176  * currently stored in the container. Use this function with great
177  * care since it only changes the m_Size member and not the actual size
178  * of the import pointer m_ImportPointer. It should typically
179  * be used only to override AllocateElements and
180  * DeallocateManagedMemory. */
181  itkSetMacro(Size, TElementIdentifier);
182 
183  /* Set the m_Capacity member that represents the capacity of
184  * the current container. Use this function with great care
185  * since it only changes the m_Capacity member and not the actual
186  * capacity of the import pointer m_ImportPointer. It should typically
187  * be used only to override AllocateElements and
188  * DeallocateManagedMemory. */
189  itkSetMacro(Capacity, TElementIdentifier);
190 
191  /* Set the m_ImportPointer member. Use this function with great care
192  * since it only changes the m_ImportPointer member but not the m_Size
193  * and m_Capacity members. It should typically be used only to override
194  * AllocateElements and DeallocateManagedMemory. */
195  void
196  SetImportPointer(TElement * ptr)
197  {
198  m_ImportPointer = ptr;
199  }
200 
201 private:
202  TElement * m_ImportPointer{};
203  TElementIdentifier m_Size{};
204  TElementIdentifier m_Capacity{};
205  bool m_ContainerManageMemory{ true };
206 };
207 } // end namespace itk
208 
209 #ifndef ITK_MANUAL_INSTANTIATION
210 # include "itkImportImageContainer.hxx"
211 #endif
212 
213 #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:98
itk::ImportImageContainer::SetImportPointer
void SetImportPointer(TElement *ptr)
Definition: itkImportImageContainer.h:196
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:86
itk::ImportImageContainer::GetBufferPointer
TElement * GetBufferPointer()
Definition: itkImportImageContainer.h:91
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:105
itk::ImportImageContainer::operator[]
TElement & operator[](const ElementIdentifier id)
Definition: itkImportImageContainer.h:83
itk::ImportImageContainer::Element
TElement Element
Definition: itkImportImageContainer.h:58