ITK  6.0.0
Insight Toolkit
itkWatershedBoundary.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 itkWatershedBoundary_h
19#define itkWatershedBoundary_h
20
21
22#include <list>
23#include <vector>
24#include "itkImage.h"
25#include "itkProcessObject.h"
26#include <unordered_map>
27
28namespace itk
29{
30namespace watershed
31{
53template <typename TScalar, unsigned int TDimension>
54class ITK_TEMPLATE_EXPORT Boundary : public DataObject
55{
56public:
57 ITK_DISALLOW_COPY_AND_MOVE(Boundary);
58
63 static constexpr unsigned int Dimension = TDimension;
64
69 using IndexType = std::pair<unsigned int, unsigned int>;
72 using ScalarType = TScalar;
73
76 {
77
92 short flow;
93
96 };
97
100 {
101
105 std::list<IdentifierType> offset_list;
106
110
114
117 };
118
122
124 using flat_hash_t = std::unordered_map<IdentifierType, flat_region_t>;
125 using FlatHashValueType = typename flat_hash_t::value_type;
126
129 using Self = Boundary;
133 itkNewMacro(Self);
134 itkOverrideGetNameOfClassMacro(Boundary);
138 using FacePointer = typename face_t::Pointer;
139
142 GetFace(const IndexType & idx)
143 {
144 return this->GetFace(idx.first, idx.second);
145 }
146
150 FacePointer
151 GetFace(unsigned int dimension, unsigned int highlow)
152 {
153 if (highlow == 0)
154 {
155 return m_Faces[dimension].first;
156 }
157
158 return m_Faces[dimension].second;
159 }
160
161 void
163 {
164 this->SetFace(f, idx.first, idx.second);
165 }
166
167 void
168 SetFace(FacePointer f, unsigned int dimension, unsigned int highlow)
169 {
170 if (highlow == 0)
171 {
172 m_Faces[dimension].first = f;
173 }
174 else
175 {
176 m_Faces[dimension].second = f;
177 }
178 this->Modified();
179 }
180
182 flat_hash_t *
184 {
185 return this->GetFlatHash(idx.first, idx.second);
186 }
187 flat_hash_t *
188 GetFlatHash(unsigned int dimension, unsigned int highlow)
189 {
190 if (highlow == 0)
191 {
192 return &(m_FlatHashes[dimension].first);
193 }
196 return &(m_FlatHashes[dimension].second);
197 }
198
199 void
201 {
202 this->SetFlatHash(l, idx.first, idx.second);
203 }
204 void
205 SetFlatHash(flat_hash_t & l, unsigned int dimension, unsigned int highlow)
206 {
207 if (highlow == 0)
208 {
209 m_FlatHashes[dimension].first = l;
210 }
211 else
212 {
213 m_FlatHashes[dimension].second = l;
214 }
215 this->Modified();
216 }
217
222 void
223 SetValid(bool & l, const IndexType & idx)
224 {
225 this->SetValid(l, idx.first, idx.second);
226 }
227 void
228 SetValid(bool b, unsigned int dimension, unsigned int highlow)
229 {
230 if (highlow == 0)
231 {
232 m_Valid[dimension].first = b;
233 }
234 else
235 {
236 m_Valid[dimension].second = b;
237 }
238 this->Modified();
239 }
242 bool
243 GetValid(const IndexType & idx) const
244 {
245 return this->GetValid(idx.first, idx.second);
246 }
247 bool
248 GetValid(unsigned int dimension, unsigned int highlow) const
249 {
250 if (highlow == 0)
251 {
252 return m_Valid[dimension].first;
253 }
254
255 return m_Valid[dimension].second;
256 }
257
258protected:
260 ~Boundary() override = default;
261
262 void
263 PrintSelf(std::ostream & os, Indent indent) const override;
264
266 std::vector<std::pair<FacePointer, FacePointer>> m_Faces{};
267
270 std::vector<std::pair<flat_hash_t, flat_hash_t>> m_FlatHashes{};
271
274 std::vector<std::pair<bool, bool>> m_Valid{};
275};
276} // end namespace watershed
277} // end namespace itk
278
279#ifndef ITK_MANUAL_INSTANTIATION
280# include "itkWatershedBoundary.hxx"
281#endif
282
283#endif
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for most ITK classes.
Definition: itkObject.h:62
flat_hash_t * GetFlatHash(const IndexType &idx)
flat_hash_t * GetFlatHash(unsigned int dimension, unsigned int highlow)
bool GetValid(unsigned int dimension, unsigned int highlow) const
void SetValid(bool b, unsigned int dimension, unsigned int highlow)
void SetFace(FacePointer f, unsigned int dimension, unsigned int highlow)
void SetValid(bool &l, const IndexType &idx)
~Boundary() override=default
void SetFace(FacePointer f, const IndexType &idx)
typename face_t::Pointer FacePointer
bool GetValid(const IndexType &idx) const
std::pair< unsigned int, unsigned int > IndexType
std::unordered_map< IdentifierType, flat_region_t > flat_hash_t
void PrintSelf(std::ostream &os, Indent indent) const override
typename flat_hash_t::value_type FlatHashValueType
void SetFlatHash(flat_hash_t &l, const IndexType &idx)
typename ImageType::IndexType ImageIndexType
FacePointer GetFace(const IndexType &idx)
FacePointer GetFace(unsigned int dimension, unsigned int highlow)
void SetFlatHash(flat_hash_t &l, unsigned int dimension, unsigned int highlow)
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
class ITK_FORWARD_EXPORT DataObject
Definition: itkDataObject.h:42
SizeValueType IdentifierType
Definition: itkIntTypes.h:90