ITK  6.0.0
Insight Toolkit
itkPointSetBase.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/*=========================================================================
19 *
20 * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21 *
22 * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23 *
24 * For complete copyright, license and disclaimer of warranty information
25 * please refer to the NOTICE file at the top of the ITK source tree.
26 *
27 *=========================================================================*/
28#ifndef itkPointSetBase_h
29#define itkPointSetBase_h
30
31#include "itkDataObject.h"
32#include "itkVectorContainer.h"
33#include <vector>
34
35namespace itk
36{
37
55template <typename TPointsContainer>
56class ITK_TEMPLATE_EXPORT PointSetBase : public DataObject
57{
58public:
59 ITK_DISALLOW_COPY_AND_MOVE(PointSetBase);
60
66
68 itkOverrideGetNameOfClassMacro(PointSetBase);
69
71
73 using PointType = typename TPointsContainer::Element;
74 using CoordRepType = typename PointType::CoordRepType;
75 using PointIdentifier = typename TPointsContainer::ElementIdentifier;
76 using PointsContainer = TPointsContainer;
77
81
83 static constexpr unsigned int PointDimension = PointType::PointDimension;
84
88
90 using PointsContainerConstIterator = typename PointsContainer::ConstIterator;
91 using PointsContainerIterator = typename PointsContainer::Iterator;
92
94 using RegionType = long;
95
98 itkGetConstMacro(MaximumNumberOfRegions, RegionType);
99
100protected:
103 PointsContainerPointer m_PointsContainer{};
104
105public:
109 void
110 PassStructure(Self * inputPointSet);
111
115 void
116 Initialize() override;
117
121
123 void
125
128 void
130
132 void
133 SetPointsByCoordinates(const std::vector<CoordRepType> & coordinates);
134
138
140 const PointsContainer *
141 GetPoints() const;
142
147
154 bool
156
159
161 void
163
164 void
166
167 void
168 CopyInformation(const DataObject * data) override;
169
170 bool
172
173 bool
175
180 void
181 SetRequestedRegion(const DataObject * data) override;
182
184 virtual void
186
187 itkGetConstMacro(RequestedRegion, RegionType);
188
190 virtual void
192
193 itkGetConstMacro(BufferedRegion, RegionType);
194
195protected:
197 PointSetBase() = default;
198
200 ~PointSetBase() override = 0;
201
202 void
203 PrintSelf(std::ostream & os, Indent indent) const override;
204
206 InternalClone() const override;
207
208 // If the RegionType is ITK_UNSTRUCTURED_REGION, then the following
209 // variables represent the maximum number of region that the data
210 // object can be broken into, which region out of how many is
211 // currently in the buffered region, and the number of regions and
212 // the specific region requested for the update. Data objects that
213 // do not support any division of the data can simply leave the
214 // MaximumNumberOfRegions as 1. The RequestedNumberOfRegions and
215 // RequestedRegion are used to define the currently requested
216 // region. The LargestPossibleRegion is always requested region = 0
217 // and number of regions = 1;
218 RegionType m_MaximumNumberOfRegions{ 1 };
219 RegionType m_NumberOfRegions{ 1 };
220 RegionType m_RequestedNumberOfRegions{};
221 RegionType m_BufferedRegion{ -1 };
222 RegionType m_RequestedRegion{ -1 };
223};
224} // end namespace itk
225
226#ifndef ITK_MANUAL_INSTANTIATION
227# include "itkPointSetBase.hxx"
228#endif
229
230#endif
Base class for all data objects in ITK.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for most ITK classes.
Definition: itkObject.h:62
A superclass of PointSet supports point (geometric coordinate and attribute) definition.
virtual void SetBufferedRegion(const RegionType &region)
void SetRequestedRegionToLargestPossibleRegion() override
void Initialize() override
void SetPoints(PointsVectorContainer *)
void CopyInformation(const DataObject *data) override
void SetPointsByCoordinates(const std::vector< CoordRepType > &coordinates)
typename PointsVectorContainer::Pointer PointsVectorContainerPointer
typename PointsContainer::Iterator PointsContainerIterator
typename PointsContainer::ConstIterator PointsContainerConstIterator
void UpdateOutputInformation() override
typename itk::VectorContainer< PointIdentifier, CoordRepType > PointsVectorContainer
void SetRequestedRegion(const DataObject *data) override
typename TPointsContainer::ElementIdentifier PointIdentifier
typename PointsContainer::Pointer PointsContainerPointer
LightObject::Pointer InternalClone() const override
typename TPointsContainer::Element PointType
TPointsContainer PointsContainer
bool GetPoint(PointIdentifier, PointType *) const
virtual void SetRequestedRegion(const RegionType &region)
const PointsContainer * GetPoints() const
PointSetBase()=default
bool VerifyRequestedRegion() override
PointsContainer * GetPoints()
void PrintSelf(std::ostream &os, Indent indent) const override
PointIdentifier GetNumberOfPoints() const
void SetPoint(PointIdentifier, PointType)
void SetPoints(PointsContainer *)
~PointSetBase() override=0
bool RequestedRegionIsOutsideOfTheBufferedRegion() override
void PassStructure(Self *inputPointSet)
typename PointType::CoordRepType CoordRepType
PointType GetPoint(PointIdentifier) const
typename PointsContainer::ConstPointer PointsContainerConstPointer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
class ITK_FORWARD_EXPORT DataObject
Definition: itkDataObject.h:42