QGIS API Documentation 3.41.0-Master (af5edcb665c)
Loading...
Searching...
No Matches
qgsleastsquares.cpp
Go to the documentation of this file.
1/***************************************************************************
2 qgsleastsquares.cpp
3 --------------------------------------
4 Date : Sun Sep 16 12:03:37 AKDT 2007
5 Copyright : (C) 2007 by Gary E. Sherman
6 Email : sherman at mrcc dot com
7 ***************************************************************************
8 * *
9 * This program is free software; you can redistribute it and/or modify *
10 * it under the terms of the GNU General Public License as published by *
11 * the Free Software Foundation; either version 2 of the License, or *
12 * (at your option) any later version. *
13 * *
14 ***************************************************************************/
15
16#include "qgsleastsquares.h"
17#include "qgsconfig.h"
18#include "qgsexception.h"
19
20#include <QObject>
21
22#include <cmath>
23#include <stdexcept>
24
25#ifdef HAVE_GSL
26#include <gsl/gsl_linalg.h>
27#include <gsl/gsl_blas.h>
28#endif
29
30void QgsLeastSquares::linear( const QVector<QgsPointXY> &sourceCoordinates, const QVector<QgsPointXY> &destinationCoordinates, QgsPointXY &origin, double &pixelXSize, double &pixelYSize )
31{
32 const int n = destinationCoordinates.size();
33 if ( n < 2 )
34 {
35 throw std::domain_error( QObject::tr( "Fit to a linear transform requires at least 2 points." ).toLocal8Bit().constData() );
36 }
37
38 double sumPx( 0 ), sumPy( 0 ), sumPx2( 0 ), sumPy2( 0 ), sumPxMx( 0 ), sumPyMy( 0 ), sumMx( 0 ), sumMy( 0 );
39 for ( int i = 0; i < n; ++i )
40 {
41 sumPx += sourceCoordinates.at( i ).x();
42 sumPy += sourceCoordinates.at( i ).y();
43 sumPx2 += std::pow( sourceCoordinates.at( i ).x(), 2 );
44 sumPy2 += std::pow( sourceCoordinates.at( i ).y(), 2 );
45 sumPxMx += sourceCoordinates.at( i ).x() * destinationCoordinates.at( i ).x();
46 sumPyMy += sourceCoordinates.at( i ).y() * destinationCoordinates.at( i ).y();
47 sumMx += destinationCoordinates.at( i ).x();
48 sumMy += destinationCoordinates.at( i ).y();
49 }
50
51 const double deltaX = n * sumPx2 - std::pow( sumPx, 2 );
52 const double deltaY = n * sumPy2 - std::pow( sumPy, 2 );
53
54 const double aX = ( sumPx2 * sumMx - sumPx * sumPxMx ) / deltaX;
55 const double aY = ( sumPy2 * sumMy - sumPy * sumPyMy ) / deltaY;
56 const double bX = ( n * sumPxMx - sumPx * sumMx ) / deltaX;
57 const double bY = ( n * sumPyMy - sumPy * sumMy ) / deltaY;
58
59 origin.setX( aX );
60 origin.setY( aY );
61
62 pixelXSize = std::fabs( bX );
63 pixelYSize = std::fabs( bY );
64}
65
66
67void QgsLeastSquares::helmert( const QVector<QgsPointXY> &sourceCoordinates, const QVector<QgsPointXY> &destinationCoordinates, QgsPointXY &origin, double &pixelSize, double &rotation )
68{
69#ifndef HAVE_GSL
70 ( void ) sourceCoordinates;
71 ( void ) destinationCoordinates;
72 ( void ) origin;
73 ( void ) pixelSize;
74 ( void ) rotation;
75 throw QgsNotSupportedException( QObject::tr( "Calculating a helmert transformation requires a QGIS build based GSL" ) );
76#else
77 const int n = destinationCoordinates.size();
78 if ( n < 2 )
79 {
80 throw std::domain_error( QObject::tr( "Fit to a Helmert transform requires at least 2 points." ).toLocal8Bit().constData() );
81 }
82
83 double A = 0;
84 double B = 0;
85 double C = 0;
86 double D = 0;
87 double E = 0;
88 double F = 0;
89 double G = 0;
90 double H = 0;
91 double I = 0;
92 double J = 0;
93 for ( int i = 0; i < n; ++i )
94 {
95 A += sourceCoordinates.at( i ).x();
96 B += sourceCoordinates.at( i ).y();
97 C += destinationCoordinates.at( i ).x();
98 D += destinationCoordinates.at( i ).y();
99 E += destinationCoordinates.at( i ).x() * sourceCoordinates.at( i ).x();
100 F += destinationCoordinates.at( i ).y() * sourceCoordinates.at( i ).y();
101 G += std::pow( sourceCoordinates.at( i ).x(), 2 );
102 H += std::pow( sourceCoordinates.at( i ).y(), 2 );
103 I += destinationCoordinates.at( i ).x() * sourceCoordinates.at( i ).y();
104 J += sourceCoordinates.at( i ).x() * destinationCoordinates.at( i ).y();
105 }
106
107 /* The least squares fit for the parameters { a, b, x0, y0 } is the solution
108 to the matrix equation Mx = b, where M and b is given below. I *think*
109 that this is correct but I derived it myself late at night. Look at
110 helmert.jpg if you suspect bugs. */
111
112 double MData[] = { A, -B, ( double ) n, 0., B, A, 0., ( double ) n, G + H, 0., A, B, 0., G + H, -B, A };
113
114 double bData[] = { C, D, E + F, J - I };
115
116 // we want to solve the equation M*x = b, where x = [a b x0 y0]
117 gsl_matrix_view M = gsl_matrix_view_array( MData, 4, 4 );
118 const gsl_vector_view b = gsl_vector_view_array( bData, 4 );
119 gsl_vector *x = gsl_vector_alloc( 4 );
120 gsl_permutation *p = gsl_permutation_alloc( 4 );
121 int s;
122 gsl_linalg_LU_decomp( &M.matrix, p, &s );
123 gsl_linalg_LU_solve( &M.matrix, p, &b.vector, x );
124 gsl_permutation_free( p );
125
126 origin.setX( gsl_vector_get( x, 2 ) );
127 origin.setY( gsl_vector_get( x, 3 ) );
128 pixelSize = std::sqrt( std::pow( gsl_vector_get( x, 0 ), 2 ) + std::pow( gsl_vector_get( x, 1 ), 2 ) );
129 rotation = std::atan2( gsl_vector_get( x, 1 ), gsl_vector_get( x, 0 ) );
130
131 gsl_vector_free( x );
132#endif
133}
134
135#if 0
136void QgsLeastSquares::affine( QVector<QgsPointXY> mapCoords,
137 QVector<QgsPointXY> pixelCoords )
138{
139 int n = mapCoords.size();
140 if ( n < 4 )
141 {
142 throw std::domain_error( QObject::tr( "Fit to an affine transform requires at least 4 points." ).toLocal8Bit().constData() );
143 }
144
145 double A = 0, B = 0, C = 0, D = 0, E = 0, F = 0,
146 G = 0, H = 0, I = 0, J = 0, K = 0;
147 for ( int i = 0; i < n; ++i )
148 {
149 A += pixelCoords[i].x();
150 B += pixelCoords[i].y();
151 C += mapCoords[i].x();
152 D += mapCoords[i].y();
153 E += std::pow( pixelCoords[i].x(), 2 );
154 F += std::pow( pixelCoords[i].y(), 2 );
155 G += pixelCoords[i].x() * pixelCoords[i].y();
156 H += pixelCoords[i].x() * mapCoords[i].x();
157 I += pixelCoords[i].y() * mapCoords[i].y();
158 J += pixelCoords[i].x() * mapCoords[i].y();
159 K += mapCoords[i].x() * pixelCoords[i].y();
160 }
161
162 /* The least squares fit for the parameters { a, b, c, d, x0, y0 } is the
163 solution to the matrix equation Mx = b, where M and b is given below.
164 I *think* that this is correct but I derived it myself late at night.
165 Look at affine.jpg if you suspect bugs. */
166
167 double MData[] = { A, B, 0, 0, ( double ) n, 0,
168 0, 0, A, B, 0, ( double ) n,
169 E, G, 0, 0, A, 0,
170 G, F, 0, 0, B, 0,
171 0, 0, E, G, 0, A,
172 0, 0, G, F, 0, B
173 };
174
175 double bData[] = { C, D, H, K, J, I };
176
177 // we want to solve the equation M*x = b, where x = [a b c d x0 y0]
178 gsl_matrix_view M = gsl_matrix_view_array( MData, 6, 6 );
179 gsl_vector_view b = gsl_vector_view_array( bData, 6 );
180 gsl_vector *x = gsl_vector_alloc( 6 );
181 gsl_permutation *p = gsl_permutation_alloc( 6 );
182 int s;
183 gsl_linalg_LU_decomp( &M.matrix, p, &s );
184 gsl_linalg_LU_solve( &M.matrix, p, &b.vector, x );
185 gsl_permutation_free( p );
186
187}
188#endif
189
195void normalizeCoordinates( const QVector<QgsPointXY> &coords, QVector<QgsPointXY> &normalizedCoords, double normalizeMatrix[9], double denormalizeMatrix[9] )
196{
197 // Calculate center of gravity
198 double cogX = 0.0, cogY = 0.0;
199 for ( int i = 0; i < coords.size(); i++ )
200 {
201 cogX += coords[i].x();
202 cogY += coords[i].y();
203 }
204 cogX *= 1.0 / coords.size();
205 cogY *= 1.0 / coords.size();
206
207 // Calculate mean distance to origin
208 double meanDist = 0.0;
209 for ( int i = 0; i < coords.size(); i++ )
210 {
211 const double X = ( coords[i].x() - cogX );
212 const double Y = ( coords[i].y() - cogY );
213 meanDist += std::sqrt( X * X + Y * Y );
214 }
215 meanDist *= 1.0 / coords.size();
216
217 const double OOD = meanDist * M_SQRT1_2;
218 const double D = 1.0 / OOD;
219 normalizedCoords.resize( coords.size() );
220 for ( int i = 0; i < coords.size(); i++ )
221 {
222 normalizedCoords[i] = QgsPointXY( ( coords[i].x() - cogX ) * D, ( coords[i].y() - cogY ) * D );
223 }
224
225 normalizeMatrix[0] = D;
226 normalizeMatrix[1] = 0.0;
227 normalizeMatrix[2] = -cogX * D;
228 normalizeMatrix[3] = 0.0;
229 normalizeMatrix[4] = D;
230 normalizeMatrix[5] = -cogY * D;
231 normalizeMatrix[6] = 0.0;
232 normalizeMatrix[7] = 0.0;
233 normalizeMatrix[8] = 1.0;
234
235 denormalizeMatrix[0] = OOD;
236 denormalizeMatrix[1] = 0.0;
237 denormalizeMatrix[2] = cogX;
238 denormalizeMatrix[3] = 0.0;
239 denormalizeMatrix[4] = OOD;
240 denormalizeMatrix[5] = cogY;
241 denormalizeMatrix[6] = 0.0;
242 denormalizeMatrix[7] = 0.0;
243 denormalizeMatrix[8] = 1.0;
244}
245
246// Fits a homography to the given corresponding points, and
247// return it in H (row-major format).
248void QgsLeastSquares::projective( const QVector<QgsPointXY> &sourceCoordinates, const QVector<QgsPointXY> &destinationCoordinates, double H[9] )
249{
250#ifndef HAVE_GSL
251 ( void ) sourceCoordinates;
252 ( void ) destinationCoordinates;
253 ( void ) H;
254 throw QgsNotSupportedException( QObject::tr( "Calculating a projective transformation requires a QGIS build based GSL" ) );
255#else
256 Q_ASSERT( sourceCoordinates.size() == destinationCoordinates.size() );
257
258 if ( destinationCoordinates.size() < 4 )
259 {
260 throw std::domain_error( QObject::tr( "Fitting a projective transform requires at least 4 corresponding points." ).toLocal8Bit().constData() );
261 }
262
263 QVector<QgsPointXY> sourceCoordinatesNormalized;
264 QVector<QgsPointXY> destinationCoordinatesNormalized;
265
266 double normSource[9], denormSource[9];
267 double normDest[9], denormDest[9];
268 normalizeCoordinates( sourceCoordinates, sourceCoordinatesNormalized, normSource, denormSource );
269 normalizeCoordinates( destinationCoordinates, destinationCoordinatesNormalized, normDest, denormDest );
270
271 // GSL does not support a full SVD, so we artificially add a linear dependent row
272 // to the matrix in case the system is underconstrained.
273 const uint m = std::max( 9u, ( uint ) destinationCoordinatesNormalized.size() * 2u );
274 const uint n = 9;
275 gsl_matrix *S = gsl_matrix_alloc( m, n );
276
277 for ( int i = 0; i < destinationCoordinatesNormalized.size(); i++ )
278 {
279 gsl_matrix_set( S, i * 2, 0, sourceCoordinatesNormalized[i].x() );
280 gsl_matrix_set( S, i * 2, 1, sourceCoordinatesNormalized[i].y() );
281 gsl_matrix_set( S, i * 2, 2, 1.0 );
282
283 gsl_matrix_set( S, i * 2, 3, 0.0 );
284 gsl_matrix_set( S, i * 2, 4, 0.0 );
285 gsl_matrix_set( S, i * 2, 5, 0.0 );
286
287 gsl_matrix_set( S, i * 2, 6, -destinationCoordinatesNormalized[i].x() * sourceCoordinatesNormalized[i].x() );
288 gsl_matrix_set( S, i * 2, 7, -destinationCoordinatesNormalized[i].x() * sourceCoordinatesNormalized[i].y() );
289 gsl_matrix_set( S, i * 2, 8, -destinationCoordinatesNormalized[i].x() * 1.0 );
290
291 gsl_matrix_set( S, i * 2 + 1, 0, 0.0 );
292 gsl_matrix_set( S, i * 2 + 1, 1, 0.0 );
293 gsl_matrix_set( S, i * 2 + 1, 2, 0.0 );
294
295 gsl_matrix_set( S, i * 2 + 1, 3, sourceCoordinatesNormalized[i].x() );
296 gsl_matrix_set( S, i * 2 + 1, 4, sourceCoordinatesNormalized[i].y() );
297 gsl_matrix_set( S, i * 2 + 1, 5, 1.0 );
298
299 gsl_matrix_set( S, i * 2 + 1, 6, -destinationCoordinatesNormalized[i].y() * sourceCoordinatesNormalized[i].x() );
300 gsl_matrix_set( S, i * 2 + 1, 7, -destinationCoordinatesNormalized[i].y() * sourceCoordinatesNormalized[i].y() );
301 gsl_matrix_set( S, i * 2 + 1, 8, -destinationCoordinatesNormalized[i].y() * 1.0 );
302 }
303
304 if ( destinationCoordinatesNormalized.size() == 4 )
305 {
306 // The GSL SVD routine only supports matrices with rows >= columns (m >= n)
307 // Unfortunately, we can't use the SVD of the transpose (i.e. S^T = (U D V^T)^T = V D U^T)
308 // to work around this, because the solution lies in the right nullspace of S, and
309 // gsl only supports a thin SVD of S^T, which does not return these vectors.
310
311 // HACK: duplicate last row to get a 9x9 equation system
312 for ( int j = 0; j < 9; j++ )
313 {
314 gsl_matrix_set( S, 8, j, gsl_matrix_get( S, 7, j ) );
315 }
316 }
317
318 // Solve Sh = 0 in the total least squares sense, i.e.
319 // with Sh = min and |h|=1. The solution "h" is given by the
320 // right singular eigenvector of S corresponding, to the smallest
321 // singular value (via SVD).
322 gsl_matrix *V = gsl_matrix_alloc( n, n );
323 gsl_vector *singular_values = gsl_vector_alloc( n );
324 gsl_vector *work = gsl_vector_alloc( n );
325
326 // V = n x n
327 // U = m x n (thin SVD) U D V^T
328 gsl_linalg_SV_decomp( S, V, singular_values, work );
329
330 // Columns of V store the right singular vectors of S
331 for ( unsigned int i = 0; i < n; i++ )
332 {
333 H[i] = gsl_matrix_get( V, i, n - 1 );
334 }
335
336 gsl_matrix *prodMatrix = gsl_matrix_alloc( 3, 3 );
337
338 gsl_matrix_view Hmatrix = gsl_matrix_view_array( H, 3, 3 );
339 const gsl_matrix_view normSourceMatrix = gsl_matrix_view_array( normSource, 3, 3 );
340 const gsl_matrix_view denormDestMatrix = gsl_matrix_view_array( denormDest, 3, 3 );
341
342 // Change coordinate frame of image and pre-image from normalized to destination and source coordinates.
343 // H' = denormalizeMapCoords*H*normalizePixelCoords
344 gsl_blas_dgemm( CblasNoTrans, CblasNoTrans, 1.0, &Hmatrix.matrix, &normSourceMatrix.matrix, 0.0, prodMatrix );
345 gsl_blas_dgemm( CblasNoTrans, CblasNoTrans, 1.0, &denormDestMatrix.matrix, prodMatrix, 0.0, &Hmatrix.matrix );
346
347 gsl_matrix_free( prodMatrix );
348 gsl_matrix_free( S );
349 gsl_matrix_free( V );
350 gsl_vector_free( singular_values );
351 gsl_vector_free( work );
352#endif
353}
static void helmert(const QVector< QgsPointXY > &sourceCoordinates, const QVector< QgsPointXY > &destinationCoordinates, QgsPointXY &origin, double &pixelSize, double &rotation)
Transforms the point at origin in-place, using a helmert transformation calculated from the list of s...
static void projective(const QVector< QgsPointXY > &sourceCoordinates, const QVector< QgsPointXY > &destinationCoordinates, double H[9])
Calculates projective parameters from the list of source and destination Ground Control Points (GCPs)...
static void linear(const QVector< QgsPointXY > &sourceCoordinates, const QVector< QgsPointXY > &destinationCoordinates, QgsPointXY &origin, double &pixelXSize, double &pixelYSize)
Transforms the point at origin in-place, using a linear transformation calculated from the list of so...
Custom exception class which is raised when an operation is not supported.
A class to represent a 2D point.
Definition qgspointxy.h:60
void setY(double y)
Sets the y value of the point.
Definition qgspointxy.h:129
void setX(double x)
Sets the x value of the point.
Definition qgspointxy.h:119
void normalizeCoordinates(const QVector< QgsPointXY > &coords, QVector< QgsPointXY > &normalizedCoords, double normalizeMatrix[9], double denormalizeMatrix[9])
Scales the given coordinates so that the center of gravity is at the origin and the mean distance to ...