ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMagnitudeAndPhaseToComplexImageFilter.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 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  MagnitudeAndPhaseToComplex() = default;
62  ~MagnitudeAndPhaseToComplex() = default;
64  {
65  return false;
66  }
67 
68  bool operator==(const MagnitudeAndPhaseToComplex & other) const
69  {
70  return !( *this != other );
71  }
72 
73  inline std::complex< TOutput > operator()(const TInput1 & A, const TInput2 & B) const
74  {
75  return std::complex< TOutput >( std::polar( static_cast< TOutput >( A ), static_cast< TOutput >( B ) ) );
76  }
77 };
78 }
79 
80 template< typename TInputImage1,
81  typename TInputImage2 = TInputImage1,
83  TInputImage1::ImageDimension > >
85  public BinaryGeneratorImageFilter< TInputImage1,
86  TInputImage2,
87  TOutputImage>
88 {
89 public:
90  ITK_DISALLOW_COPY_AND_ASSIGN(MagnitudeAndPhaseToComplexImageFilter);
91 
94 
95  using Superclass = BinaryGeneratorImageFilter< TInputImage1,
96  TInputImage2,
97  TOutputImage >;
98 
99  using InputPixel1Type = typename TInputImage1::PixelType;
100  using InputPixel2Type = typename TInputImage2::PixelType;
101  using OutputPixelType = typename TOutputImage::PixelType;
102 
105 
106  using FunctorType = Functor::MagnitudeAndPhaseToComplex< typename TInputImage1::PixelType,
107  typename TInputImage2::PixelType,
108  typename TOutputImage::PixelType::value_type >;
109 
111  itkNewMacro(Self);
112 
115 
116 #ifdef ITK_USE_CONCEPT_CHECKING
117  // Begin concept checking
118  itkConceptMacro( Input1ConvertibleToDoubleCheck,
120  itkConceptMacro( Input2ConvertibleToDoubleCheck,
122  itkConceptMacro( DoubleConvertibleToOutputCheck,
124  // End concept checking
125 #endif
126 
127 protected:
129  {
130 #if !defined( ITK_WRAPPING_PARSER )
132 #endif
133  }
134 
135  ~MagnitudeAndPhaseToComplexImageFilter() override = default;
136 };
137 } // end namespace itk
138 
139 #endif
~MagnitudeAndPhaseToComplexImageFilter() override=default
std::complex< TOutput > operator()(const TInput1 &A, const TInput2 &B) const
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
bool operator!=(const MagnitudeAndPhaseToComplex &) const
Implements pixel-wise generic operation of two images, or of an image and a constant.
void SetFunctor(const std::function< ConstRefFunctionType > &f)
Implements pixel-wise conversion of magnitude and phase data into complex voxels. ...
Functor::MagnitudeAndPhaseToComplex< typename TInputImage1::PixelType, typename TInputImage2::PixelType, typename TOutputImage::PixelType::value_type > FunctorType
bool operator==(const MagnitudeAndPhaseToComplex &other) const
#define itkConceptMacro(name, concept)
Templated n-dimensional image class.
Definition: itkImage.h:75