ITK  6.0.0
Insight Toolkit
itkMagnitudeAndPhaseToComplexImageFilter.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 itkMagnitudeAndPhaseToComplexImageFilter_h
19 #define itkMagnitudeAndPhaseToComplexImageFilter_h
20 
22 
23 namespace itk
24 {
55 namespace Functor
56 {
57 template <typename TInput1, typename TInput2, typename TOutput>
59 {
60 public:
61  bool
63  {
64  return true;
65  }
66 
68 
69  inline std::complex<TOutput>
70  operator()(const TInput1 & A, const TInput2 & B) const
71  {
72  return std::complex<TOutput>(std::polar(static_cast<TOutput>(A), static_cast<TOutput>(B)));
73  }
74 };
75 } // namespace Functor
76 
77 template <typename TInputImage1,
78  typename TInputImage2 = TInputImage1,
79  typename TOutputImage =
80  itk::Image<std::complex<typename TInputImage1::PixelType>, TInputImage1::ImageDimension>>
82  : public BinaryGeneratorImageFilter<TInputImage1, TInputImage2, TOutputImage>
83 {
84 public:
85  ITK_DISALLOW_COPY_AND_MOVE(MagnitudeAndPhaseToComplexImageFilter);
86 
89 
91 
92  using InputPixel1Type = typename TInputImage1::PixelType;
93  using InputPixel2Type = typename TInputImage2::PixelType;
94  using OutputPixelType = typename TOutputImage::PixelType;
95 
98 
99  using FunctorType = Functor::MagnitudeAndPhaseToComplex<typename TInputImage1::PixelType,
100  typename TInputImage2::PixelType,
101  typename TOutputImage::PixelType::value_type>;
102 
104  itkNewMacro(Self);
105 
107  itkOverrideGetNameOfClassMacro(MagnitudeAndPhaseToComplexImageFilter);
108 
109 #ifdef ITK_USE_CONCEPT_CHECKING
110  // Begin concept checking
111  itkConceptMacro(Input1ConvertibleToDoubleCheck, (Concept::Convertible<InputPixel1Type, double>));
112  itkConceptMacro(Input2ConvertibleToDoubleCheck, (Concept::Convertible<InputPixel2Type, double>));
113  itkConceptMacro(DoubleConvertibleToOutputCheck, (Concept::Convertible<double, OutputPixelType>));
114  // End concept checking
115 #endif
116 
117 protected:
119  {
120 #if !defined(ITK_WRAPPING_PARSER)
122 #endif
123  }
124 
125  ~MagnitudeAndPhaseToComplexImageFilter() override = default;
126 };
127 } // end namespace itk
128 
129 #endif
itk::MagnitudeAndPhaseToComplexImageFilter::~MagnitudeAndPhaseToComplexImageFilter
~MagnitudeAndPhaseToComplexImageFilter() override=default
itk::MagnitudeAndPhaseToComplexImageFilter::MagnitudeAndPhaseToComplexImageFilter
MagnitudeAndPhaseToComplexImageFilter()
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:118
itk::BinaryGeneratorImageFilter
Implements pixel-wise generic operation of two images, or of an image and a constant.
Definition: itkBinaryGeneratorImageFilter.h:56
itk::SmartPointer< Self >
itk::MagnitudeAndPhaseToComplexImageFilter
Implements pixel-wise conversion of magnitude and phase data into complex voxels.
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:81
itk::Functor::MagnitudeAndPhaseToComplex
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:58
itk::MagnitudeAndPhaseToComplexImageFilter::InputPixel2Type
typename TInputImage2::PixelType InputPixel2Type
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:93
itkBinaryGeneratorImageFilter.h
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::Functor::MagnitudeAndPhaseToComplex::operator==
bool operator==(const MagnitudeAndPhaseToComplex &) const
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:62
itk::Functor::MagnitudeAndPhaseToComplex::operator()
std::complex< TOutput > operator()(const TInput1 &A, const TInput2 &B) const
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:70
itk::Functor::MagnitudeAndPhaseToComplex::ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION
ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION(MagnitudeAndPhaseToComplex)
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::MagnitudeAndPhaseToComplexImageFilter::FunctorType
Functor::MagnitudeAndPhaseToComplex< typename TInputImage1::PixelType, typename TInputImage2::PixelType, typename TOutputImage::PixelType::value_type > FunctorType
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:101
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::MagnitudeAndPhaseToComplexImageFilter::InputPixel1Type
typename TInputImage1::PixelType InputPixel1Type
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:92
itk::MagnitudeAndPhaseToComplexImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkMagnitudeAndPhaseToComplexImageFilter.h:94
itk::BinaryGeneratorImageFilter::SetFunctor
void SetFunctor(const std::function< ConstRefFunctionType > &f)
Definition: itkBinaryGeneratorImageFilter.h:150