ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkBinaryGrindPeakImageFilter.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 itkBinaryGrindPeakImageFilter_h
19 #define itkBinaryGrindPeakImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk {
24 
44 template<typename TInputImage>
46  public ImageToImageFilter<TInputImage, TInputImage>
47 {
48 public:
54 
56  typedef TInputImage InputImageType;
57  typedef TInputImage OutputImageType;
58  typedef typename InputImageType::Pointer InputImagePointer;
59  typedef typename InputImageType::ConstPointer InputImageConstPointer;
60  typedef typename InputImageType::RegionType InputImageRegionType;
61  typedef typename InputImageType::PixelType InputImagePixelType;
62  typedef typename OutputImageType::Pointer OutputImagePointer;
63  typedef typename OutputImageType::ConstPointer OutputImageConstPointer;
64  typedef typename OutputImageType::RegionType OutputImageRegionType;
65  typedef typename OutputImageType::PixelType OutputImagePixelType;
66 
68  itkStaticConstMacro(InputImageDimension, unsigned int,
69  TInputImage::ImageDimension);
70  itkStaticConstMacro(OutputImageDimension, unsigned int,
71  OutputImageType::ImageDimension);
73 
75  itkNewMacro(Self);
76 
78  itkTypeMacro(BinaryGrindPeakImageFilter,
80 
87  itkSetMacro(FullyConnected, bool);
88  itkGetConstReferenceMacro(FullyConnected, bool);
89  itkBooleanMacro(FullyConnected);
91 
92 #ifdef ITK_USE_CONCEPT_CHECKING
93  // Begin concept checking
94  itkConceptMacro(InputOStreamWritableCheck,
96  // End concept checking
97 #endif
98 
101  itkSetMacro(ForegroundValue, InputImagePixelType);
102 
105  itkGetMacro(ForegroundValue, InputImagePixelType);
106 
108  itkSetMacro(BackgroundValue, InputImagePixelType);
109 
111  itkGetMacro(BackgroundValue, InputImagePixelType);
112 
113 protected:
116  void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
117 
121  void GenerateInputRequestedRegion() ITK_OVERRIDE;
122 
124  void EnlargeOutputRequestedRegion(DataObject *itkNotUsed(output)) ITK_OVERRIDE;
125 
128  void GenerateData() ITK_OVERRIDE;
129 
130 private:
131  BinaryGrindPeakImageFilter(const Self&); //purposely not implemented
132  void operator=(const Self&); //purposely not implemented
133 
135 
137 
139 
140 }; // end of class
141 
142 } // end namespace itk
143 
144 #ifndef ITK_MANUAL_INSTANTIATION
145 #include "itkBinaryGrindPeakImageFilter.hxx"
146 #endif
147 
148 #endif
ImageToImageFilter< TInputImage, TInputImage > Superclass
Light weight base class for most itk classes.
static const unsigned int OutputImageDimension
InputImageType::ConstPointer InputImageConstPointer
OutputImageType::ConstPointer OutputImageConstPointer
Remove the objects not connected to the boundary of the image.
void PrintSelf(std::ostream &os, Indent indent) const override
void EnlargeOutputRequestedRegion(DataObject *) override
void GenerateInputRequestedRegion() override
OutputImageType::PixelType OutputImagePixelType
Base class for filters that take an image as input and produce an image as output.
OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.
InputImageType::RegionType InputImageRegionType