isce3 0.25.0
Loading...
Searching...
No Matches
ResampSlc.icc
1// Since fmod(a,b) in C++ != MODULO(a,b) in Fortran for all a,b, define a C++
2// equivalent
3#define modulo_f(a, b) fmod(fmod(a, b) + (b), (b))
4
5namespace isce3 { namespace image {
6
7// Constructor from an isce3::product::RadarGridProduct
9 char frequency, const std::complex<float> invalid_value)
10 : ResampSlc(product.swath(frequency), invalid_value)
11{
12 // Set the doppler
13 _dopplerLUT = product.metadata().procInfo().dopplerCentroid(frequency);
14 // Save the filename
15 _filename = product.filename();
16}
17
18// Constructor from an isce3::product::RadarGridProduct and reference product
19// (flattening)
21 const isce3::product::RadarGridProduct& refProduct,
22 char frequency, const std::complex<float> invalid_value)
23 : ResampSlc(product.swath(frequency), refProduct.swath(frequency),
24 invalid_value)
25{
26 // Set the doppler
27 _dopplerLUT = product.metadata().procInfo().dopplerCentroid(frequency);
28 // Save the filename
29 _filename = product.filename();
30}
31
32// Constructor from an isce3::product::Swath (no flattening)
34 const std::complex<float> invalid_value)
35 : _haveRefData(false), _invalid_value(invalid_value)
36{
37 // Set radar parameters
38 _setDataFromSwath(swath);
39}
40
41// Constructor from an isce3::product::Swath and reference swath (flattening)
43 const isce3::product::Swath& refSwath,
44 const std::complex<float> invalid_value)
45 : _haveRefData(true), _invalid_value(invalid_value)
46{
47 // Set radar parameters
48 _setDataFromSwath(swath);
49 // Set reference radar parameters
50 _setRefDataFromSwath(refSwath);
51}
52
53// Constructor from an isce3::product::RadarGridParameters (no flattening)
56 const std::complex<float> invalid_value)
57 : _haveRefData(false), _dopplerLUT(doppler),
58 _startingRange(rdr_grid.startingRange()),
59 _rangePixelSpacing(rdr_grid.rangePixelSpacing()),
60 _sensingStart(rdr_grid.sensingStart()), _prf(rdr_grid.prf()),
61 _wavelength(rdr_grid.wavelength()), _invalid_value(invalid_value)
62{}
63
64// Constructor from an isce3::product::RadarGridParameters and reference radar
65// grid (flattening)
68 const isce3::product::RadarGridParameters& ref_rdr_grid,
70 const std::complex<float> invalid_value)
71 : _haveRefData(true), _dopplerLUT(doppler),
72 _startingRange(rdr_grid.startingRange()),
73 _rangePixelSpacing(rdr_grid.rangePixelSpacing()),
74 _sensingStart(rdr_grid.sensingStart()), _prf(rdr_grid.prf()),
75 _wavelength(rdr_grid.wavelength()),
76 _refStartingRange(ref_rdr_grid.startingRange()),
77 _refRangePixelSpacing(ref_rdr_grid.rangePixelSpacing()),
78 _refWavelength(ref_rdr_grid.wavelength()), _invalid_value(invalid_value)
79{}
80
81// Constructor from individual components (no flattening)
83 double startingRange, double rangePixelSpacing,
84 double sensingStart, double prf, double wvl,
85 const std::complex<float> invalid_value)
86 : _haveRefData(false), _dopplerLUT(doppler), _startingRange(startingRange),
87 _rangePixelSpacing(rangePixelSpacing), _sensingStart(sensingStart),
88 _prf(prf), _wavelength(wvl), _invalid_value(invalid_value)
89{}
90
91// Constructor from individual components (flattening)
93 double startingRange, double rangePixelSpacing,
94 double sensingStart, double prf, double wvl,
95 double refStartingRange,
96 double refRangePixelSpacing, double refWvl,
97 const std::complex<float> invalid_value)
98 : _haveRefData(true), _dopplerLUT(doppler), _startingRange(startingRange),
99 _rangePixelSpacing(rangePixelSpacing), _sensingStart(sensingStart),
100 _prf(prf), _wavelength(wvl), _refStartingRange(refStartingRange),
101 _refRangePixelSpacing(refRangePixelSpacing), _refWavelength(refWvl),
102 _invalid_value(invalid_value)
103{}
104
105// Set radar parameters from an isce3::product::Swath
106inline void ResampSlc::_setDataFromSwath(const isce3::product::Swath& swath)
107{
108 _startingRange = swath.slantRange()[0];
109 _rangePixelSpacing = swath.rangePixelSpacing();
110 _sensingStart = swath.zeroDopplerTime()[0];
111 _prf = 1.0 / swath.zeroDopplerTimeSpacing();
112 _wavelength = swath.processedWavelength();
113}
114
115// Set reference radar parameters from an isce3::product::Swath (for flattening)
116inline void ResampSlc::_setRefDataFromSwath(const isce3::product::Swath& swath)
117{
118 _refStartingRange = swath.slantRange()[0];
119 _refRangePixelSpacing = swath.rangePixelSpacing();
120 _refWavelength = swath.processedWavelength();
121}
122
123// Get range carrier polynomial
124inline isce3::core::Poly2d ResampSlc::rgCarrier() const { return _rgCarrier; }
125
126// Set range carrier polynomial
127inline void ResampSlc::rgCarrier(const isce3::core::Poly2d& lut)
128{
129 _rgCarrier = lut;
130}
131
132// Get azimuth carrier polynomial
133inline isce3::core::Poly2d ResampSlc::azCarrier() const { return _azCarrier; }
134
135// Set azimuth carrier polynomial
136inline void ResampSlc::azCarrier(const isce3::core::Poly2d& lut)
137{
138 _azCarrier = lut;
139}
140
141// Get read-only reference to Doppler LUT
143{
144 return _dopplerLUT;
145}
146
147// Get reference to Doppler LUT
148inline isce3::core::LUT2d<double>& ResampSlc::doppler() { return _dopplerLUT; }
149
150// Set Doppler LUT
152{
153 _dopplerLUT = lut;
154}
155
156// Set reference product
157inline void
158ResampSlc::referenceProduct(const isce3::product::RadarGridProduct& refProduct,
159 char frequency)
160{
161 _setRefDataFromSwath(refProduct.swath(frequency));
162 _haveRefData = true;
163}
164
165// Announce my properties to the world
166inline void ResampSlc::declare(size_t inLength, size_t inWidth, size_t outLength,
167 size_t outWidth) const
168{
169 // Make info channel
170 pyre::journal::info_t channel("isce.core.ResampSlc");
171 // Basic info
172 channel << pyre::journal::newline
173 << "Resample one image to another image coordinates >>"
174 << pyre::journal::newline << pyre::journal::newline;
175 channel << "Input Image Dimensions: " << inLength << " lines, " << inWidth
176 << " pixels" << pyre::journal::newline;
177 channel << "Output Image Dimensions: " << outLength << " lines, "
178 << outWidth << " pixels" << pyre::journal::newline
179 << pyre::journal::newline;
180 channel << "Complex data interpolation" << pyre::journal::newline;
181}
182
183// Get the number of lines per tile
184inline size_t ResampSlc::linesPerTile() const { return _linesPerTile; }
185
186// Set the number of lines per tile
187inline void ResampSlc::linesPerTile(size_t value) { _linesPerTile = value; }
188
189// Compute number of tiles given a specified nominal tile size
190inline size_t ResampSlc::_computeNumberOfTiles(size_t outLength, size_t linesPerTile)
191{
192 // Compute floor(nTiles)
193 size_t nTiles = outLength / linesPerTile;
194 // See if there are any leftover lines to add an extra tile
195 size_t extraLines = outLength - nTiles * linesPerTile;
196 if (extraLines > 0) {
197 nTiles += 1;
198 }
199 return nTiles;
200}
201
202// Prepare interpolation pointer
203inline void ResampSlc::_prepareInterpMethods(isce3::core::dataInterpMethod,
204 int sinc_len)
205{
206 _interp = new isce3::core::Sinc2dInterpolator<std::complex<float>>(
207 sinc_len, isce3::core::SINC_SUB);
208}
209
210}} // namespace isce3::image
Data structure to store 2D Lookup table.
Definition LUT2d.h:20
Data structure for representing 2D polynomials.
Definition Poly2d.h:20
ResampSlc(const isce3::product::RadarGridProduct &product, char frequency='A', const std::complex< float > invalid_value=std::complex< float >(0.0))
Constructor from an isce3::product::RadarGridProduct (no flattening)
Definition ResampSlc.icc:8
const isce3::core::LUT2d< double > & doppler() const
Get read-only reference to Doppler LUT2d.
Definition ResampSlc.icc:142
Definition RadarGridParameters.h:16
RadarGridProduct class declaration.
Definition RadarGridProduct.h:71
const Swath & swath(char freq) const
Get a read-only reference to a swath.
Definition RadarGridProduct.h:85
Definition Swath.h:23
double rangePixelSpacing() const
Get the range pixel spacing.
Definition Swath.h:38
const std::valarray< double > & zeroDopplerTime() const
Get zero Doppler time array.
Definition Swath.h:41
double processedWavelength() const
Get processed wavelength.
Definition Swath.h:65
double zeroDopplerTimeSpacing() const
Get time spacing of raster grid.
Definition Swath.h:85
const std::valarray< double > & slantRange() const
Get slant range array.
Definition Swath.h:30
The isce3::product namespace.
Definition forward.h:3
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5

Generated for ISCE3.0 by doxygen 1.13.2.