13 #include <isce3/core/LUT1d.h>
14 #include <isce3/io/forward.h>
55 std::valarray<std::complex<float>> &shiftImpact);
59 std::valarray<std::complex<float>> &secSlc,
60 std::valarray<std::complex<float>> geometryIfgram,
61 std::valarray<std::complex<float>> geometryIfgramConj,
62 std::valarray<std::complex<float>> &refSpectrum,
63 std::valarray<std::complex<float>> &secSpectrum,
64 std::valarray<double> &rangeFrequencies,
75 inline void prf(
double);
94 inline void beta(
double);
111 std::valarray<std::complex<float>> &secAvgSpectrum,
112 std::valarray<double> &rangeFrequencies,
113 size_t blockRowsData,
115 double &frequencyShift);
132 double _rangeSamplingFrequency;
135 double _rangeBandwidth;
138 double _rangePixelSpacing;
144 double _commonAzimuthBandwidth;
153 int _azimuthLooks = 1;
155 bool _doMultiLook =
false;
158 bool _doCommonAzimuthbandFilter =
false;
161 bool _doCommonRangebandFilter =
false;
164 bool _computeCoherence =
true;
167 size_t blockRows = 8192;
170 size_t oversample = 1;
176 #define ISCE_SIGNAL_CROSSMUL_ICC
177 #include "Crossmul.icc"
178 #undef ISCE_SIGNAL_CROSSMUL_ICC
void rangeLooks(int)
Set number of range looks.
Definition: Crossmul.icc:90
void rangeCommonBandFilter(std::valarray< std::complex< float >> &refSlc, std::valarray< std::complex< float >> &secSlc, std::valarray< std::complex< float >> geometryIfgram, std::valarray< std::complex< float >> geometryIfgramConj, std::valarray< std::complex< float >> &refSpectrum, std::valarray< std::complex< float >> &secSpectrum, std::valarray< double > &rangeFrequencies, isce3::signal::Filter< float > &rngFilter, size_t blockRows, size_t ncols)
Range common band filtering.
Definition: Crossmul.cpp:446
Intereferogram generation by cross-multiplication of reference and secondary SLCs.
Definition: Crossmul.h:20
void lookdownShiftImpact(size_t oversample, size_t fft_size, size_t blockRows, std::valarray< std::complex< float >> &shiftImpact)
Compute the frequency response due to a subpixel shift introduced by upsampling and downsampling...
Definition: Crossmul.cpp:391
void doCommonAzimuthbandFiltering(bool)
Set common azimuth band filtering flag.
Definition: Crossmul.icc:108
void getPeakIndex(std::valarray< float > data, size_t &peakIndex)
estimate the index of the maximum of a vector of data
Definition: Crossmul.icc:171
void rangeSamplingFrequency(double)
Set range sampling frequency.
Definition: Crossmul.icc:39
void azimuthLooks(int)
Set number of azimuth looks.
Definition: Crossmul.icc:99
void commonAzimuthBandwidth(double)
Set azimuth common bandwidth.
Definition: Crossmul.icc:73
void crossmul(isce3::io::Raster &referenceSLC, isce3::io::Raster &secondarySLC, isce3::io::Raster &rngOffset, isce3::io::Raster &interferogram, isce3::io::Raster &coherence)
Run crossmul.
Definition: Crossmul.cpp:87
void doCommonRangebandFiltering(bool)
Set common range band filtering flag.
Definition: Crossmul.icc:115
void rangeFrequencyShift(std::valarray< std::complex< float >> &refAvgSpectrum, std::valarray< std::complex< float >> &secAvgSpectrum, std::valarray< double > &rangeFrequencies, size_t blockRowsData, size_t fft_size, double &frequencyShift)
Compute the avergae frequency shift in range direction between two SLCs.
Definition: Crossmul.icc:128
void beta(double)
Set beta parameter for the azimuth common band filter.
Definition: Crossmul.icc:81
void rangeBandwidth(double)
Set the range bandwidth.
Definition: Crossmul.icc:47
void doppler(isce3::core::LUT1d< double >, isce3::core::LUT1d< double >)
Set doppler LUTs for reference and secondary SLCs.
Definition: Crossmul.icc:18
void prf(double)
Set pulse repetition frequency (PRF)
Definition: Crossmul.icc:31
void wavelength(double)
Set Wavelength.
Definition: Crossmul.icc:65
void rangePixelSpacing(double)
Range pixel spacing.
Definition: Crossmul.icc:56
Data structure meant to handle Raster I/O operations.
Definition: Raster.h:34