ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkContourMeanDistanceImageFilter.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 itkContourMeanDistanceImageFilter_h
19 #define itkContourMeanDistanceImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkNumericTraits.h"
23 
24 namespace itk
25 {
68 template< typename TInputImage1, typename TInputImage2 >
69 class ITK_TEMPLATE_EXPORT ContourMeanDistanceImageFilter:
70  public ImageToImageFilter< TInputImage1, TInputImage1 >
71 {
72 public:
73  ITK_DISALLOW_COPY_AND_ASSIGN(ContourMeanDistanceImageFilter);
74 
80 
82  itkNewMacro(Self);
83 
86 
88  using InputImage1Type = TInputImage1;
89  using InputImage2Type = TInputImage2;
90  using InputImage1Pointer = typename TInputImage1::Pointer;
91  using InputImage2Pointer = typename TInputImage2::Pointer;
92  using InputImage1ConstPointer = typename TInputImage1::ConstPointer;
93  using InputImage2ConstPointer = typename TInputImage2::ConstPointer;
94 
96  using SizeType = typename TInputImage1::SizeType;
98 
99  using InputImage1PixelType = typename TInputImage1::PixelType;
100  using InputImage2PixelType = typename TInputImage2::PixelType;
101 
103  static constexpr unsigned int ImageDimension = TInputImage1::ImageDimension;
104 
107 
109  void SetInput1(const InputImage1Type *image);
110 
112  void SetInput2(const InputImage2Type *image);
113 
115  const InputImage1Type * GetInput1();
116 
118  const InputImage2Type * GetInput2();
119 
121  itkGetConstMacro(MeanDistance, RealType);
122 
124  itkSetMacro( UseImageSpacing, bool );
125  itkGetConstMacro( UseImageSpacing, bool );
127 
128 #ifdef ITK_USE_CONCEPT_CHECKING
129  // Begin concept checking
130  itkConceptMacro( InputHasNumericTraitsCheck,
132  // End concept checking
133 #endif
134 
135 protected:
137  ~ContourMeanDistanceImageFilter() override = default;
138  void PrintSelf(std::ostream & os, Indent indent) const override;
139 
141  void GenerateData() override;
142 
143  // Override since the filter needs all the data for the algorithm
144  void GenerateInputRequestedRegion() override;
145 
146  // Override since the filter produces all of its output
147  void EnlargeOutputRequestedRegion(DataObject *data) override;
148 
149 private:
152 }; // end of class
153 } // end namespace itk
154 
155 #ifndef ITK_MANUAL_INSTANTIATION
156 #include "itkContourMeanDistanceImageFilter.hxx"
157 #endif
158 
159 #endif
Light weight base class for most itk classes.
Define numeric traits for std::vector.
Computes the Mean distance between the boundaries of non-zero regions of two images.
typename NumericTraits< InputImage1PixelType >::RealType RealType
typename TInputImage1::ConstPointer InputImage1ConstPointer
typename TInputImage2::ConstPointer InputImage2ConstPointer
typename TInputImage1::PixelType InputImage1PixelType
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)
typename TInputImage2::PixelType InputImage2PixelType
Base class for all data objects in ITK.