ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkRecursiveSeparableImageFilter.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 itkRecursiveSeparableImageFilter_h
19 #define itkRecursiveSeparableImageFilter_h
20 
21 #include "itkInPlaceImageFilter.h"
22 #include "itkNumericTraits.h"
24 
25 namespace itk
26 {
50 template< typename TInputImage, typename TOutputImage = TInputImage >
51 class ITK_TEMPLATE_EXPORT RecursiveSeparableImageFilter:
52  public InPlaceImageFilter< TInputImage, TOutputImage >
53 {
54 public:
55  ITK_DISALLOW_COPY_AND_ASSIGN(RecursiveSeparableImageFilter);
56 
62 
65 
67  using InputImagePointer = typename TInputImage::Pointer;
68  using InputImageConstPointer = typename TInputImage::ConstPointer;
69 
75  using InputPixelType = typename TInputImage::PixelType;
78 
80 
82  using InputImageType = TInputImage;
83 
85  using OutputImageType = TOutputImage;
86 
88  itkGetConstMacro(Direction, unsigned int);
89 
91  itkSetMacro(Direction, unsigned int);
92 
94  void SetInputImage(const TInputImage *);
95 
97  const TInputImage * GetInputImage();
98 
99 protected:
101  ~RecursiveSeparableImageFilter() override = default;
102  void PrintSelf(std::ostream & os, Indent indent) const override;
103 
104  void BeforeThreadedGenerateData() override;
105 
106  void GenerateData() override;
107 
108  void DynamicThreadedGenerateData( const OutputImageRegionType & ) override;
109 
118  void EnlargeOutputRequestedRegion(DataObject *output) override;
119 
124  virtual void SetUp(ScalarRealType spacing) = 0;
125 
132  void FilterDataArray(RealType *outs, const RealType *data, RealType *scratch,
133  SizeValueType ln);
134 
135 protected:
141 
149 
155 
162 
167 
168 
169  template <typename T1, typename T2>
170  inline void MathEMAMAMAM(T1 &out,
171  const T1 &a1, const T2 &b1,
172  const T1 &a2, const T2 &b2,
173  const T1 &a3, const T2 &b3,
174  const T1 &a4, const T2 &b4 )
175  {
176  out = a1*b1 + a2*b2 + a3*b3 + a4*b4;
177  }
178 
179 
180  template <typename T1, typename T2>
182  const VariableLengthVector<T1> &a1, const T2 &b1,
183  const VariableLengthVector<T1> &a2, const T2 &b2,
184  const VariableLengthVector<T1> &a3, const T2 &b3,
185  const VariableLengthVector<T1> &a4, const T2 &b4 )
186  {
187  const unsigned int sz = a1.GetSize();
188  if (sz != out.GetSize() )
189  {
190  out.SetSize(sz);
191  }
192  for ( unsigned int i = 0; i < sz; ++i)
193  {
194  out[i] = a1[i]*b1 + a2[i]*b2 + a3[i]*b3 + a4[i]*b4;
195  }
196  }
197 
198  template <typename T1, typename T2>
199  inline void MathSMAMAMAM(T1 &out,
200  const T1 &a1, const T2 &b1,
201  const T1 &a2, const T2 &b2,
202  const T1 &a3, const T2 &b3,
203  const T1 &a4, const T2 &b4 )
204  {
205  out -= a1*b1 + a2*b2 + a3*b3 + a4*b4;
206  }
207 
208  template <typename T1, typename T2>
210  const VariableLengthVector<T1> &a1, const T2 &b1,
211  const VariableLengthVector<T1> &a2, const T2 &b2,
212  const VariableLengthVector<T1> &a3, const T2 &b3,
213  const VariableLengthVector<T1> &a4, const T2 &b4 )
214  {
215  const unsigned int sz = a1.GetSize();
216  if (sz != out.GetSize() )
217  {
218  out.SetSize(sz);
219  }
220  for ( unsigned int i = 0; i < sz; ++i)
221  {
222  out[i] -= a1[i]*b1 + a2[i]*b2 + a3[i]*b3 + a4[i]*b4;
223  }
224  }
225 
226 private:
229  unsigned int m_Direction{ 0 };
230 };
231 } // end namespace itk
232 
233 #ifndef ITK_MANUAL_INSTANTIATION
234 #include "itkRecursiveSeparableImageFilter.hxx"
235 #endif
236 
237 #endif
Base class for recursive convolution with a kernel.
Define numeric traits for std::vector.
unsigned long SizeValueType
Definition: itkIntTypes.h:83
void MathSMAMAMAM(VariableLengthVector< T1 > &out, const VariableLengthVector< T1 > &a1, const T2 &b1, const VariableLengthVector< T1 > &a2, const T2 &b2, const VariableLengthVector< T1 > &a3, const T2 &b3, const VariableLengthVector< T1 > &a4, const T2 &b4)
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.
void MathSMAMAMAM(T1 &out, const T1 &a1, const T2 &b1, const T1 &a2, const T2 &b2, const T1 &a3, const T2 &b3, const T1 &a4, const T2 &b4)
typename InputImageType::Pointer InputImagePointer
Represents an array whose length can be defined at run-time.
typename OutputImageType::RegionType OutputImageRegionType
typename NumericTraits< InputPixelType >::RealType RealType
TOutputImage OutputImageType
void MathEMAMAMAM(VariableLengthVector< T1 > &out, const VariableLengthVector< T1 > &a1, const T2 &b1, const VariableLengthVector< T1 > &a2, const T2 &b2, const VariableLengthVector< T1 > &a3, const T2 &b3, const VariableLengthVector< T1 > &a4, const T2 &b4)
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename NumericTraits< InputPixelType >::ScalarRealType ScalarRealType
void MathEMAMAMAM(T1 &out, const T1 &a1, const T2 &b1, const T1 &a2, const T2 &b2, const T1 &a3, const T2 &b3, const T1 &a4, const T2 &b4)
void SetSize(unsigned int sz, TReallocatePolicy reallocatePolicy, TKeepValuesPolicy keepValues)
Base class for filters that take an image as input and overwrite that image as the output...
typename InputImageType::ConstPointer InputImageConstPointer
Base class for all data objects in ITK.
typename TInputImage::PixelType InputPixelType