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 else
158 {
159 return m_Faces[dimension].second;
160 }
161 }
164 void
166 {
167 this->SetFace(f, idx.first, idx.second);
168 }
169
170 void
171 SetFace(FacePointer f, unsigned int dimension, unsigned int highlow)
172 {
173 if (highlow == 0)
174 {
175 m_Faces[dimension].first = f;
176 }
177 else
178 {
179 m_Faces[dimension].second = f;
180 }
181 this->Modified();
182 }
183
185 flat_hash_t *
187 {
188 return this->GetFlatHash(idx.first, idx.second);
189 }
190 flat_hash_t *
191 GetFlatHash(unsigned int dimension, unsigned int highlow)
192 {
193 if (highlow == 0)
194 {
195 return &(m_FlatHashes[dimension].first);
196 }
197 else
198 {
199 return &(m_FlatHashes[dimension].second);
200 }
201 }
204 void
206 {
207 this->SetFlatHash(l, idx.first, idx.second);
208 }
209 void
210 SetFlatHash(flat_hash_t & l, unsigned int dimension, unsigned int highlow)
211 {
212 if (highlow == 0)
213 {
214 m_FlatHashes[dimension].first = l;
215 }
216 else
217 {
218 m_FlatHashes[dimension].second = l;
219 }
220 this->Modified();
221 }
222
227 void
228 SetValid(bool & l, const IndexType & idx)
229 {
230 this->SetValid(l, idx.first, idx.second);
231 }
232 void
233 SetValid(bool b, unsigned int dimension, unsigned int highlow)
234 {
235 if (highlow == 0)
236 {
237 m_Valid[dimension].first = b;
238 }
239 else
240 {
241 m_Valid[dimension].second = b;
242 }
243 this->Modified();
244 }
247 bool
248 GetValid(const IndexType & idx) const
249 {
250 return this->GetValid(idx.first, idx.second);
251 }
252 bool
253 GetValid(unsigned int dimension, unsigned int highlow) const
254 {
255 if (highlow == 0)
256 {
257 return m_Valid[dimension].first;
258 }
259 else
260 {
261 return m_Valid[dimension].second;
262 }
263 }
264
265protected:
267 ~Boundary() override = default;
268
269 void
270 PrintSelf(std::ostream & os, Indent indent) const override;
271
273 std::vector<std::pair<FacePointer, FacePointer>> m_Faces{};
274
277 std::vector<std::pair<flat_hash_t, flat_hash_t>> m_FlatHashes{};
278
281 std::vector<std::pair<bool, bool>> m_Valid{};
282};
283} // end namespace watershed
284} // end namespace itk
285
286#ifndef ITK_MANUAL_INSTANTIATION
287# include "itkWatershedBoundary.hxx"
288#endif
289
290#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