isce3  0.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
Signal.h
1 // -*- C++ -*-
2 // -*- coding: utf-8 -*-
3 //
4 // Author: Heresh Fattahi
5 // Copyright 2018-
6 //
7 
8 #pragma once
9 
10 #include "forward.h"
11 
12 #include <cmath>
13 #include <memory>
14 #include <valarray>
15 
16 #include <isce3/core/Constants.h>
17 
20 template<class T>
22  public:
24  Signal();
25 
27  Signal(int nthreads);
28 
29  ~Signal() {};
30 
34  void fftPlanForward(std::valarray<std::complex<T>> &input,
35  std::valarray<std::complex<T>> &output,
36  int rank, int* n, int howmany,
37  int* inembed, int istride, int idist,
38  int* onembed, int ostride, int odist, int sign);
39 
43  void fftPlanForward(std::complex<T>* input,
44  std::complex<T>* output,
45  int rank, int* n, int howmany,
46  int* inembed, int istride, int idist,
47  int* onembed, int ostride, int odist, int sign);
48 
52  void fftPlanForward(std::valarray<T> &input,
53  std::valarray<std::complex<T>> &output,
54  int rank, int* n, int howmany,
55  int* inembed, int istride, int idist,
56  int* onembed, int ostride, int odist);
57 
61  void fftPlanForward(T* input,
62  std::complex<T>* output,
63  int rank, int* n, int howmany,
64  int* inembed, int istride, int idist,
65  int* onembed, int ostride, int odist);
66 
70  void fftPlanBackward(std::valarray<std::complex<T>> &input,
71  std::valarray<std::complex<T>> &output,
72  int rank, int* n, int howmany,
73  int* inembed, int istride, int idist,
74  int* onembed, int ostride, int odist, int sign);
75 
79  void fftPlanBackward(std::complex<T>* input,
80  std::complex<T>* output,
81  int rank, int* n, int howmany,
82  int* inembed, int istride, int idist,
83  int* onembed, int ostride, int odist, int sign);
84 
89  void fftPlanBackward(std::valarray<std::complex<T>>& input,
90  std::valarray<T>& output,
91  int rank, int* n, int howmany,
92  int* inembed, int istride, int idist,
93  int* onembed, int ostride, int odist);
94 
99  void fftPlanBackward(std::complex<T>* input,
100  T* output,
101  int rank, int* n, int howmany,
102  int* inembed, int istride, int idist,
103  int* onembed, int ostride, int odist);
104 
106  void forward(std::valarray<std::complex<T>> &input,
107  std::valarray<std::complex<T>> &output);
108 
110  void forward(std::complex<T>* input,
111  std::complex<T>* output);
112 
114  void forward(std::valarray<T> &input,
115  std::valarray<std::complex<T>> &output);
116 
118  void forward(T* input,
119  std::complex<T>* output);
120 
121 
122 
124  void inverse(std::valarray<std::complex<T>> &input,
125  std::valarray<std::complex<T>> &output);
126 
128  void inverse(std::complex<T>* input,
129  std::complex<T>* output);
130 
132  void inverse(std::valarray<std::complex<T>>& input,
133  std::valarray<T>& output);
134 
136  void inverse(std::complex<T>* input,
137  T* output);
138 
144  void forwardRangeFFT(std::valarray<std::complex<T>>& signal,
145  std::valarray<std::complex<T>>& spectrum,
146  int ncolumns, int nrows);
147 
153  void forwardRangeFFT(std::complex<T>* signal,
154  std::complex<T>* spectrum,
155  int ncolumns, int nrows);
156 
162  void forwardRangeFFT(std::valarray<T>& signal,
163  std::valarray<std::complex<T>>& spectrum,
164  int ncolumns, int nrows);
165 
171  void forwardRangeFFT(T* signal,
172  std::complex<T>* spectrum,
173  int ncolumns, int nrows);
174 
180  void forwardAzimuthFFT(std::valarray<std::complex<T>> &signal,
181  std::valarray<std::complex<T>> &spectrum,
182  int ncolumns, int nrows);
183 
189  void forwardAzimuthFFT(std::complex<T>* signal,
190  std::complex<T>* spectrum,
191  int ncolumns, int nrows);
192 
198  void forwardAzimuthFFT(std::valarray<T> &signal,
199  std::valarray<std::complex<T>> &spectrum,
200  int ncolumns, int nrows);
201 
207  void forwardAzimuthFFT(T* signal,
208  std::complex<T>* spectrum,
209  int ncolumns, int nrows);
210 
213  void forward2DFFT(std::valarray<std::complex<T>> &signal,
214  std::valarray<std::complex<T>> &spectrum,
215  int ncolumns, int nrows);
216 
219  void forward2DFFT(std::valarray<std::complex<T>> &signal,
220  std::valarray<std::complex<T>> &spectrum,
221  int incolumns, int inrows,
222  int oncolumns, int onrows);
223 
226  void forward2DFFT(std::complex<T>* signal,
227  std::complex<T>* spectrum,
228  int ncolumns, int nrows);
229 
232  void forward2DFFT(std::complex<T>* signal,
233  std::complex<T>* spectrum,
234  int incolumns, int inrows,
235  int oncolumns, int onrows);
236 
239  void forward2DFFT(std::valarray<T>& signal,
240  std::valarray<std::complex<T>>& spectrum,
241  int ncolumns, int nrows);
242 
245  void forward2DFFT(std::valarray<T>& signal,
246  std::valarray<std::complex<T>>& spectrum,
247  int incolumns, int inrows,
248  int oncolumns, int onrows);
249 
252  void forward2DFFT(T* signal,
253  std::complex<T>* spectrum,
254  int ncolumns, int nrows);
255 
258  void forward2DFFT(T* signal,
259  std::complex<T>* spectrum,
260  int incolumns, int inrows,
261  int oncolumns, int onrows);
262 
265  void inverseRangeFFT(std::valarray<std::complex<T>> &spectrum,
266  std::valarray<std::complex<T>> &signal,
267  int ncolumns, int nrows);
268 
271  void inverseRangeFFT(std::complex<T>* spectrum,
272  std::complex<T>* signal,
273  int ncolumns, int nrows);
274 
277  void inverseRangeFFT(std::valarray<std::complex<T>> &spectrum,
278  std::valarray<T> &signal,
279  int ncolumns, int nrows);
280 
283  void inverseRangeFFT(std::complex<T>* spectrum,
284  T* signal,
285  int ncolumns, int nrows);
286 
289  void inverseAzimuthFFT(std::valarray<std::complex<T>> &spectrum,
290  std::valarray<std::complex<T>> &signal,
291  int ncolumns, int nrows);
292 
295  void inverseAzimuthFFT(std::complex<T>* spectrum,
296  std::complex<T>* signal,
297  int ncolumns, int nrows);
300  void inverseAzimuthFFT(std::valarray<std::complex<T>> &spectrum,
301  std::valarray<T> &signal,
302  int ncolumns, int nrows);
303 
306  void inverseAzimuthFFT(std::complex<T>* spectrum,
307  T* signal,
308  int ncolumns, int nrows);
309 
312  void inverse2DFFT(std::valarray<std::complex<T>> &spectrum,
313  std::valarray<std::complex<T>> &signal,
314  int ncolumns, int nrows);
315 
318  void inverse2DFFT(std::complex<T>* spectrum,
319  std::complex<T>* signal,
320  int ncolumns, int nrows);
321 
324  void inverse2DFFT(std::valarray<std::complex<T>> &spectrum,
325  std::valarray<T> &signal,
326  int ncolumns, int nrows);
327 
330  void inverse2DFFT(std::complex<T>* spectrum,
331  T* signal,
332  int ncolumns, int nrows);
333 
334 
336  void upsample(std::valarray<std::complex<T>> &signal,
337  std::valarray<std::complex<T>> &signalOversampled,
338  int rows, int fft_size, int oversampleFactor);
339 
344  void upsample(std::valarray<std::complex<T>> &signal,
345  std::valarray<std::complex<T>> &signalOversampled,
346  int rows, int fft_size, int oversampleFactor,
347  std::valarray<std::complex<T>> shiftImpact);
348 
349 
350 
351 
353  void upsample2D(std::valarray<std::complex<T>> &signal,
354  std::valarray<std::complex<T>> &signalOversampled,
355  int oversampleFactor);
356 
361  void upsample2D(std::valarray<std::complex<T>> &signal,
362  std::valarray<std::complex<T>> &signalOversampled,
363  int oversampleFactor,
364  std::valarray<std::complex<T>> &shiftImpact);
365 
366 
368  void upsample2D(std::complex<T> *signal,
369  std::complex<T> *signalOversampled,
370  int oversampleFactor);
371 
376  void upsample2D(std::complex<T> *signal,
377  std::complex<T> *signalOversampled,
378  int oversampleFactor,
379  std::complex<T> *shiftImpact);
380 
381 
382 
383 
385  inline void nextPowerOfTwo(size_t N, size_t& fftLength);
386 
387 
389  inline void _fwd_configure(int rank, int* n, int howmany,
390  int* inembed, int istride, int idist,
391  int* onembed, int ostride, int odist);
392 
394  inline void _rev_configure(int rank, int* n, int howmany,
395  int* inembed, int istride, int idist,
396  int* onembed, int ostride, int odist);
397 
398 
400  inline void _fwd_configureRangeFFT(int ncolumns, int nrows);
401 
403  inline void _fwd_configureAzimuthFFT(int ncolumns, int nrows);
404 
406  inline void _fwd_configure2DFFT(int incolumns, int inrows, int oncolumns, int onrows);
407 
409  inline void _rev_configureRangeFFT(int ncolumns, int nrows);
410 
412  inline void _rev_configureAzimuthFFT(int ncolumns, int nrows);
413 
415  inline void _rev_configure2DFFT(int ncolumns, int nrows);
416 
417 
418  private:
419  int _fwd_rank;
420  int* _fwd_n;
421  int _fwd_howmany;
422  int* _fwd_inembed;
423  int _fwd_istride;
424  int _fwd_idist;
425  int* _fwd_onembed;
426  int _fwd_ostride;
427  int _fwd_odist;
428 
429  int _rev_rank;
430  int* _rev_n;
431  int _rev_howmany;
432  int* _rev_inembed;
433  int _rev_istride;
434  int _rev_idist;
435  int* _rev_onembed;
436  int _rev_ostride;
437  int _rev_odist;
438 
439  struct impl;
440  std::unique_ptr<impl, void(*)(impl*)> pimpl;
441 };
442 
443 #define ISCE_SIGNAL_SIGNAL_ICC
444 #include "Signal.icc"
445 #undef ISCE_SIGNAL_SIGNAL_ICC
void inverse2DFFT(std::valarray< std::complex< T >> &spectrum, std::valarray< std::complex< T >> &signal, int ncolumns, int nrows)
initiate plan for inverse two dimensional FFT for a block of data
Definition: Signal.cpp:871
Signal()
Default constructor.
Definition: Signal.cpp:20
A class to handle 2D FFT or 1D FFT in range or azimuth directions.
Definition: forward.h:10
void forwardAzimuthFFT(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &spectrum, int ncolumns, int nrows)
initiate plan for forward FFT in azimuth direction for a block of complex data.
Definition: Signal.cpp:459
void _fwd_configure(int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist)
save FFT plan parameters
Definition: Signal.icc:221
void inverse(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output)
inverse transform
Definition: Signal.cpp:341
void inverseAzimuthFFT(std::valarray< std::complex< T >> &spectrum, std::valarray< std::complex< T >> &signal, int ncolumns, int nrows)
initiate plan for inverse FFT in azimuth direction for a block of data
Definition: Signal.cpp:790
void _rev_configure(int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist)
save FFT plan parameters
Definition: Signal.icc:271
void inverseRangeFFT(std::valarray< std::complex< T >> &spectrum, std::valarray< std::complex< T >> &signal, int ncolumns, int nrows)
initiate plan for backward FFT in range direction for a block of data
Definition: Signal.cpp:733
void _fwd_configure2DFFT(int incolumns, int inrows, int oncolumns, int onrows)
determine the required parameters for setting 2D FFT plans
Definition: Signal.icc:150
void forward(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output)
forward transform
Definition: Signal.cpp:288
void _rev_configure2DFFT(int ncolumns, int nrows)
determine the required parameters for setting 2D FFT plans
Definition: Signal.icc:182
void _fwd_configureRangeFFT(int ncolumns, int nrows)
determine the required parameters for setting range FFT plans
Definition: Signal.icc:35
void upsample(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &signalOversampled, int rows, int fft_size, int oversampleFactor)
upsampling a block of data in range direction
Definition: Signal.cpp:952
void _rev_configureRangeFFT(int ncolumns, int nrows)
determine the required parameters for setting range FFT plans
Definition: Signal.icc:62
void fftPlanBackward(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output, int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist, int sign)
initiate iverse FFTW3 plan for a block of spectrum input parameters follow FFTW3 interface for fftw_p...
Definition: Signal.cpp:171
void forwardRangeFFT(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &spectrum, int ncolumns, int nrows)
initiate plan for forward FFT in range direction for a block of complex data.
Definition: Signal.cpp:381
void nextPowerOfTwo(size_t N, size_t &fftLength)
next power of two
Definition: Signal.icc:19
void forward2DFFT(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &spectrum, int ncolumns, int nrows)
initiate plan for forward two imensional FFT for a block of complex data
Definition: Signal.cpp:535
void upsample2D(std::valarray< std::complex< T >> &signal, std::valarray< std::complex< T >> &signalOversampled, int oversampleFactor)
2D upsampling a block of 2D data
Definition: Signal.cpp:1051
void _fwd_configureAzimuthFFT(int ncolumns, int nrows)
determine the required parameters for setting azimuth FFT plans
Definition: Signal.icc:92
void _rev_configureAzimuthFFT(int ncolumns, int nrows)
determine the required parameters for setting azimuth FFT plans
Definition: Signal.icc:120
void fftPlanForward(std::valarray< std::complex< T >> &input, std::valarray< std::complex< T >> &output, int rank, int *n, int howmany, int *inembed, int istride, int idist, int *onembed, int ostride, int odist, int sign)
initiate forward FFTW3 plan for a block of complex data input parameters follow FFTW3 interface for f...
Definition: Signal.cpp:46

Generated for ISCE3.0 by doxygen 1.8.5.