1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143
|
/*=========================================================================
Program: Insight Segmentation & Registration Toolkit
Module: $RCSfile: itkNaryMaximumImageFilter.h,v $
Language: C++
Date: $Date: 2007-09-27 11:36:41 $
Version: $Revision: 1.8 $
Copyright (c) Insight Software Consortium. All rights reserved.
See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#ifndef __itkNaryMaximumImageFilter_h
#define __itkNaryMaximumImageFilter_h
#include "itkNaryFunctorImageFilter.h"
#include "itkNumericTraits.h"
namespace itk
{
/** \class NaryMaximumImageFilter
* \brief Implements an operator computing the pixel-wise maximum of several images.
*
* This class is parametrized over the types of the input images and the type
* of the output image. Numeric conversions (castings) are done by the C++
* defaults.
*
* The pixel type of the output images must have a valid defintion of the
* operator<. This condition is required because internally this filter will
* perform an operation similar to:
*
* const OutputPixelType query_value = static_cast<OutputPixelType>(pixel_from_input_n);
* if(current_maximum < query_value)
* {
* current_maximum = query_value;
* }
* (where current_maximum is also of type OutputPixelType)
*
* for each of the n input images.
*
* For example, this filter could be used directly to find a "maximum projection"
* of a series of images, often used in preliminary analysis of time-series data.
*
* \author Zachary Pincus
*
* This filter was contributed by Zachary Pincus from the Department of
* Biochemistry and Program in Biomedical Informatics at Stanford University
* School of Medicine
*
* \ingroup IntensityImageFilters Multithreaded
*/
namespace Functor {
template< class TInput, class TOutput >
class Maximum1
{
public:
typedef typename NumericTraits< TOutput >::ValueType OutputValueType;
// not sure if this typedef really makes things more clear... could just use TOutput?
Maximum1() {}
~Maximum1() {}
inline TOutput operator()( const std::vector< TInput > & B)
{
OutputValueType A = NumericTraits< TOutput >::NonpositiveMin();
for( unsigned int i=0; i< B.size(); i++ )
{
if( A < static_cast<OutputValueType>(B[i]) )
{
A = static_cast< OutputValueType > (B[i]);
}
}
return A;
}
bool operator== (const Maximum1&) const
{
return true;
}
bool operator!= (const Maximum1&) const
{
return false;
}
};
}
template <class TInputImage, class TOutputImage>
class ITK_EXPORT NaryMaximumImageFilter :
public
NaryFunctorImageFilter<TInputImage,TOutputImage,
Functor::Maximum1< typename TInputImage::PixelType,
typename TInputImage::PixelType > >
{
public:
/** Standard class typedefs. */
typedef NaryMaximumImageFilter Self;
typedef NaryFunctorImageFilter<TInputImage,TOutputImage,
Functor::Maximum1< typename TInputImage::PixelType,
typename TInputImage::PixelType > > Superclass;
typedef SmartPointer<Self> Pointer;
typedef SmartPointer<const Self> ConstPointer;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** Runtime information support. */
itkTypeMacro(NaryMaximumImageFilter,
NaryFunctorImageFilter);
#ifdef ITK_USE_CONCEPT_CHECKING
/** Begin concept checking */
itkConceptMacro(InputConvertibleToOutputCheck,
(Concept::Convertible<typename TInputImage::PixelType,
typename TOutputImage::PixelType>));
itkConceptMacro(InputLessThanComparableCheck,
(Concept::LessThanComparable<typename TInputImage::PixelType>));
itkConceptMacro(InputHasNumericTraitsCheck,
(Concept::HasNumericTraits<typename TInputImage::PixelType>));
/** End concept checking */
#endif
protected:
NaryMaximumImageFilter() {}
virtual ~NaryMaximumImageFilter() {}
private:
NaryMaximumImageFilter(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
};
} // end namespace itk
#endif
|