ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkMaximumProjectionImageFilter.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
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 * https://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 itkMaximumProjectionImageFilter_h
19#define itkMaximumProjectionImageFilter_h
20
22#include "itkConceptChecking.h"
23
24namespace itk
25{
44
45namespace Functor
46{
47template <typename TInputPixel>
49{
50public:
53
54 inline void
56 {
57
58 // check if scalar or fixed length array type
59 if constexpr (std::is_same<TInputPixel, typename NumericTraits<TInputPixel>::ValueType>::value)
60 {
62 }
63 else
64 {
65 m_Maximum = TInputPixel();
67 }
68 }
69
70 inline void
71 operator()(const TInputPixel & input)
72 {
73 if constexpr (std::is_same<TInputPixel, typename NumericTraits<TInputPixel>::ValueType>::value)
74 {
75 m_Maximum = std::max(m_Maximum, input);
76 }
77 else
78 {
80 {
81 m_Maximum = input;
82 }
83 else
84 {
85 for (unsigned int i = 0; i < itk::NumericTraits<TInputPixel>::GetLength(m_Maximum); ++i)
86 {
87 m_Maximum[i] = std::max(m_Maximum[i], input[i]);
88 }
89 }
90 }
91 }
92
93 inline TInputPixel
95 {
96 return m_Maximum;
97 }
98
99 TInputPixel m_Maximum;
100};
101} // namespace Functor
102
103template <typename TInputImage, typename TOutputImage>
105 : public ProjectionImageFilter<TInputImage,
106 TOutputImage,
107 Functor::MaximumAccumulator<typename TInputImage::PixelType>>
108{
109public:
110 ITK_DISALLOW_COPY_AND_MOVE(MaximumProjectionImageFilter);
111
115
116 using InputImageType = TInputImage;
117 using InputPixelType = typename InputImageType::PixelType;
118
121
123 itkOverrideGetNameOfClassMacro(MaximumProjectionImageFilter);
124
126 itkNewMacro(Self);
127
128 itkConceptMacro(InputPixelTypeGreaterThanComparable,
130 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
131
132protected:
134 ~MaximumProjectionImageFilter() override = default;
135}; // end
136 // MaximumProjectionImageFilter
137} // end namespace itk
138
139#endif
~MaximumProjectionImageFilter() override=default
ProjectionImageFilter< TInputImage, TOutputImage, Functor::MaximumAccumulator< typename TInputImage::PixelType > > Superclass
typename InputImageType::PixelType InputPixelType
static constexpr T NonpositiveMin()
static unsigned int GetLength(const T &)
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86