ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
130 using Self = Boundary;
134 itkNewMacro(Self);
135 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
183 flat_hash_t *
185 {
186 return this->GetFlatHash(idx.first, idx.second);
187 }
188 flat_hash_t *
189 GetFlatHash(unsigned int dimension, unsigned int highlow)
190 {
191 if (highlow == 0)
192 {
193 return &(m_FlatHashes[dimension].first);
194 }
195
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
218
224 void
225 SetValid(bool & l, const IndexType & idx)
226 {
227 this->SetValid(l, idx.first, idx.second);
228 }
229 void
230 SetValid(bool b, unsigned int dimension, unsigned int highlow)
231 {
232 if (highlow == 0)
233 {
234 m_Valid[dimension].first = b;
235 }
236 else
237 {
238 m_Valid[dimension].second = b;
239 }
240 this->Modified();
241 }
242
243 bool
244 GetValid(const IndexType & idx) const
245 {
246 return this->GetValid(idx.first, idx.second);
247 }
248 bool
249 GetValid(unsigned int dimension, unsigned int highlow) const
250 {
251 if (highlow == 0)
252 {
253 return m_Valid[dimension].first;
254 }
255
256 return m_Valid[dimension].second;
257 }
258
259
260protected:
262 ~Boundary() override = default;
263
264 void
265 PrintSelf(std::ostream & os, Indent indent) const override;
266
268 std::vector<std::pair<FacePointer, FacePointer>> m_Faces{};
269
272 std::vector<std::pair<flat_hash_t, flat_hash_t>> m_FlatHashes{};
273
276 std::vector<std::pair<bool, bool>> m_Valid{};
277};
278} // end namespace watershed
279} // end namespace itk
280
281#ifndef ITK_MANUAL_INSTANTIATION
282# include "itkWatershedBoundary.hxx"
283#endif
284
285#endif
Templated n-dimensional image class.
Definition itkImage.h:89
Index< VImageDimension > IndexType
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void Modified() const
Implements transparent reference counting.
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)
Image< IdentifierType, TDimension > ImageType
void SetFace(FacePointer f, unsigned int dimension, unsigned int highlow)
std::vector< std::pair< FacePointer, FacePointer > > m_Faces
SmartPointer< const Self > ConstPointer
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)
Image< face_pixel_t, TDimension > face_t
std::vector< std::pair< flat_hash_t, flat_hash_t > > m_FlatHashes
typename ImageType::IndexType ImageIndexType
std::vector< std::pair< bool, bool > > m_Valid
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)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType IdentifierType
Definition itkIntTypes.h:90