ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkMinimumMaximumImageFilter.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 __itkMinimumMaximumImageFilter_h
19 #define __itkMinimumMaximumImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 
24 #include <vector>
25 
26 #include "itkNumericTraits.h"
27 
28 namespace itk
29 {
43 template< typename TInputImage >
45  public ImageToImageFilter< TInputImage, TInputImage >
46 {
47 public:
49  itkStaticConstMacro(InputImageDimension, unsigned int,
50  TInputImage::ImageDimension);
51  itkStaticConstMacro(OutputImageDimension, unsigned int,
52  TInputImage::ImageDimension);
54 
60 
62  typedef typename TInputImage::Pointer InputImagePointer;
63 
64  typedef typename TInputImage::RegionType RegionType;
65  typedef typename TInputImage::SizeType SizeType;
66  typedef typename TInputImage::IndexType IndexType;
67  typedef typename TInputImage::PixelType PixelType;
68 
71 
73  itkNewMacro(Self);
74 
77 
79  typedef TInputImage InputImageType;
80 
83 
86  { return this->GetMinimumOutput()->Get(); }
89 
90  const PixelObjectType * GetMinimumOutput() const;
91 
94  { return this->GetMaximumOutput()->Get(); }
97 
98  const PixelObjectType * GetMaximumOutput() const;
99 
105 
106 #ifdef ITK_USE_CONCEPT_CHECKING
107  // Begin concept checking
108  itkConceptMacro( LessThanComparableCheck,
110  itkConceptMacro( GreaterThanComparableCheck,
112  itkConceptMacro( OStreamWritableCheck,
114  // End concept checking
115 #endif
116 
117 protected:
120  void PrintSelf(std::ostream & os, Indent indent) const;
121 
124  void AllocateOutputs();
125 
128 
132 
134  void ThreadedGenerateData(const RegionType &
135  outputRegionForThread,
136  ThreadIdType threadId);
137 
138  // Override since the filter needs all the data for the algorithm
140 
141  // Override since the filter produces all of its output
143 
144 private:
145  MinimumMaximumImageFilter(const Self &); //purposely not implemented
146  void operator=(const Self &); //purposely not implemented
147 
148  std::vector< PixelType > m_ThreadMin;
149  std::vector< PixelType > m_ThreadMax;
150 };
151 } // end namespace itk
152 
153 #ifndef ITK_MANUAL_INSTANTIATION
154 #include "itkMinimumMaximumImageFilter.hxx"
155 #endif
156 
157 #endif
SimpleDataObjectDecorator< PixelType > PixelObjectType
Light weight base class for most itk classes.
PixelObjectType * GetMinimumOutput()
Computes the minimum and the maximum intensity values of an image.
virtual ProcessObject::DataObjectPointer MakeOutput(ProcessObject::DataObjectPointerArraySizeType idx) ITK_OVERRIDE
PixelObjectType * GetMaximumOutput()
Decorates any &quot;simple&quot; data type (data types without smart pointers) with a DataObject API...
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
static const unsigned int OutputImageDimension
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx)
ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
void operator=(const Self &)
void ThreadedGenerateData(const RegionType &outputRegionForThread, ThreadIdType threadId)
void EnlargeOutputRequestedRegion(DataObject *data)
ImageToImageFilter< TInputImage, TInputImage > Superclass
static const unsigned int InputImageDimension
Base class for filters that take an image as input and produce an image as output.
void PrintSelf(std::ostream &os, Indent indent) const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.
unsigned int ThreadIdType
Definition: itkIntTypes.h:159