postprocessing

displacement filters

functions of use for the extraction of information of displacement fields, or filtering thereof.

dhdt.postprocessing.displacement_filters.local_infilling_filter(Z, tsize=5)

apply local inverse distance weighting, on no-data values in an array

Parameters
  • Z (numpy.array, size=(m,n), dtype=float) – intensity array with no-data values

  • tsize (integer, {x ∈ ℕ | x ≥ 1}, unit=pixel) – size of the neighborhood

Returns

Z_new – intensity array with interpolated values at the data gap locations

Return type

numpy.array, size=(m,n), dtype=float

References

Jo02

Joughin “Ice-sheet velocity mapping: a combined interferometric and speckle-tracking approach”, Annuals of glaciology vol.34 pp.195-201, 2002.

Jo17

Joughin et al. “Greenland ice mapping project 2 (GIMP-2) algorithm theoretical basis document”, Making earth system data records for use in research environment (MEaSUREs) documentation, 2017.

dhdt.postprocessing.displacement_filters.local_variance(V, tsize=5)

local non-linear variance calculation

Parameters
  • V (numpy.array, size=(m,n), dtype=float) – array with one velocity component, all algorithms are indepent of their axis.

  • sig_V (numpy.array, size=(m,n), dtype=float) – statistical local variance, based on the procedure described in [JoXX], [JoYY]

References

Jo02

Joughin “Ice-sheet velocity mapping: a combined interferometric and speckle-tracking approach”, Annuals of glaciology vol.34 pp.195-201, 2002.

Jo17

Joughin et al. “Greenland ice mapping project 2 (GIMP-2) algorithm theoretical basis document”, Making earth system data records for use in research environment (MEaSUREs) documentation, 2017.

dhdt.postprocessing.displacement_filters.nan_resistant_filter(Z, kernel)

operate zonal filter on image with missing values. Energy within the kernel that fals on missing instances is distributed over other elements. Hence, regions with no data are shrinking, instead of growing. For a more detailed describtion see [Al22]

Parameters
  • Z ({numpy.ndarray, numpy.masked.array}, size=(m,n)) – grid with values

  • kernel (numpy.ndarray, size=(k,l)) – convolution kernel

  • Returns

  • Z_new (numpy.ndarray, size=(m,n)) – convolved estimate

References

Al22

Altena et al. “Correlation dispersion as a measure to better estimate uncertainty of remotely sensed glacier displacements”, The Crysophere, vol.16(6) pp.2285–2300, 2022.

solar tools

functions of use for illumination calculations

dhdt.postprocessing.solar_tools.az_to_sun_vector(az, indexing='ij')

transform azimuth angle to 2D-unit vector

Parameters
  • az (float, unit=degrees) – azimuth of sun.

  • indexing ({‘xy’, ‘ij’}) –

    • “xy” : using map coordinates

    • ”ij” : using local image coordinates

Returns

sun – unit vector in the direction of the sun.

Return type

numpy.array, size=(2,1), range=0…1

Notes

The azimuth angle declared in the following coordinate frame:

     ^ North & y
     |
- <--┼--> +
     |
     ┼----> East & x

The angles related to the suns’ heading are as follows:

surface normal              * sun
^                     ^    /
|                     |   /
├-- zenith angle      |  /
| /                   | /|
|/                    |/ | elevation angle
└----                 └--┴---

Two different coordinate system are used here:

indexing   |           indexing    ^ y
system 'ij'|           system 'xy' |
           |                       |
           |       i               |       x
   --------┼-------->      --------┼-------->
           |                       |
           |                       |
image      | j         map         |
based      v           based       |
dhdt.postprocessing.solar_tools.make_doppler_range(Z, az, zn, Lambertian=True, spac=10)
Parameters
  • Z (numpy.array, unit=meters) – array with elevation values

  • az (float, unit=degrees, range=-180...+180) – flight orientation of the satellite

  • zn ({float,array}, unit=degrees, range=0...+90) – illumination angle from the satellite

