31void QgsLeastSquares::linear(
const QVector<QgsPointXY> &sourceCoordinates,
const QVector<QgsPointXY> &destinationCoordinates,
QgsPointXY &origin,
double &pixelXSize,
double &pixelYSize )
33 const int n = destinationCoordinates.size();
36 throw std::domain_error( QObject::tr(
"Fit to a linear transform requires at least 2 points." ).toLocal8Bit().constData() );
39 double sumPx( 0 ), sumPy( 0 ), sumPx2( 0 ), sumPy2( 0 ), sumPxMx( 0 ), sumPyMy( 0 ), sumMx( 0 ), sumMy( 0 );
40 for (
int i = 0; i < n; ++i )
42 sumPx += sourceCoordinates.at( i ).x();
43 sumPy += sourceCoordinates.at( i ).y();
44 sumPx2 += std::pow( sourceCoordinates.at( i ).x(), 2 );
45 sumPy2 += std::pow( sourceCoordinates.at( i ).y(), 2 );
46 sumPxMx += sourceCoordinates.at( i ).x() * destinationCoordinates.at( i ).x();
47 sumPyMy += sourceCoordinates.at( i ).y() * destinationCoordinates.at( i ).y();
48 sumMx += destinationCoordinates.at( i ).x();
49 sumMy += destinationCoordinates.at( i ).y();
52 const double deltaX = n * sumPx2 - std::pow( sumPx, 2 );
53 const double deltaY = n * sumPy2 - std::pow( sumPy, 2 );
55 const double aX = ( sumPx2 * sumMx - sumPx * sumPxMx ) / deltaX;
56 const double aY = ( sumPy2 * sumMy - sumPy * sumPyMy ) / deltaY;
57 const double bX = ( n * sumPxMx - sumPx * sumMx ) / deltaX;
58 const double bY = ( n * sumPyMy - sumPy * sumMy ) / deltaY;
63 pixelXSize = std::fabs( bX );
64 pixelYSize = std::fabs( bY );
68void QgsLeastSquares::helmert(
const QVector<QgsPointXY> &sourceCoordinates,
const QVector<QgsPointXY> &destinationCoordinates,
QgsPointXY &origin,
double &pixelSize,
double &rotation )
71 ( void ) sourceCoordinates;
72 ( void ) destinationCoordinates;
76 throw QgsNotSupportedException( QObject::tr(
"Calculating a helmert transformation requires a QGIS build based GSL" ) );
78 const int n = destinationCoordinates.size();
81 throw std::domain_error( QObject::tr(
"Fit to a Helmert transform requires at least 2 points." ).toLocal8Bit().constData() );
94 for (
int i = 0; i < n; ++i )
96 A += sourceCoordinates.at( i ).x();
97 B += sourceCoordinates.at( i ).y();
98 C += destinationCoordinates.at( i ).x();
99 D += destinationCoordinates.at( i ).y();
100 E += destinationCoordinates.at( i ).x() * sourceCoordinates.at( i ).x();
101 F += destinationCoordinates.at( i ).y() * sourceCoordinates.at( i ).y();
102 G += std::pow( sourceCoordinates.at( i ).x(), 2 );
103 H += std::pow( sourceCoordinates.at( i ).y(), 2 );
104 I += destinationCoordinates.at( i ).x() * sourceCoordinates.at( i ).y();
105 J += sourceCoordinates.at( i ).x() * destinationCoordinates.at( i ).y();
113 double MData[] = { A, -B, ( double ) n, 0., B, A, 0., (
double ) n, G + H, 0., A, B, 0., G + H, -B, A };
115 double bData[] = { C, D, E + F, J - I };
118 gsl_matrix_view M = gsl_matrix_view_array( MData, 4, 4 );
119 const gsl_vector_view b = gsl_vector_view_array( bData, 4 );
120 gsl_vector *x = gsl_vector_alloc( 4 );
121 gsl_permutation *p = gsl_permutation_alloc( 4 );
123 gsl_linalg_LU_decomp( &M.matrix, p, &s );
124 gsl_linalg_LU_solve( &M.matrix, p, &b.vector, x );
125 gsl_permutation_free( p );
127 origin.
setX( gsl_vector_get( x, 2 ) );
128 origin.
setY( gsl_vector_get( x, 3 ) );
129 pixelSize = std::sqrt( std::pow( gsl_vector_get( x, 0 ), 2 ) + std::pow( gsl_vector_get( x, 1 ), 2 ) );
130 rotation = std::atan2( gsl_vector_get( x, 1 ), gsl_vector_get( x, 0 ) );
132 gsl_vector_free( x );
196void normalizeCoordinates(
const QVector<QgsPointXY> &coords, QVector<QgsPointXY> &normalizedCoords,
double normalizeMatrix[9],
double denormalizeMatrix[9] )
199 double cogX = 0.0, cogY = 0.0;
200 for (
int i = 0; i < coords.size(); i++ )
202 cogX += coords[i].x();
203 cogY += coords[i].y();
205 cogX *= 1.0 / coords.size();
206 cogY *= 1.0 / coords.size();
209 double meanDist = 0.0;
210 for (
int i = 0; i < coords.size(); i++ )
212 const double X = ( coords[i].x() - cogX );
213 const double Y = ( coords[i].y() - cogY );
214 meanDist += std::sqrt( X * X + Y * Y );
216 meanDist *= 1.0 / coords.size();
218 const double OOD = meanDist * M_SQRT1_2;
219 const double D = 1.0 / OOD;
220 normalizedCoords.resize( coords.size() );
221 for (
int i = 0; i < coords.size(); i++ )
223 normalizedCoords[i] =
QgsPointXY( ( coords[i].x() - cogX ) * D, ( coords[i].y() - cogY ) * D );
226 normalizeMatrix[0] = D;
227 normalizeMatrix[1] = 0.0;
228 normalizeMatrix[2] = -cogX * D;
229 normalizeMatrix[3] = 0.0;
230 normalizeMatrix[4] = D;
231 normalizeMatrix[5] = -cogY * D;
232 normalizeMatrix[6] = 0.0;
233 normalizeMatrix[7] = 0.0;
234 normalizeMatrix[8] = 1.0;
236 denormalizeMatrix[0] = OOD;
237 denormalizeMatrix[1] = 0.0;
238 denormalizeMatrix[2] = cogX;
239 denormalizeMatrix[3] = 0.0;
240 denormalizeMatrix[4] = OOD;
241 denormalizeMatrix[5] = cogY;
242 denormalizeMatrix[6] = 0.0;
243 denormalizeMatrix[7] = 0.0;
244 denormalizeMatrix[8] = 1.0;
252 ( void ) sourceCoordinates;
253 ( void ) destinationCoordinates;
255 throw QgsNotSupportedException( QObject::tr(
"Calculating a projective transformation requires a QGIS build based GSL" ) );
257 Q_ASSERT( sourceCoordinates.size() == destinationCoordinates.size() );
259 if ( destinationCoordinates.size() < 4 )
261 throw std::domain_error( QObject::tr(
"Fitting a projective transform requires at least 4 corresponding points." ).toLocal8Bit().constData() );
264 QVector<QgsPointXY> sourceCoordinatesNormalized;
265 QVector<QgsPointXY> destinationCoordinatesNormalized;
267 double normSource[9], denormSource[9];
268 double normDest[9], denormDest[9];
269 normalizeCoordinates( sourceCoordinates, sourceCoordinatesNormalized, normSource, denormSource );
270 normalizeCoordinates( destinationCoordinates, destinationCoordinatesNormalized, normDest, denormDest );
274 const uint m = std::max( 9u, ( uint ) destinationCoordinatesNormalized.size() * 2u );
276 gsl_matrix *S = gsl_matrix_alloc( m, n );
278 for (
int i = 0; i < destinationCoordinatesNormalized.size(); i++ )
280 gsl_matrix_set( S, i * 2, 0, sourceCoordinatesNormalized[i].x() );
281 gsl_matrix_set( S, i * 2, 1, sourceCoordinatesNormalized[i].y() );
282 gsl_matrix_set( S, i * 2, 2, 1.0 );
284 gsl_matrix_set( S, i * 2, 3, 0.0 );
285 gsl_matrix_set( S, i * 2, 4, 0.0 );
286 gsl_matrix_set( S, i * 2, 5, 0.0 );
288 gsl_matrix_set( S, i * 2, 6, -destinationCoordinatesNormalized[i].x() * sourceCoordinatesNormalized[i].x() );
289 gsl_matrix_set( S, i * 2, 7, -destinationCoordinatesNormalized[i].x() * sourceCoordinatesNormalized[i].y() );
290 gsl_matrix_set( S, i * 2, 8, -destinationCoordinatesNormalized[i].x() * 1.0 );
292 gsl_matrix_set( S, i * 2 + 1, 0, 0.0 );
293 gsl_matrix_set( S, i * 2 + 1, 1, 0.0 );
294 gsl_matrix_set( S, i * 2 + 1, 2, 0.0 );
296 gsl_matrix_set( S, i * 2 + 1, 3, sourceCoordinatesNormalized[i].x() );
297 gsl_matrix_set( S, i * 2 + 1, 4, sourceCoordinatesNormalized[i].y() );
298 gsl_matrix_set( S, i * 2 + 1, 5, 1.0 );
300 gsl_matrix_set( S, i * 2 + 1, 6, -destinationCoordinatesNormalized[i].y() * sourceCoordinatesNormalized[i].x() );
301 gsl_matrix_set( S, i * 2 + 1, 7, -destinationCoordinatesNormalized[i].y() * sourceCoordinatesNormalized[i].y() );
302 gsl_matrix_set( S, i * 2 + 1, 8, -destinationCoordinatesNormalized[i].y() * 1.0 );
305 if ( destinationCoordinatesNormalized.size() == 4 )
313 for (
int j = 0; j < 9; j++ )
315 gsl_matrix_set( S, 8, j, gsl_matrix_get( S, 7, j ) );
323 gsl_matrix *V = gsl_matrix_alloc( n, n );
324 gsl_vector *singular_values = gsl_vector_alloc( n );
325 gsl_vector *work = gsl_vector_alloc( n );
329 gsl_linalg_SV_decomp( S, V, singular_values, work );
332 for (
unsigned int i = 0; i < n; i++ )
334 H[i] = gsl_matrix_get( V, i, n - 1 );
337 gsl_matrix *prodMatrix = gsl_matrix_alloc( 3, 3 );
339 gsl_matrix_view Hmatrix = gsl_matrix_view_array( H, 3, 3 );
340 const gsl_matrix_view normSourceMatrix = gsl_matrix_view_array( normSource, 3, 3 );
341 const gsl_matrix_view denormDestMatrix = gsl_matrix_view_array( denormDest, 3, 3 );
345 gsl_blas_dgemm( CblasNoTrans, CblasNoTrans, 1.0, &Hmatrix.matrix, &normSourceMatrix.matrix, 0.0, prodMatrix );
346 gsl_blas_dgemm( CblasNoTrans, CblasNoTrans, 1.0, &denormDestMatrix.matrix, prodMatrix, 0.0, &Hmatrix.matrix );
348 gsl_matrix_free( prodMatrix );
349 gsl_matrix_free( S );
350 gsl_matrix_free( V );
351 gsl_vector_free( singular_values );
352 gsl_vector_free( work );