ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkAttributeMorphologyBaseImageFilter.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 itkAttributeMorphologyBaseImageFilter_h
19#define itkAttributeMorphologyBaseImageFilter_h
20
22#include <vector>
23#include <memory> // For unique_ptr.
24
25namespace itk
26{
50
51template <typename TInputImage, typename TOutputImage, typename TAttribute, typename TFunction>
52class ITK_TEMPLATE_EXPORT AttributeMorphologyBaseImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
53{
54public:
60
64 using typename Superclass::InputImagePointer;
65
70 using OutputPixelType = typename TOutputImage::PixelType;
71 using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
72 using InputPixelType = typename TInputImage::PixelType;
73 using InputInternalPixelType = typename TInputImage::InternalPixelType;
74 using IndexType = typename TInputImage::IndexType;
75 using OffsetType = typename TInputImage::OffsetType;
76 using SizeType = typename TInputImage::SizeType;
77
78 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
79
83 using InputImageType = TInputImage;
84 using OutputImageType = TOutputImage;
85 // using IndexType = typename TInputImage::IndexType;
86 // using SizeType = typename TInputImage::SizeType;
87 using RegionType = typename TOutputImage::RegionType;
88 using ListType = std::list<IndexType>;
89 using AttributeType = TAttribute;
90
96
100 itkOverrideGetNameOfClassMacro(AttributeMorphologyBaseImageFilter);
101
105 itkNewMacro(Self);
106
113 itkSetMacro(FullyConnected, bool);
114 itkGetConstReferenceMacro(FullyConnected, bool);
115 itkBooleanMacro(FullyConnected);
117
123 itkSetMacro(Lambda, AttributeType);
124 itkGetConstMacro(Lambda, AttributeType);
126
127protected:
134
137 void
138 PrintSelf(std::ostream & os, Indent indent) const override;
139
143 void
144 GenerateData() override;
145
149 void
151
156 void
157 EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
158
160
161private:
164
165 // some constants used several times in the code
166 static constexpr OffsetValueType INACTIVE = -1;
167 static constexpr OffsetValueType ACTIVE = -2;
168
169 // Just used for area/volume openings at the moment
170 std::unique_ptr<AttributeType[]> m_AuxData;
171
172 using OffsetVecType = std::vector<OffsetType>;
173 // offset in the linear array.
174 using OffsetDirectVecType = std::vector<OffsetValueType>;
175
176 void
178
179 // m_SortPixels contains offsets into the raw image
180 // it is sorted with a stable sort by grey level as the
181 // first step in the algorithm. The sorting step avoids
182 // the need to explicitly locate regional extrema.
183 std::unique_ptr<OffsetValueType[]> m_SortPixels;
184 std::unique_ptr<OffsetValueType[]> m_Parent;
185
186 // This is a bit ugly, but I can't see an easy way around
187 std::unique_ptr<InputPixelType[]> m_Raw;
188
190 {
191 public:
192 TFunction m_TFunction;
193 // buf contains the raw data, which is what
194 // we want to sort by. i.e. the first value in
195 // the sorted buffer will be the location of the
196 // largest or smallest pixel.
198 bool
199 operator()(const OffsetValueType & l, const OffsetValueType & r) const
200 {
201 return (m_TFunction(buf[l], buf[r]));
202 }
203 };
204
205 CompareOffsetType m_CompareOffset{};
206 // version from PAMI. Note - using the AuxData array rather than the
207 // parent array to store area
208 void
214
217 {
218 if (m_Parent[x] >= 0)
219 {
220 m_Parent[x] = FindRoot(m_Parent[x]);
221 return (m_Parent[x]);
222 }
223 else
224 {
225 return (x);
226 }
227 }
228
229 bool
231 {
232 return ((m_Raw[x] == m_Raw[y]) || (m_AuxData[x] < m_Lambda));
233 }
234
235 void
237 {
239
240 if (r != p)
241 {
242 if (Criterion(r, p))
243 {
244 m_AuxData[p] += m_AuxData[r];
245 m_Parent[r] = p;
246 }
247 else
248 {
249 m_AuxData[p] = m_Lambda;
250 }
251 }
252 }
253};
254} // end namespace itk
255
256#ifndef ITK_MANUAL_INSTANTIATION
257# include "itkAttributeMorphologyBaseImageFilter.hxx"
258#endif
259
260#endif
bool operator()(const OffsetValueType &l, const OffsetValueType &r) const
void EnlargeOutputRequestedRegion(DataObject *output) override
void PrintSelf(std::ostream &os, Indent indent) const override
void SetupOffsetVec(OffsetDirectVecType &PosOffsets, OffsetVecType &Offsets)
bool Criterion(OffsetValueType x, OffsetValueType y)
~AttributeMorphologyBaseImageFilter() override=default
Base class for all data objects in ITK.
typename InputImageType::Pointer InputImagePointer
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long OffsetValueType
Definition itkIntTypes.h:97