ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkContourDirectedMeanDistanceImageFilter.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 itkContourDirectedMeanDistanceImageFilter_h
19 #define itkContourDirectedMeanDistanceImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNumericTraits.h"
23 #include "itkArray.h"
24 #include "itkImage.h"
25 
26 namespace itk
27 {
62 template< typename TInputImage1, typename TInputImage2 >
63 class ITK_TEMPLATE_EXPORT ContourDirectedMeanDistanceImageFilter:
64  public ImageToImageFilter< TInputImage1, TInputImage1 >
65 {
66 public:
67  ITK_DISALLOW_COPY_AND_ASSIGN(ContourDirectedMeanDistanceImageFilter);
68 
74 
76  itkNewMacro(Self);
77 
80 
82  using InputImage1Type = TInputImage1;
83  using InputImage2Type = TInputImage2;
84  using InputImage1Pointer = typename TInputImage1::Pointer;
85  using InputImage2Pointer = typename TInputImage2::Pointer;
86  using InputImage1ConstPointer = typename TInputImage1::ConstPointer;
87  using InputImage2ConstPointer = typename TInputImage2::ConstPointer;
88 
90  using SizeType = typename TInputImage1::SizeType;
92 
93  using InputImage1PixelType = typename TInputImage1::PixelType;
94  using InputImage2PixelType = typename TInputImage2::PixelType;
95 
97  static constexpr unsigned int ImageDimension = TInputImage1::ImageDimension;
98 
101 
103  void SetInput1(const InputImage1Type *image);
104 
106  void SetInput2(const InputImage2Type *image);
107 
109  const InputImage1Type * GetInput1();
110 
112  const InputImage2Type * GetInput2();
113 
115  itkGetConstMacro(ContourDirectedMeanDistance, RealType);
116 
118  itkSetMacro( UseImageSpacing, bool );
119  itkGetConstMacro( UseImageSpacing, bool );
121 
122 #ifdef ITK_USE_CONCEPT_CHECKING
123  // Begin concept checking
124  itkConceptMacro( InputHasNumericTraitsCheck,
126  // End concept checking
127 #endif
128 
129 protected:
131  ~ContourDirectedMeanDistanceImageFilter() override = default;
132  void PrintSelf(std::ostream & os, Indent indent) const override;
133 
136  void AllocateOutputs() override;
137 
139  void BeforeThreadedGenerateData() override;
140 
143  void AfterThreadedGenerateData() override;
144 
146  void ThreadedGenerateData(const RegionType &
147  outputRegionForThread,
148  ThreadIdType threadId) override;
149 
150  void DynamicThreadedGenerateData( const RegionType & ) override
151  {
152  itkExceptionMacro("This class requires threadId so it must use classic multi-threading model");
153  }
154 
155  // Override since the filter needs all the data for the algorithm
156  void GenerateInputRequestedRegion() override;
157 
158  // Override since the filter produces all of its output
159  void EnlargeOutputRequestedRegion(DataObject *data) override;
160 
161 private:
163 
165 
170 }; // end of class
171 } // end namespace itk
172 
173 #ifndef ITK_MANUAL_INSTANTIATION
174 #include "itkContourDirectedMeanDistanceImageFilter.hxx"
175 #endif
176 
177 #endif
Light weight base class for most itk classes.
Define numeric traits for std::vector.
typename NumericTraits< InputImage1PixelType >::RealType RealType
Computes the directed Mean distance between the boundaries of non-zero pixel regions of two images...
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75