isce3 0.25.0
Loading...
Searching...
No Matches
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#include <isce3/core/EMatrix.h>
18
21template<class T>
23 public:
25 Signal();
26
28 Signal(int nthreads);
29
30 ~Signal() {};
31
35 void fftPlanForward(std::valarray<std::complex<T>> &input,
36 std::valarray<std::complex<T>> &output,
37 int rank, int* n, int howmany,
38 int* inembed, int istride, int idist,
39 int* onembed, int ostride, int odist, int sign);
40
44 void fftPlanForward(std::complex<T>* input,
45 std::complex<T>* output,
46 int rank, int* n, int howmany,
47 int* inembed, int istride, int idist,
48 int* onembed, int ostride, int odist, int sign);
49
53 void fftPlanForward(std::valarray<T> &input,
54 std::valarray<std::complex<T>> &output,
55 int rank, int* n, int howmany,
56 int* inembed, int istride, int idist,
57 int* onembed, int ostride, int odist);
58
62 void fftPlanForward(T* input,
63 std::complex<T>* output,
64 int rank, int* n, int howmany,
65 int* inembed, int istride, int idist,
66 int* onembed, int ostride, int odist);
67
71 void fftPlanBackward(std::valarray<std::complex<T>> &input,
72 std::valarray<std::complex<T>> &output,
73 int rank, int* n, int howmany,
74 int* inembed, int istride, int idist,
75 int* onembed, int ostride, int odist, int sign);
76
80 void fftPlanBackward(std::complex<T>* input,
81 std::complex<T>* output,
82 int rank, int* n, int howmany,
83 int* inembed, int istride, int idist,
84 int* onembed, int ostride, int odist, int sign);
85
90 void fftPlanBackward(std::valarray<std::complex<T>>& input,
91 std::valarray<T>& output,
92 int rank, int* n, int howmany,
93 int* inembed, int istride, int idist,
94 int* onembed, int ostride, int odist);
95
100 void fftPlanBackward(std::complex<T>* input,
101 T* output,
102 int rank, int* n, int howmany,
103 int* inembed, int istride, int idist,
104 int* onembed, int ostride, int odist);
105
107 void forward(std::valarray<std::complex<T>> &input,
108 std::valarray<std::complex<T>> &output);
109
111 void forward(std::complex<T>* input,
112 std::complex<T>* output);
113
115 void forward(std::valarray<T> &input,
116 std::valarray<std::complex<T>> &output);
117
119 void forward(T* input,
120 std::complex<T>* output);
121
122
123
125 void inverse(std::valarray<std::complex<T>> &input,
126 std::valarray<std::complex<T>> &output);
127
129 void inverse(std::complex<T>* input,
130 std::complex<T>* output);
131
133 void inverse(std::valarray<std::complex<T>>& input,
134 std::valarray<T>& output);
135
137 void inverse(std::complex<T>* input,
138 T* output);
139
145 void forwardRangeFFT(std::valarray<std::complex<T>>& signal,
146 std::valarray<std::complex<T>>& spectrum,
147 int ncolumns, int nrows);
148
154 void forwardRangeFFT(std::complex<T>* signal,
155 std::complex<T>* spectrum,
156 int ncolumns, int nrows);
157
163 void forwardRangeFFT(std::valarray<T>& signal,
164 std::valarray<std::complex<T>>& spectrum,
165 int ncolumns, int nrows);
166
172 void forwardRangeFFT(T* signal,
173 std::complex<T>* spectrum,
174 int ncolumns, int nrows);
175
181 void forwardAzimuthFFT(std::valarray<std::complex<T>> &signal,
182 std::valarray<std::complex<T>> &spectrum,
183 int ncolumns, int nrows);
184
190 void forwardAzimuthFFT(std::complex<T>* signal,
191 std::complex<T>* spectrum,
192 int ncolumns, int nrows);
193
199 void forwardAzimuthFFT(std::valarray<T> &signal,
200 std::valarray<std::complex<T>> &spectrum,
201 int ncolumns, int nrows);
202
208 void forwardAzimuthFFT(T* signal,
209 std::complex<T>* spectrum,
210 int ncolumns, int nrows);
211
214 void forward2DFFT(std::valarray<std::complex<T>> &signal,
215 std::valarray<std::complex<T>> &spectrum,
216 int ncolumns, int nrows);
217
220 void forward2DFFT(std::valarray<std::complex<T>> &signal,
221 std::valarray<std::complex<T>> &spectrum,
222 int incolumns, int inrows,
223 int oncolumns, int onrows);
224
227 void forward2DFFT(std::complex<T>* signal,
228 std::complex<T>* spectrum,
229 int ncolumns, int nrows);
230
233 void forward2DFFT(std::complex<T>* signal,
234 std::complex<T>* spectrum,
235 int incolumns, int inrows,
236 int oncolumns, int onrows);
237
240 void forward2DFFT(std::valarray<T>& signal,
241 std::valarray<std::complex<T>>& spectrum,
242 int ncolumns, int nrows);
243
246 void forward2DFFT(std::valarray<T>& signal,
247 std::valarray<std::complex<T>>& spectrum,
248 int incolumns, int inrows,
249 int oncolumns, int onrows);
250
253 void forward2DFFT(T* signal,
254 std::complex<T>* spectrum,
255 int ncolumns, int nrows);
256
259 void forward2DFFT(T* signal,
260 std::complex<T>* spectrum,
261 int incolumns, int inrows,
262 int oncolumns, int onrows);
263
266 void inverseRangeFFT(std::valarray<std::complex<T>> &spectrum,
267 std::valarray<std::complex<T>> &signal,
268 int ncolumns, int nrows);
269
272 void inverseRangeFFT(std::complex<T>* spectrum,
273 std::complex<T>* signal,
274 int ncolumns, int nrows);
275
278 void inverseRangeFFT(std::valarray<std::complex<T>> &spectrum,
279 std::valarray<T> &signal,
280 int ncolumns, int nrows);
281
284 void inverseRangeFFT(std::complex<T>* spectrum,
285 T* signal,
286 int ncolumns, int nrows);
287
290 void inverseAzimuthFFT(std::valarray<std::complex<T>> &spectrum,
291 std::valarray<std::complex<T>> &signal,
292 int ncolumns, int nrows);
293
296 void inverseAzimuthFFT(std::complex<T>* spectrum,
297 std::complex<T>* signal,
298 int ncolumns, int nrows);
301 void inverseAzimuthFFT(std::valarray<std::complex<T>> &spectrum,
302 std::valarray<T> &signal,
303 int ncolumns, int nrows);
304
307 void inverseAzimuthFFT(std::complex<T>* spectrum,
308 T* signal,
309 int ncolumns, int nrows);
310
313 void inverse2DFFT(std::valarray<std::complex<T>> &spectrum,
314 std::valarray<std::complex<T>> &signal,
315 int ncolumns, int nrows);
316
319 void inverse2DFFT(std::complex<T>* spectrum,
320 std::complex<T>* signal,
321 int ncolumns, int nrows);
322
325 void inverse2DFFT(std::valarray<std::complex<T>> &spectrum,
326 std::valarray<T> &signal,
327 int ncolumns, int nrows);
328
331 void inverse2DFFT(std::complex<T>* spectrum,
332 T* signal,
333 int ncolumns, int nrows);
334
335
337 void upsample(std::valarray<std::complex<T>> &signal,
338 std::valarray<std::complex<T>> &signalOversampled,
339 int rows, int fft_size, int oversampleFactor);
340
345 void upsample(std::valarray<std::complex<T>> &signal,
346 std::valarray<std::complex<T>> &signalOversampled,
347 int rows, int fft_size, int oversampleFactor,
348 std::valarray<std::complex<T>> shiftImpact);
349
354 void upsample(
355 isce3::core::EArray2D<std::complex<T>>& signal,
356 isce3::core::EArray2D<std::complex<T>>& signalUpsampled,
357 const isce3::core::EArray2D<std::complex<T>>& shiftImpact = {});
358
360 void upsample2D(std::valarray<std::complex<T>> &signal,
361 std::valarray<std::complex<T>> &signalOversampled,
362 int oversampleFactor);
363
368 void upsample2D(std::valarray<std::complex<T>> &signal,
369 std::valarray<std::complex<T>> &signalOversampled,
370 int oversampleFactor,
371 std::valarray<std::complex<T>> &shiftImpact);
372
373
375 void upsample2D(std::complex<T> *signal,
376 std::complex<T> *signalOversampled,
377 int oversampleFactor);
378
383 void upsample2D(std::complex<T> *signal,
384 std::complex<T> *signalOversampled,
385 int oversampleFactor,
386 std::complex<T> *shiftImpact);
387
388
389
390
392 inline void nextPowerOfTwo(size_t N, size_t& fftLength);
393
394
396 inline void _fwd_configure(int rank, int* n, int howmany,
397 int* inembed, int istride, int idist,
398 int* onembed, int ostride, int odist);
399
401 inline void _rev_configure(int rank, int* n, int howmany,
402 int* inembed, int istride, int idist,
403 int* onembed, int ostride, int odist);
404
405
407 inline void _fwd_configureRangeFFT(int ncolumns, int nrows);
408
410 inline void _fwd_configureAzimuthFFT(int ncolumns, int nrows);
411
413 inline void _fwd_configure2DFFT(int incolumns, int inrows, int oncolumns, int onrows);
414
416 inline void _rev_configureRangeFFT(int ncolumns, int nrows);
417
419 inline void _rev_configureAzimuthFFT(int ncolumns, int nrows);
420
422 inline void _rev_configure2DFFT(int ncolumns, int nrows);
423
424
425 private:
426 int _fwd_rank;
427 int* _fwd_n;
428 int _fwd_howmany;
429 int* _fwd_inembed;
430 int _fwd_istride;
431 int _fwd_idist;
432 int* _fwd_onembed;
433 int _fwd_ostride;
434 int _fwd_odist;
435
436 int _rev_rank;
437 int* _rev_n;
438 int _rev_howmany;
439 int* _rev_inembed;
440 int _rev_istride;
441 int _rev_idist;
442 int* _rev_onembed;
443 int _rev_ostride;
444 int _rev_odist;
445
446 struct impl;
447 std::unique_ptr<impl, void(*)(impl*)> pimpl;
448};
449
450#define ISCE_SIGNAL_SIGNAL_ICC
451#include "Signal.icc"
452#undef ISCE_SIGNAL_SIGNAL_ICC
A class to handle 2D FFT or 1D FFT in range or azimuth directions.
Definition Signal.h:22
void inverse(std::valarray< std::complex< T > > &input, std::valarray< std::complex< T > > &output)
inverse transform
Definition Signal.cpp:334
void _rev_configureAzimuthFFT(int ncolumns, int nrows)
determine the required parameters for setting azimuth FFT plans
Definition Signal.icc:120
void _rev_configureRangeFFT(int ncolumns, int nrows)
determine the required parameters for setting range FFT plans
Definition Signal.icc:62
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:39
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:945
void _fwd_configureAzimuthFFT(int ncolumns, int nrows)
determine the required parameters for setting azimuth FFT plans
Definition Signal.icc:92
void _fwd_configureRangeFFT(int ncolumns, int nrows)
determine the required parameters for setting range FFT plans
Definition Signal.icc:35
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:864
Signal()
Default constructor.
Definition Signal.cpp:13
void _rev_configure2DFFT(int ncolumns, int nrows)
determine the required parameters for setting 2D FFT plans
Definition Signal.icc:182
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:164
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 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:528
void nextPowerOfTwo(size_t N, size_t &fftLength)
next power of two
Definition Signal.icc:19
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 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:783
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:452
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:726
void forward(std::valarray< std::complex< T > > &input, std::valarray< std::complex< T > > &output)
forward transform
Definition Signal.cpp:281
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:374
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:1097
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

Generated for ISCE3.0 by doxygen 1.13.2.