Notes

dhdt.postprocessing.solar_tools.make_shading(Z, az, zn, spac=10)

create synthetic shading image from given sun angles

A simple Lambertian reflection model is used here.

Parameters
  • Z (numpy.array, size=(m,n), dtype={integer,float}, unit=meter) – grid with elevation data

  • az (float, unit=degrees) – azimuth angle

  • zn (float, unit=degrees) – zenith angle

  • spac ({float, tuple}, default=10, unit=meter) – resolution of the square grid. Optionally the geoTransform is given.

Returns

Sh – estimated shading grid

Return type

numpy.array, size=(m,n), dtype=float, range=0…1

Notes

The azimuth angle declared in the following coordinate frame:

     ^ North & y
     |
- <--┼--> +
     |
     ┼----> East & x

The angles related to the sun are as follows:

surface normal              * sun
^                     ^    /
|                     |   /
├-- zenith angle      |  /
| /                   | /|
|/                    |/ | elevation angle
└----                 └--┴---
dhdt.postprocessing.solar_tools.make_shading_minnaert(Z, az, zn, k=1, spac=10)

create synthetic shading image from given sun angles

A simple Minnaert reflection model is used here.

Parameters
  • Z (numpy.array, size=(m,n), dtype={integer,float}, unit=meter) – grid with elevation data

  • az (float, unit=degrees) – azimuth angle

  • zn (float, unit=degrees) – zenith angle

  • spac (float, default=10, unit=meter) – resolution of the square grid.

Returns

Sh – estimated shading grid

Return type

numpy.array, size=(m,n), dtype=float, range=0…1

Notes

The azimuth angle declared in the following coordinate frame:

     ^ North & y
     |
- <--┼--> +
     |
     ┼----> East & x

The angles related to the sun are as follows:

surface normal              * sun
^                     ^    /
|                     |   /
├-- zenith angle      |  /
| /                   | /|
|/                    |/ | elevation angle
└----                 └--┴---
dhdt.postprocessing.solar_tools.make_shadowing(Z, az, zn, spac=10, weights=None, radiation=False, border='copy')

create synthetic shadow image from given sun angles

Parameters
  • Z (numpy.array, size=(m,n), dtype={integer,float}) – grid with elevation data

  • az (float, unit=degrees) – azimuth angle

  • zn ({float, numpy.array}, unit=degrees) – zenith angle

  • spac (float, optional) – resolution of the square grid. The default is 10.

  • weights (numpy.array) – if the zenith angle is a vector, it can be used to get a weighted sum

  • border (string, {'copy','zeros'}) – specify what to do at the border, since it is a zonal operation

Returns

Sw – estimated shadow grid

Return type

numpy.array, size=(m,n), dtype=bool

Notes

The azimuth angle declared in the following coordinate frame:

     ^ North & y
     |
- <--┼--> +
     |
     ┼----> East & x

The angles related to the sun are as follows:

surface normal              * sun
^                     ^    /
|                     |   /
├-- zenith angle      |  /
| /                   | /|
|/                    |/ | elevation angle
└----                 └--┴---
dhdt.postprocessing.solar_tools.sun_angles_to_vector(az, zn, indexing='ij')

transform azimuth and zenith angle to 3D-unit vector

Parameters
  • az ({float, numpy.ndarray}, unit=degrees) – azimuth angle of sun.

  • zn ({float, numpy.ndarray}, unit=degrees) – zenith angle of sun.

  • indexing ({‘xy’, ‘ij’}) –

    • “xy” : using map coordinates

    • ”ij” : using local image coordinates

Returns

sun – unit vector in the direction of the sun.

Return type

numpy.array, size=(3,1), dtype=float, range=0…1

Notes

The azimuth angle declared in the following coordinate frame:

     ^ North & y
     |
- <--┼--> +
     |
     ┼----> East & x

The angles related to the sun are as follows:

