Skip to content

Commit

Permalink
Merge pull request OSGeo#8573 from rouault/fix_8572
Browse files Browse the repository at this point in the history
TPS transformer: use an iterative method to refine the inverse transformation (fixes OSGeo#8572)
  • Loading branch information
rouault authored Oct 27, 2023
2 parents b5a6279 + 01eefef commit 35ce51c
Show file tree
Hide file tree
Showing 6 changed files with 232 additions and 8 deletions.
4 changes: 3 additions & 1 deletion alg/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -35,7 +35,9 @@ add_library(
rasterfill.cpp
thinplatespline.cpp
gdal_simplesurf.cpp
viewshed.cpp)
viewshed.cpp
gdalgenericinverse.cpp
)

add_dependencies(alg generate_gdal_version_h)

Expand Down
20 changes: 20 additions & 0 deletions alg/gdal_tps.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -45,6 +45,7 @@
#include "gdal_alg.h"
#include "gdal_alg_priv.h"
#include "gdal_priv.h"
#include "gdalgenericinverse.h"

CPL_CVSID("$Id$")

Expand Down Expand Up @@ -363,7 +364,26 @@ int GDALTPSTransform(void *pTransformArg, int bDstToSrc, int nPointCount,

if (bDstToSrc)
{
// Compute initial guess
psInfo->poReverse->get_point(x[i], y[i], xy_out);

const auto ForwardTransformer = [](double xIn, double yIn,
double &xOut, double &yOut,
void *pUserData)
{
double xyOut[2] = {0.0, 0.0};
TPSTransformInfo *l_psInfo =
static_cast<TPSTransformInfo *>(pUserData);
l_psInfo->poForward->get_point(xIn, yIn, xyOut);
xOut = xyOut[0];
yOut = xyOut[1];
return true;
};

// Refine the initial guess
GDALGenericInverse2D(x[i], y[i], xy_out[0], xy_out[1],
ForwardTransformer, psInfo, xy_out[0],
xy_out[1]);
x[i] = xy_out[0];
y[i] = xy_out[1];
}
Expand Down
117 changes: 117 additions & 0 deletions alg/gdalgenericinverse.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,117 @@
/******************************************************************************
*
* Project: GDAL
* Purpose: Generic method to compute inverse coordinate transformation from
* forward method
* Author: Even Rouault <even dot rouault at spatialys dot com>
*
******************************************************************************
* Copyright (c) 2023, Even Rouault <even dot rouault at spatialys dot com>
*
* Permission is hereby granted, free of charge, to any person obtaining a
* copy of this software and associated documentation files (the "Software"),
* to deal in the Software without restriction, including without limitation
* the rights to use, copy, modify, merge, publish, distribute, sublicense,
* and/or sell copies of the Software, and to permit persons to whom the
* Software is furnished to do so, subject to the following conditions:
*
* The above copyright notice and this permission notice shall be included
* in all copies or substantial portions of the Software.
*
* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
* OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
* FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
* THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
* FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
* DEALINGS IN THE SOFTWARE.
****************************************************************************/

#include <algorithm>
#include <cmath>

#include "gdalgenericinverse.h"
#include <cstdio>

/** Compute (xOut, yOut) corresponding to input (xIn, yIn) using
* the provided forward transformation to emulate the reverse direction.
*
* Uses Newton-Raphson method, extended to 2D variables, that is using
* inversion of the Jacobian 2D matrix of partial derivatives. The derivatives
* are estimated numerically from the forward method evaluated at close points.
*
* Starts with initial guess provided by user in (guessedXOut, guessedYOut)
*
* It iterates at most for 15 iterations or as soon as we get below the specified
* tolerance (on input coordinates)
*/
bool GDALGenericInverse2D(double xIn, double yIn, double guessedXOut,
double guessedYOut,
GDALForwardCoordTransformer pfnForwardTranformer,
void *pfnForwardTranformerUserData, double &xOut,
double &yOut, double toleranceOnInputCoordinates)
{
const double dfAbsValOut = std::max(fabs(guessedXOut), fabs(guessedYOut));
const double dfEps = dfAbsValOut > 0 ? dfAbsValOut * 1e-6 : 1e-6;
if (toleranceOnInputCoordinates == 0)
{
const double dfAbsValIn = std::max(fabs(xIn), fabs(yIn));
toleranceOnInputCoordinates =
dfAbsValIn > 0 ? dfAbsValIn * 1e-12 : 1e-12;
}
xOut = guessedXOut;
yOut = guessedYOut;
double deriv_lam_X = 0;
double deriv_lam_Y = 0;
double deriv_phi_X = 0;
double deriv_phi_Y = 0;
for (int i = 0; i < 15; i++)
{
double xApprox;
double yApprox;
if (!pfnForwardTranformer(xOut, yOut, xApprox, yApprox,
pfnForwardTranformerUserData))
return false;
const double deltaX = xApprox - xIn;
const double deltaY = yApprox - yIn;
if (fabs(deltaX) < toleranceOnInputCoordinates &&
fabs(deltaY) < toleranceOnInputCoordinates)
{
return true;
}

// Compute Jacobian matrix
double xTmp = xOut + dfEps;
double yTmp = yOut;
double xTmpOut;
double yTmpOut;
if (!pfnForwardTranformer(xTmp, yTmp, xTmpOut, yTmpOut,
pfnForwardTranformerUserData))
return false;
const double deriv_X_lam = (xTmpOut - xApprox) / dfEps;
const double deriv_Y_lam = (yTmpOut - yApprox) / dfEps;

xTmp = xOut;
yTmp = yOut + dfEps;
if (!pfnForwardTranformer(xTmp, yTmp, xTmpOut, yTmpOut,
pfnForwardTranformerUserData))
return false;
const double deriv_X_phi = (xTmpOut - xApprox) / dfEps;
const double deriv_Y_phi = (yTmpOut - yApprox) / dfEps;

// Inverse of Jacobian matrix
const double det =
deriv_X_lam * deriv_Y_phi - deriv_X_phi * deriv_Y_lam;
if (det != 0)
{
deriv_lam_X = deriv_Y_phi / det;
deriv_lam_Y = -deriv_X_phi / det;
deriv_phi_X = -deriv_Y_lam / det;
deriv_phi_Y = deriv_X_lam / det;
}

xOut -= deltaX * deriv_lam_X + deltaY * deriv_lam_Y;
yOut -= deltaX * deriv_phi_X + deltaY * deriv_phi_Y;
}
return false;
}
44 changes: 44 additions & 0 deletions alg/gdalgenericinverse.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,44 @@
/******************************************************************************
*
* Project: GDAL
* Purpose: Generic method to compute inverse coordinate transformation from
* forward method
* Author: Even Rouault <even dot rouault at spatialys dot com>
*
******************************************************************************
* Copyright (c) 2023, Even Rouault <even dot rouault at spatialys dot com>
*
* Permission is hereby granted, free of charge, to any person obtaining a
* copy of this software and associated documentation files (the "Software"),
* to deal in the Software without restriction, including without limitation
* the rights to use, copy, modify, merge, publish, distribute, sublicense,
* and/or sell copies of the Software, and to permit persons to whom the
* Software is furnished to do so, subject to the following conditions:
*
* The above copyright notice and this permission notice shall be included
* in all copies or substantial portions of the Software.
*
* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
* OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
* FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
* THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
* FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
* DEALINGS IN THE SOFTWARE.
****************************************************************************/

#ifndef GDALGENERICINVERSE_H
#define GDALGENERICINVERSE_H
#include <stdbool.h>

typedef bool (*GDALForwardCoordTransformer)(double xIn, double yIn,
double &xOut, double &yOt,
void *pUserData);

bool GDALGenericInverse2D(double xIn, double yIn, double guessedXOut,
double guessedYOut,
GDALForwardCoordTransformer pfnForwardTranformer,
void *pfnForwardTranformerUserData, double &xOut,
double &yOut, double toleranceOnInputCoordinates = 0);

#endif // GDALGENERICINVERSE_H
2 changes: 1 addition & 1 deletion autotest/alg/warp.py
Original file line number Diff line number Diff line change
Expand Up @@ -1515,7 +1515,7 @@ def test_warp_55():

ds = gdal.Open("data/warpedvrt_with_ovr.vrt")
cs = ds.GetRasterBand(1).GetOverview(0).Checksum()
assert cs == 25478
assert cs == 51966
ds = None


Expand Down
53 changes: 47 additions & 6 deletions autotest/gcore/transformer.py
Original file line number Diff line number Diff line change
Expand Up @@ -643,11 +643,11 @@ def test_transformer_12():
"""
<VRTDataset rasterXSize="20" rasterYSize="20">
<GCPList Projection="PROJCS[&quot;NAD27 / UTM zone 11N&quot;,GEOGCS[&quot;NAD27&quot;,DATUM[&quot;North_American_Datum_1927&quot;,SPHEROID[&quot;Clarke 1866&quot;,6378206.4,294.9786982139006,AUTHORITY[&quot;EPSG&quot;,&quot;7008&quot;]],AUTHORITY[&quot;EPSG&quot;,&quot;6267&quot;]],PRIMEM[&quot;Greenwich&quot;,0],UNIT[&quot;degree&quot;,0.0174532925199433],AUTHORITY[&quot;EPSG&quot;,&quot;4267&quot;]],PROJECTION[&quot;Transverse_Mercator&quot;],PARAMETER[&quot;latitude_of_origin&quot;,0],PARAMETER[&quot;central_meridian&quot;,-117],PARAMETER[&quot;scale_factor&quot;,0.9996],PARAMETER[&quot;false_easting&quot;,500000],PARAMETER[&quot;false_northing&quot;,0],UNIT[&quot;metre&quot;,1,AUTHORITY[&quot;EPSG&quot;,&quot;9001&quot;]],AUTHORITY[&quot;EPSG&quot;,&quot;26711&quot;]]">
<GCP Id="" Pixel="0" Line="0" X="0" Y="0"/>
<GCP Id="" Pixel="20" Line="0" X="20" Y="0"/>
<GCP Id="" Pixel="0" Line="20" X="0" Y="20"/>
<GCP Id="" Pixel="20" Line="20" X="20" Y="20"/>
<GCP Id="" Pixel="0" Line="0" X="0" Y="0"/> <!-- duplicate entry -->
<GCP Id="" Pixel="0" Line="0" X="-1" Y="-1"/>
<GCP Id="" Pixel="20" Line="0" X="1" Y="-1"/>
<GCP Id="" Pixel="0" Line="20" X="-1" Y="1"/>
<GCP Id="" Pixel="20" Line="20" X="1" Y="1"/>
<GCP Id="" Pixel="0" Line="0" X="-1" Y="-1"/> <!-- duplicate entry -->
</GCPList>
<VRTRasterBand dataType="Byte" band="1">
<ColorInterp>Gray</ColorInterp>
Expand All @@ -661,12 +661,26 @@ def test_transformer_12():
tr = gdal.Transformer(ds, None, ["METHOD=GCP_TPS"])
assert tr is not None
(success, pnt) = tr.TransformPoint(0, 0, 0)
assert (
success
and pnt[0] == pytest.approx(-1, abs=1e-7)
and pnt[1] == pytest.approx(-1, abs=1e-7)
)

(success, pnt) = tr.TransformPoint(1, -1, -1)
assert (
success
and pnt[0] == pytest.approx(0, abs=1e-7)
and pnt[1] == pytest.approx(0, abs=1e-7)
)

(success, pnt) = tr.TransformPoint(1, 0.2, 0.8)
assert (
success
and pnt[0] == pytest.approx(12, abs=1e-7)
and pnt[1] == pytest.approx(18, abs=1e-7)
)

ds = gdal.Open(
"""
<VRTDataset rasterXSize="20" rasterYSize="20">
Expand Down Expand Up @@ -1009,13 +1023,40 @@ def test_transformer_tps_precision():

success = True
maxDiffResult = 0.0
for gcp in ds.GetGCPs():
gcps = ds.GetGCPs()
for i, gcp in enumerate(gcps):
(s, result) = tr.TransformPoint(0, gcp.GCPPixel, gcp.GCPLine)
success &= s
diffResult = math.sqrt(
(gcp.GCPX - result[0]) ** 2 + (gcp.GCPY - result[1]) ** 2
)
maxDiffResult = max(maxDiffResult, diffResult)

(s, result) = tr.TransformPoint(1, result[0], result[1])
assert s
assert result[0] == pytest.approx(gcp.GCPPixel), (i, result)
assert result[1] == pytest.approx(gcp.GCPLine), (i, result)

if i + 1 < len(gcps):
# Try transforming "random" points
xIn = 0.5 * (gcps[i].GCPPixel + gcps[i + 1].GCPPixel)
yIn = 0.5 * (gcps[i].GCPLine + gcps[i + 1].GCPLine)
(s, result) = tr.TransformPoint(0, xIn, yIn)
assert s
(s, result) = tr.TransformPoint(1, result[0], result[1])
assert s
# A few points fail to converge with reasonable accuracy
if i not in (
775,
776,
787,
1010,
1634,
1638,
):
assert result[0] == pytest.approx(xIn, abs=1e-3), (i, xIn, yIn, result)
assert result[1] == pytest.approx(yIn, abs=1e-3), (i, xIn, yIn, result)

assert success, "at least one point could not be transformed"
assert maxDiffResult < 1e-3, "at least one transformation exceeds the error bound"

Expand Down

0 comments on commit 35ce51c

Please sign in to comment.