ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkMorphologicalContourInterpolator.h
Go to the documentation of this file.
1 /*=========================================================================
2 *
3 * Copyright Insight Software Consortium
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 * http://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 itkMorphologicalContourInterpolator_h
19 #define itkMorphologicalContourInterpolator_h
20 
23 #include "itkExtractImageFilter.h"
24 #include "itkImageToImageFilter.h"
25 #include "itksys/hash_map.hxx"
26 
27 namespace itk
28 {
64 template< typename TImage >
66  public ImageToImageFilter< TImage, TImage >
67 {
68  template< typename T >
70 
71 public:
76  typedef Image< typename TImage::PixelType, TImage::ImageDimension - 1 > SliceType;
77 
79  itkNewMacro( Self );
80 
83 
85  itkSetMacro( Label, typename TImage::PixelType );
86 
88  itkGetMacro( Label, typename TImage::PixelType );
89 
91  itkGetConstMacro( Label, typename TImage::PixelType );
92 
94  itkSetMacro( Axis, int );
95 
97  itkGetMacro( Axis, int );
98 
100  itkGetConstMacro( Axis, int );
101 
104  itkSetMacro( HeuristicAlignment, bool );
105 
108  itkGetMacro( HeuristicAlignment, bool );
109 
112  itkGetConstMacro( HeuristicAlignment, bool );
113 
117  itkSetMacro( UseDistanceTransform, bool );
118 
122  itkGetMacro( UseDistanceTransform, bool );
123 
127  itkGetConstMacro( UseDistanceTransform, bool );
128 
131  itkSetMacro( UseCustomSlicePositions, bool );
132 
135  itkGetMacro( UseCustomSlicePositions, bool );
136 
139  itkGetConstMacro( UseCustomSlicePositions, bool );
140 
142  void
144  {
145  if ( useBall != m_UseBallStructuringElement )
146  {
147  m_UseBallStructuringElement = useBall;
148  m_ConnectedComponents->SetFullyConnected( useBall );
149  this->Modified();
150  }
151  }
153 
155  itkGetMacro( UseBallStructuringElement, bool );
156 
158  itkGetConstMacro( UseBallStructuringElement, bool );
159 
164  void
166 
168  typedef std::set< typename TImage::IndexValueType > SliceSetType;
169 
171  void
173  {
174  m_LabeledSlices.clear();
175  m_LabeledSlices.resize( TImage::ImageDimension );
176  this->Modified();
177  }
179 
182  void
183  SetLabeledSliceIndices( unsigned int axis,
184  typename TImage::PixelType label,
185  const std::vector< typename TImage::IndexValueType >& indices )
186  {
187  SliceSetType sliceSet;
188  sliceSet.insert( indices.begin(), indices.end() );
189  m_LabeledSlices[axis][label] = sliceSet;
190  this->Modified();
191  }
193 
196  void
197  SetLabeledSliceIndices( unsigned int axis, typename TImage::PixelType label, const SliceSetType& indices )
198  {
199  m_LabeledSlices[axis][label] = indices;
200  this->Modified();
201  }
203 
206  GetLabeledSliceIndices( unsigned int axis, typename TImage::PixelType label )
207  {
208  return m_LabeledSlices[axis][label];
209  }
210 
211  // each label gets a set of slices in which it is present
212  typedef itksys::hash_map< typename TImage::PixelType, SliceSetType > LabeledSlicesType;
213  typedef std::vector< LabeledSlicesType > SliceIndicesType;
214 
218  {
219  return m_LabeledSlices;
220  }
221 
222 protected:
225  typename TImage::PixelType m_Label;
226  int m_Axis;
231  IdentifierType m_MinAlignIters; // minimum number of iterations in align method
232  IdentifierType m_MaxAlignIters; // maximum number of iterations in align method
233  IdentifierType m_ThreadCount; // for thread local instances
234  SliceIndicesType m_LabeledSlices; // one for each axis
235 
238  typedef Image< float, TImage::ImageDimension - 1 > FloatSliceType;
239  typedef Image< bool, TImage::ImageDimension - 1 > BoolSliceType;
240 
242  bool
243  ImagesEqual( typename BoolSliceType::Pointer& a, typename BoolSliceType::Pointer& b );
244 
246  virtual void
247  GenerateData() ITK_OVERRIDE;
248 
250  void
251  InterpolateBetweenTwo( int axis,
252  TImage* out,
253  typename TImage::PixelType label,
254  typename TImage::IndexValueType i,
255  typename TImage::IndexValueType j,
256  typename SliceType::Pointer& iconn,
257  typename SliceType::Pointer& jconn,
258  ThreadIdType threadId );
259 
265  void
266  InterpolateAlong( int axis, TImage* out );
267 
269  void
270  Extrapolate( int axis,
271  TImage* out,
272  typename TImage::PixelType label,
273  typename TImage::IndexValueType i,
274  typename TImage::IndexValueType j,
275  typename SliceType::Pointer& iConn,
276  typename TImage::PixelType iRegionId,
277  ThreadIdType threadId );
278 
280  typename FloatSliceType::Pointer
281  MaurerDM( typename BoolSliceType::Pointer& inImage, ThreadIdType threadId );
282 
285  std::vector< typename BoolSliceType::Pointer >
287  typename BoolSliceType::Pointer& end,
288  ThreadIdType threadId );
289 
291  typename BoolSliceType::Pointer
292  FindMedianImageDilations( typename BoolSliceType::Pointer& intersection,
293  typename BoolSliceType::Pointer& iMask,
294  typename BoolSliceType::Pointer& jMask,
295  ThreadIdType threadId );
296 
298  typename BoolSliceType::Pointer
299  FindMedianImageDistances( typename BoolSliceType::Pointer& intersection,
300  typename BoolSliceType::Pointer& iMask,
301  typename BoolSliceType::Pointer& jMask,
302  ThreadIdType threadId );
303 
305  void
306  Interpolate1to1( int axis,
307  TImage* out,
308  typename TImage::PixelType label,
309  typename TImage::IndexValueType i,
310  typename TImage::IndexValueType j,
311  typename SliceType::Pointer& iConn,
312  typename TImage::PixelType iRegionId,
313  typename SliceType::Pointer& jConn,
314  typename TImage::PixelType jRegionId,
315  const typename SliceType::IndexType& translation,
316  bool recursive,
317  ThreadIdType threadId );
318 
319  typedef std::vector< typename TImage::PixelType > PixelList;
320 
322  void
323  Interpolate1toN( int axis,
324  TImage* out,
325  typename TImage::PixelType label,
326  typename TImage::IndexValueType i,
327  typename TImage::IndexValueType j,
328  typename SliceType::Pointer& iConn,
329  typename TImage::PixelType iRegionId,
330  typename SliceType::Pointer& jConn,
331  const PixelList& jRegionIds,
332  const typename SliceType::IndexType& translation,
333  ThreadIdType threadId );
334 
336  typename SliceType::Pointer
337  TranslateImage( typename SliceType::Pointer& image,
338  const typename SliceType::IndexType& translation,
339  typename SliceType::RegionType newRegion );
340 
344  CardSymDifference( typename BoolSliceType::Pointer& shape1, typename BoolSliceType::Pointer& shape2 );
345 
347  virtual void
348  AllocateOutputs() ITK_OVERRIDE;
349 
351  typename SliceType::IndexType
352  Centroid( typename SliceType::Pointer& conn, const PixelList& regionIds );
353 
356  void
357  IntersectionRegions( const typename SliceType::IndexType& translation,
358  typename SliceType::RegionType& iRegion,
359  typename SliceType::RegionType& jRegion );
360 
363  Intersection( typename SliceType::Pointer& iConn,
364  typename TImage::PixelType iRegionId,
365  typename SliceType::Pointer& jConn,
366  const PixelList& jRegionIds,
367  const typename SliceType::IndexType& translation );
368 
370  typename SliceType::IndexType
371  Align( typename SliceType::Pointer& iConn,
372  typename TImage::PixelType iRegionId,
373  typename SliceType::Pointer& jConn,
374  const PixelList& jRegionIds );
375 
376  typedef itksys::hash_map< typename TImage::PixelType, typename TImage::RegionType > BoundingBoxesType;
377  BoundingBoxesType m_BoundingBoxes; // bounding box for each label
378 
380  typename SliceType::RegionType
381  BoundingBox( itk::SmartPointer< SliceType > image );
382 
386  template< typename T2 >
387  void
388  ExpandRegion( typename T2::RegionType& region, const typename T2::IndexType& index );
389 
391  typename SliceType::Pointer
392  RegionedConnectedComponents( const typename TImage::RegionType& region,
393  typename TImage::PixelType label,
394  IdentifierType& objectCount );
395 
397  typename BoolSliceType::Pointer
398  Dilate1( typename BoolSliceType::Pointer& seed, typename BoolSliceType::Pointer& mask, ThreadIdType threadId );
399 
401  typename RoiType::Pointer m_RoI;
402 
404  typename BinarizerType::Pointer m_Binarizer;
405 
407  typename ConnectedComponentsType::Pointer m_ConnectedComponents;
408 
409 private:
410  MorphologicalContourInterpolator( const Self & ) ITK_DELETE_FUNCTION;
411  void
412  operator=( const Self & ) ITK_DELETE_FUNCTION;
413 };
414 } // namespace itk
415 
416 #ifndef ITK_MANUAL_INSTANTIATION
417 #include "itkMorphologicalContourInterpolator.hxx"
418 #endif
419 
420 #endif // itkMorphologicalContourInterpolator_h
IdentifierType CardSymDifference(typename BoolSliceType::Pointer &shape1, typename BoolSliceType::Pointer &shape2)
SliceType::Pointer RegionedConnectedComponents(const typename TImage::RegionType &region, typename TImage::PixelType label, IdentifierType &objectCount)
Light weight base class for most itk classes.
Decrease the image size by cropping the image to the selected region bounds.
Image< float, TImage::ImageDimension-1 > FloatSliceType
void IntersectionRegions(const typename SliceType::IndexType &translation, typename SliceType::RegionType &iRegion, typename SliceType::RegionType &jRegion)
std::vector< typename TImage::PixelType > PixelList
itksys::hash_map< typename TImage::PixelType, typename TImage::RegionType > BoundingBoxesType
BoolSliceType::Pointer FindMedianImageDistances(typename BoolSliceType::Pointer &intersection, typename BoolSliceType::Pointer &iMask, typename BoolSliceType::Pointer &jMask, ThreadIdType threadId)
SliceSetType GetLabeledSliceIndices(unsigned int axis, typename TImage::PixelType label)
void Interpolate1toN(int axis, TImage *out, typename TImage::PixelType label, typename TImage::IndexValueType i, typename TImage::IndexValueType j, typename SliceType::Pointer &iConn, typename TImage::PixelType iRegionId, typename SliceType::Pointer &jConn, const PixelList &jRegionIds, const typename SliceType::IndexType &translation, ThreadIdType threadId)
signed long IndexValueType
Definition: itkIntTypes.h:150
virtual void AllocateOutputs() override
BoolSliceType::Pointer Dilate1(typename BoolSliceType::Pointer &seed, typename BoolSliceType::Pointer &mask, ThreadIdType threadId)
void Extrapolate(int axis, TImage *out, typename TImage::PixelType label, typename TImage::IndexValueType i, typename TImage::IndexValueType j, typename SliceType::Pointer &iConn, typename TImage::PixelType iRegionId, ThreadIdType threadId)
Image< bool, TImage::ImageDimension-1 > BoolSliceType
void SetLabeledSliceIndices(unsigned int axis, typename TImage::PixelType label, const std::vector< typename TImage::IndexValueType > &indices)
void Interpolate1to1(int axis, TImage *out, typename TImage::PixelType label, typename TImage::IndexValueType i, typename TImage::IndexValueType j, typename SliceType::Pointer &iConn, typename TImage::PixelType iRegionId, typename SliceType::Pointer &jConn, typename TImage::PixelType jRegionId, const typename SliceType::IndexType &translation, bool recursive, ThreadIdType threadId)
std::vector< typename BoolSliceType::Pointer > GenerateDilationSequence(typename BoolSliceType::Pointer &begin, typename BoolSliceType::Pointer &end, ThreadIdType threadId)
std::set< typename TImage::IndexValueType > SliceSetType
Binarize an input image by thresholding.
SizeValueType IdentifierType
Definition: itkIntTypes.h:147
bool ImagesEqual(typename BoolSliceType::Pointer &a, typename BoolSliceType::Pointer &b)
Image< typename TImage::PixelType, TImage::ImageDimension-1 > SliceType
BoolSliceType::Pointer FindMedianImageDilations(typename BoolSliceType::Pointer &intersection, typename BoolSliceType::Pointer &iMask, typename BoolSliceType::Pointer &jMask, ThreadIdType threadId)
SliceType::IndexType Centroid(typename SliceType::Pointer &conn, const PixelList &regionIds)
SliceType::Pointer TranslateImage(typename SliceType::Pointer &image, const typename SliceType::IndexType &translation, typename SliceType::RegionType newRegion)
FloatSliceType::Pointer MaurerDM(typename BoolSliceType::Pointer &inImage, ThreadIdType threadId)
Interpolates contours between slices. Based on a paper by Albu et al.
virtual void Modified() const
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
Image< bool, TImage::ImageDimension > BoolImageType
virtual void GenerateData() override
itksys::hash_map< typename TImage::PixelType, SliceSetType > LabeledSlicesType
void ExpandRegion(typename T2::RegionType &region, const typename T2::IndexType &index)
Base class for filters that take an image as input and produce an image as output.
void SetLabeledSliceIndices(unsigned int axis, typename TImage::PixelType label, const SliceSetType &indices)
SliceType::IndexType Align(typename SliceType::Pointer &iConn, typename TImage::PixelType iRegionId, typename SliceType::Pointer &jConn, const PixelList &jRegionIds)
void InterpolateBetweenTwo(int axis, TImage *out, typename TImage::PixelType label, typename TImage::IndexValueType i, typename TImage::IndexValueType j, typename SliceType::Pointer &iconn, typename SliceType::Pointer &jconn, ThreadIdType threadId)
Templated n-dimensional image class.
Definition: itkImage.h:75
Represent and compute information about bounding boxes.
IdentifierType Intersection(typename SliceType::Pointer &iConn, typename TImage::PixelType iRegionId, typename SliceType::Pointer &jConn, const PixelList &jRegionIds, const typename SliceType::IndexType &translation)
void InterpolateAlong(int axis, TImage *out)
Label the objects in a binary image.