surface normal              * sun
^                     ^    /
|                     |   /
├-- zenith angle      |  /
| /                   | /|
|/                    |/ | elevation angle
└---- surface -----   └--┴---

Two different coordinate system are used here:

indexing   |           indexing    ^ y
system 'ij'|           system 'xy' |
           |                       |
           |       i               |       x
   --------┼-------->      --------┼-------->
           |                       |
           |                       |
image      | j         map         |
based      v           based       |
dhdt.postprocessing.solar_tools.vector_to_sun_angles(sun, indexing='ij')

transform 3D-unit vector to azimuth and zenith angle

Parameters

sun ({float,numpy.array}, size=(m,), dtype=float, range=0...1) – the three different unit vectors in the direction of the sun.

Returns

  • az ({float, numpy.ndarray}, unit=degrees) – azimuth angle of sun.

  • zn ({float, numpy.ndarray}, unit=degrees) – zenith angle of sun.

  • The angles related to the suns’ heading are as follows

    surface normal              * sun
    ^                     ^    /
    |                     |   /
    ├-- zenith angle      |  /
    | /                   | /|
    |/                    |/ | elevation angle
    └----                 └--┴---
    

terrain tools

dhdt.postprocessing.terrain_tools.d8_catchment(V_x, V_y, geoTransform, x, y, disperse=True, flat=True)

estimate the catchment behind a seed point

Parameters
  • V_x (numpy.ndarray, size=(m,n)) – grids with horizontal and vertical directions

  • V_y (numpy.ndarray, size=(m,n)) – grids with horizontal and vertical directions

  • geoTransform

  • x ({float, numpy.ndarray}) – map locations of the seeds

  • y ({float, numpy.ndarray}) – map locations of the seeds

Returns

C – grid with the specific catchments

Return type

numpy.ndarray, size=(m,n), dtype=boolean

See also

get_d8_dir, d8_flow

dhdt.postprocessing.terrain_tools.d8_flow(Class, V_x, V_y, iter=20, analysis=True)

re-organize labelled dataset by steering towards convergent zones

Parameters
  • Class (numpy.ndarray, size=(m,n)) – labelled array

  • V_x (numpy.ndarray, size=(m,n)) – grids with horizontal and vertical directions

  • V_y (numpy.ndarray, size=(m,n)) – grids with horizontal and vertical directions

  • iter (integer, {x ∈ ℕ | x ≥ 0}, default=20) – number of iterations to be executed

  • analysis (boolean) – keeps intermediate results, to see the evolution of the algorithm

Returns

Class_new – newly evolved labelled array

Return type

numpy.ndarray, size=(m,n)

See also

get_d8_dir, d8_flow

dhdt.postprocessing.terrain_tools.get_d8_dir(V_x, V_y)

get directional data of direct neighbors, following [GM97].

Parameters
  • V_x (numpy.ndarray, size=(m,n), float) – array of displacements in horizontal direction

  • V_y (numpy.ndarray, size=(m,n), float) – array of displacements in vertical direction

Returns

d8_dir – compass directions starting from North going in clockwise direction

Return type

numpy.ndarray, size=(m,n), range=0…7

Notes

It is important to know what type of coordinate systems is used, thus:

     N             0
  NW | NE       7  |  1
W----+----E   6----8----2
  SW | SE       5  |  3
     S             4

References

GM97

Garbrecht & Martz, “The assignment of drainage direction over flat surfaces in raster digital elevation models”, Journal of hydrology, vol.193 pp.204-213, 1997.

dhdt.postprocessing.terrain_tools.get_d8_offsets(shape, order='C')

get the offsets of the direct neighbors for linear indexing

Parameters
  • shape (list) – dimensions of the array of interest

  • order ({'C','F'}) – index order in either C-like or Fortran like order

Returns

d8_offsets – relative offsets of neighbors

Return type

numpy.ndarray, size=(8)

Notes

Two different orders are used:

0 1 2 3        0 4 8
4 5 6 7        1 5 9
8 9 ...        2 6 ...

order 'C'      order 'F'