isce3  0.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
Crossmul.h
1 // -*- C++ -*-
2 // -*- coding: utf-8 -*-
3 //
4 // Author: Heresh Fattahi, Bryan Riel
5 // Copyright 2018-
6 //
7 
8 #pragma once
9 
10 #include "forward.h"
11 
12 #include <complex>
13 #include <isce3/core/LUT1d.h>
14 #include <isce3/io/forward.h>
15 
21  public:
22  // Constructor from product
23  Crossmul() {};
24 
25  ~Crossmul() {};
26 
27  /*
28  void Crossmul(const isce3::product::Product& referenceSLC,
29  const isce3::product::Product& secondarySLC,
30  const isce3::product::Product& outputInterferogram);
31  */
32 
33 
35  void crossmul(isce3::io::Raster& referenceSLC,
36  isce3::io::Raster& secondarySLC,
37  isce3::io::Raster& rngOffset,
38  isce3::io::Raster& interferogram,
39  isce3::io::Raster& coherence);
40 
42  void crossmul(isce3::io::Raster& referenceSLC,
43  isce3::io::Raster& secondarySLC,
44  isce3::io::Raster& interferogram,
45  isce3::io::Raster& coherence);
46 
48  void crossmul(isce3::io::Raster& referenceSLC,
49  isce3::io::Raster& secondarySLC,
50  isce3::io::Raster& interferogram);
51 
53  void lookdownShiftImpact(size_t oversample, size_t fft_size,
54  size_t blockRows,
55  std::valarray<std::complex<float>> &shiftImpact);
56 
58  void rangeCommonBandFilter(std::valarray<std::complex<float>> &refSlc,
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,
66  size_t blockRows,
67  size_t ncols);
68 
69 
73 
75  inline void prf(double);
76 
78  inline void rangeSamplingFrequency(double);
79 
81  inline void rangeBandwidth(double);
82 
84  inline void rangePixelSpacing(double);
85 
87  inline void wavelength(double);
88 
89 
91  inline void commonAzimuthBandwidth(double);
92 
94  inline void beta(double);
95 
96 
98  inline void rangeLooks(int);
99 
101  inline void azimuthLooks(int);
102 
104  inline void doCommonAzimuthbandFiltering(bool);
105 
107  inline void doCommonRangebandFiltering(bool);
108 
110  inline void rangeFrequencyShift(std::valarray<std::complex<float>> &refAvgSpectrum,
111  std::valarray<std::complex<float>> &secAvgSpectrum,
112  std::valarray<double> &rangeFrequencies,
113  size_t blockRowsData,
114  size_t fft_size,
115  double &frequencyShift);
116 
118  inline void getPeakIndex(std::valarray<float> data,
119  size_t &peakIndex);
120 
121  private:
122  //Doppler LUT for the refernce SLC
123  isce3::core::LUT1d<double> _refDoppler;
124 
125  //Doppler LUT for the secondary SLC
126  isce3::core::LUT1d<double> _secDoppler;
127 
128  //pulse repetition frequency
129  double _prf;
130 
131  // range samping frequency
132  double _rangeSamplingFrequency;
133 
134  // range signal bandwidth
135  double _rangeBandwidth;
136 
137  // range pixel spacing
138  double _rangePixelSpacing;
139 
140  // radar wavelength
141  double _wavelength;
142 
143  //azimuth common bandwidth
144  double _commonAzimuthBandwidth;
145 
146  // beta parameter for constructing common azimuth band filter
147  double _beta;
148 
149  // number of range looks
150  int _rangeLooks = 1;
151 
152  // number of azimuth looks
153  int _azimuthLooks = 1;
154 
155  bool _doMultiLook = false;
156 
157  // Flag for common azimuth band filtering
158  bool _doCommonAzimuthbandFilter = false;
159 
160  // Flag for common range band filtering
161  bool _doCommonRangebandFilter = false;
162 
163  // Flag for computing coherence
164  bool _computeCoherence = true;
165 
166  // number of lines per block
167  size_t blockRows = 8192;
168 
169  // upsampling factor
170  size_t oversample = 1;
171 
172 
173 };
174 
175 // Get inline implementations for Crossmul
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
Definition: Filter.h:27
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

Generated for ISCE3.0 by doxygen 1.8.5.