preprocessing¶
image transforms¶
functions adjusting the gray scale of an image
- dhdt.preprocessing.image_transforms.cum_hist(img)¶
compute the cumulative histogram of an image
- Parameters
img (numpy.array, size=(m,n), ndim=2) – grid with intensity values
- Returns
cdf – emperical cumulative distribution function, i.e. a cumulative histogram.
- Return type
numpy.array, size=(2**n,), dtype=integer
- dhdt.preprocessing.image_transforms.gamma_adjustment(Z, gamma=1.0)¶
transform intensity in non-linear way
enhances high intensities when gamma>1
- Parameters
Z (numpy.ndarray, size=(m,n), dtype=integer) – array with intensity values
gamma (float) – power law parameter
- Returns
Z_new – array with transform
- Return type
numpy.ndarray, size=(m,n)
- dhdt.preprocessing.image_transforms.general_midway_equalization(Z)¶
equalization of multiple imagery, based on [De04].
- Parameters
Z (numpy.array, dim=3, type={uint8,uint16}) – stack of imagery
- Returns
Z_new – harmonized stack of imagery
- Return type
numpy.array, dim=3, type={uint8,uint16}
References
- De04
Delon, “Midway image equalization”, Journal of mathematical imaging and Vision”, vol.21(2) pp.119-134, 2004
Notes
The following acronyms are used:
CDF - cumulative distribution function
- dhdt.preprocessing.image_transforms.get_histogram(img)¶
calculate the normalized histogram of an image
- Parameters
img (numpy.ndarray, size=(m,n)) – gray scale image array
- Returns
hist – histogram array
- Return type
numpy.ndarray, size={(2**8,2), (2**16,2)}
- dhdt.preprocessing.image_transforms.gompertz_adjustment(Z, a=None, b=None)¶
transform intensity in non-linear way, so high reflected surfaces like snow have more emphasized.
- Parameters
Z (numpy.ndarray, size=(m,n), dtype=integer) – array with intensity values
a (float) – parameters for the Gompertz function
b (float) – parameters for the Gompertz function
- Returns
Z_new – grid with transform
- Return type
numpy.ndarray, size=(m,n)
- dhdt.preprocessing.image_transforms.high_pass_im(Im, radius=10)¶
- Parameters
Im (numpy.array, ndim={2,3}, size=(m,n,_)) – data array with intensities
radius (float, unit=pixels) – hard threshold of the sphere in Fourier space
- Returns
Inew – data array with high-frequency signal
- Return type
numpy.array, ndim={2,3}, size=(m,n,_)
- dhdt.preprocessing.image_transforms.histogram_equalization(img, img_ref)¶
transform the intensity values of an array, so they have a similar histogram as the reference.
- Parameters
img (numpy.array, size=(m,n), ndim=2) – asd
img_ref (numpy.array, size=(m,n), ndim=2) –
- Returns
img – transformed image
- Return type
numpy.array, size=(m,n), ndim=2
- dhdt.preprocessing.image_transforms.hyperbolic_adjustment(Z, intercept)¶
transform intensity through an hyperbolic sine function
- Parameters
Z (numpy.ndarray, size=(m,n), dtype=float) – grid with intensity values
- Returns
Z_new – grid with transform
- Return type
numpy.ndarray, size=(m,n), dtype=float, range=-1…+1
- dhdt.preprocessing.image_transforms.inverse_tangent_transformation(x)¶
enhances high and low intensities
- Parameters
x (numpy.ndarray, size=(m,n)) – grid with intensity values
- Returns
fx – grid with transform
- Return type
numpy.ndarray, size=(_,_)
Examples
>>> import numpy as np >>> import matplotlib.pyplot as plt >>> x = np.arange(-1,1,.001) >>> y = inverse_tangent_transformation(x) >>> plt.figure(), plt.plot(x,y), plt.show() shows the mapping function
- dhdt.preprocessing.image_transforms.log_adjustment(Z)¶
transform intensity in non-linear way
enhances low intensities
- Parameters
Z (numpy.ndarray, size=(m,n)) – grid with intensity values
- Returns
Z_new – grid with transform
- Return type
numpy.ndarray, size=(m,n)
- dhdt.preprocessing.image_transforms.mat_to_gray(Z, notZ=None, vmin=None, vmax=None)¶
transform matix array to float, omitting nodata values
- Parameters
Z (numpy.ndarray, size=(m,n), dtype={integer,float}) – matrix of integers with data
notZ (numpy.ndarray, size=(m,n), dtype=bool) – matrix assigning which is a no data. The default is None.
- Returns
Znew – array with normalized intensity values
- Return type
numpy.ndarray, size=(m,n), dtype=float, range=0…1
See also
gamma_adjustment
,log_adjustment
,s_curve
,inverse_tangent_transformation
Examples
>>> import numpy as np >>> from skimage import data >>> img = data.astronaut() >>> np.max(img) 255 >>> gray = mat_to_gray(img[:,:,0], img[:,:,0]==0) >>> np.max(gray) 1.0
- dhdt.preprocessing.image_transforms.multi_spectral_grad(Ims)¶
get the dominant multi-spectral gradient of the stack, stripped down version of the fusion technique, described in [SW02].
- Parameters
Ims (numpy.array, dim=3, dtype=float) – stack of imagery
- Returns
grad_I – fused gradient of the image stack
- Return type
numpy.array, dim=2, dtype=complex
References
- SW02
Socolinsky & Wolff, “Multispectral image visualization through first-order fusion”, IEEE transactions on image fusion, vol.11(8) pp.923-931, 2002.
- dhdt.preprocessing.image_transforms.normalize_histogram(img)¶
transform image to have a uniform distribution
- Parameters
img (numpy.ndarray, size=(m,n), dtype={uint8,unit16}) – image array
- Returns
new_img – transfrormed image
- Return type
numpy.ndarray, size=(m,n), dtype={uint8,unit16}
See also
- dhdt.preprocessing.image_transforms.psuedo_inv(cumhist, val)¶
- Parameters
cumhist (numpy.array, size=(2**n,), dtype={integer, float}) – emperical cumulative distribution function, i.e. a cumulative histogram.
val ({integer,float}) –
- Returns
res – the pseudo-inverse of cumhist
- Return type
{integer,float}
Notes
The following acronyms are used:
CDF - cumulative distribution function
color transforms¶
functions transforming to and from different color spaces
- dhdt.preprocessing.color_transforms.erdas2hsi(Blue, Green, Red)¶
transform red, green, blue arrays to hue, saturation, intensity arrays
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Red (numpy.array, size=(m,n)) – red band of satellite image
- Returns
Hue (numpy.array, size=(m,n), float) – hue
Sat (numpy.array, size=(m,n), float) – saturation
Int (numpy.array, size=(m,n), float) – intensity
See also
References
- ER13
ERDAS, “User handbook”, 2013.
- dhdt.preprocessing.color_transforms.hcv2rgb(H, C, V)¶
transform a color bubble model to red green blue arrays, following [Sh95].
- Parameters
H (numpy.array, size=(m,n)) – array with amount of color
C (numpy.array, size=(m,n)) – luminance
V (numpy.array, size=(m,n)) – array with dominant frequency
- Returns
Red, Green, Blue – array with colours
- Return type
numpy.array, size=(m,n)
References
- Sh95
Shih, “The reversibility of six geometric color spaces”, Photogrammetric engineering & remote sensing, vol.61(10) pp.1223-1232, 1995.
- dhdt.preprocessing.color_transforms.lab2lch(L, a, b)¶
transform XYZ tristimulus arrays to Lab values
- Parameters
L (numpy.array, size=(m,n)) –
a (numpy.array, size=(m,n)) –
b (numpy.array, size=(m,n)) –
- Returns
C,h
- Return type
numpy.array, size=(m,n)
References
- FR98
Ford & Roberts. “Color space conversion”, pp. 1–31, 1998.
- Si18
Silva et al. “Near real-time shadow detection and removal in aerial motion imagery application” ISPRS journal of photogrammetry and remote sensing, vol.140 pp.104–121, 2018.
- dhdt.preprocessing.color_transforms.lms2lab(L, M, S)¶
transform L, M, S arrays to lab color space
- Parameters
L (numpy.array, size=(m,n)) –
M (numpy.array, size=(m,n)) –
S (numpy.array, size=(m,n)) –
- Returns
l,a,b
- Return type
numpy.array, size=(m,n)
References
- Re01
Reinhard et al. “Color transfer between images” IEEE Computer graphics and applications vol.21(5) pp.34-41, 2001.
- dhdt.preprocessing.color_transforms.rgb2hcv(Blue, Green, Red)¶
transform red green blue arrays to a color space
- Parameters
Blue (numpy.array, size=(m,n)) – Blue band of satellite image
Green (numpy.array, size=(m,n)) – Green band of satellite image
Red (numpy.array, size=(m,n)) – Red band of satellite image
- Returns
H (numpy.array, size=(m,n)) – array with amount of color
C (numpy.array, size=(m,n)) – luminance
V (numpy.array, size=(m,n)) – array with dominant frequency
References
- Sm93
Smith, “Putting colors in order”, Dr. Dobb’s Journal, pp 40, 1993.
- Ts06
Tsai, “A comparative study on shadow compensation of color aerial images in invariant color models”, IEEE transactions in geoscience and remote sensing, vol. 44(6) pp. 1661–1671, 2006.
- dhdt.preprocessing.color_transforms.rgb2hsi(Red, Green, Blue)¶
transform red, green, blue arrays to hue, saturation, intensity arrays
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Blue (numpy.array, size=(m,n)) – blue band of satellite image
- Returns
Hue (numpy.array, size=(m,n), range=0…1) – Hue
Sat (numpy.array, size=(m,n), range=0…1) – Saturation
Int (numpy.array, size=(m,n), range=0…1) – Intensity
References
- Ts06
Tsai, “A comparative study on shadow compensation of color aerial images in invariant color models”, IEEE transactions in geoscience and remote sensing, vol. 44(6) pp. 1661–1671, 2006.
- Pr91
Pratt, “Digital image processing” Wiley, 1991.
- dhdt.preprocessing.color_transforms.rgb2lms(Red, Green, Blue)¶
transform red, green, blue arrays to XYZ tristimulus values
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Blue (numpy.array, size=(m,n)) – blue band of satellite image
- Returns
L,M,S
- Return type
numpy.array, size=(m,n)
References
- Re01
Reinhard et al. “Color transfer between images” IEEE Computer graphics and applications vol.21(5) pp.34-41, 2001.
- dhdt.preprocessing.color_transforms.rgb2xyz(Red, Green, Blue, method='reinhardt')¶
transform red, green, blue arrays to XYZ tristimulus values
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Blue (numpy.array, size=(m,n)) – blue band of satellite image
method –
- ‘reinhardt’
XYZitu601-1 axis
- ’ford’
D65 illuminant
- Returns
X,Y,Z
- Return type
numpy.array, size=(m,n)
References
- Re01
Reinhard et al. “Color transfer between images” IEEE Computer graphics and applications vol.21(5) pp.34-41, 2001.
- FR98
Ford & Roberts. “Color space conversion”, pp. 1–31, 1998.
- dhdt.preprocessing.color_transforms.rgb2ycbcr(Red, Green, Blue)¶
transform red, green, blue arrays to luna and chroma values
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Blue (numpy.array, size=(m,n)) – blue band of satellite image
- Returns
Y (numpy.array, size=(m,n)) – luma
Cb (numpy.array, size=(m,n)) – chroma
Cr (numpy.array, size=(m,n)) – chroma
References
- Ts06
Tsai, “A comparative study on shadow compensation of color aerial images in invariant color models”, IEEE transactions in geoscience and remote sensing, vol. 44(6) pp. 1661–1671, 2006.
- dhdt.preprocessing.color_transforms.rgb2yiq(Red, Green, Blue)¶
transform red, green, blue to luminance, inphase, quadrature values
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Blue (numpy.array, size=(m,n)) – blue band of satellite image
- Returns
Y (numpy.array, size=(m,n)) – luminance
I (numpy.array, size=(m,n)) – inphase
Q (numpy.array, size=(m,n)) – quadrature
References
- GW92
Gonzalez & Woods “Digital image processing”, 1992.
- dhdt.preprocessing.color_transforms.xyz2lab(X, Y, Z, th=0.008856)¶
transform XYZ tristimulus arrays to Lab values
- Parameters
X (numpy.array, size=(m,n)) –
Y (numpy.array, size=(m,n)) –
Z (numpy.array, size=(m,n)) –
- Returns
L,a,b
- Return type
numpy.array, size=(m,n)
References
- FR98
Ford & Roberts. “Color space conversion”, pp. 1–31, 1998.
- Si18
Silva et al. “Near real-time shadow detection and removal in aerial motion imagery application” ISPRS journal of photogrammetry and remote sensing, vol.140 pp.104–121, 2018.
- dhdt.preprocessing.color_transforms.xyz2lms(X, Y, Z)¶
transform XYZ tristimulus arrays to LMS values
- Parameters
X (numpy.array, size=(m,n)) – modified XYZitu601-1 axis
Y (numpy.array, size=(m,n)) – modified XYZitu601-1 axis
Z (numpy.array, size=(m,n)) – modified XYZitu601-1 axis
- Returns
L,M,S
- Return type
numpy.array, size=(m,n)
References
- Re01
Reinhard et al. “Color transfer between images” IEEE Computer graphics and applications vol.21(5) pp.34-41, 2001.
- dhdt.preprocessing.color_transforms.yiq2rgb(Y, Inph, Quadr)¶
transform luminance, inphase, quadrature values to red, green, blue
- Parameters
Y (numpy.array, size=(m,n)) – luminance
Inph (numpy.array, size=(m,n)) – inphase
Quadr (numpy.array, size=(m,n)) – quadrature
- Returns
Red (numpy.array, size=(m,n)) – red band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Blue (numpy.array, size=(m,n)) – blue band of satellite image
See also
References
- GW92
Gonzalez & Woods “Digital image processing”, 1992.
acquisition geometry¶
functions of use for geometric information extraction and retrival
- dhdt.preprocessing.acquisition_geometry.compensate_ortho_offset(Img, Z, dx, dy, obs_az, obs_zn, geoTransform)¶
- Parameters
Img (numpy.array, size=(m,n), ndim=2, dtype={float,integer}) – grid with intensities to be compensated
Z (numpy.array, size=(m,n), ndim=2, dtype={float,integer}, unit=meters) – grid with elevations
dx (float, unit=meter) – mis-registration offset
dy (float, unit=meter) – mis-registration offset
obs_az (numpy.array, size=(m,n), unit=degrees) – azimuth angle of observation
obs_zn (numpy.ndarray, size=(m,n), unit=degrees) – zenith angle of observation
geoTransform1 (tuple) – affine transformation coefficients of array Z
- Returns
Img_cor – corrected intensities
- Return type
numpy.array, size=(m,n), ndim=2, dtype=float
See also
- dhdt.preprocessing.acquisition_geometry.get_ortho_offset(Z, dx, dy, obs_az, obs_zn, geoTransform)¶
get terrain displacement due to miss-registration
- Parameters
Z (numpy.array, size=(m,n), unit=meter) – elevation model
dx (float, unit=meter) – mis-registration
dy (float, unit=meter) – mis-registration
obs_az (numpy.array, size=(m,n), unit=degrees) – azimuth angle of observation
obs_zn (numpy.ndarray, size=(m,n), unit=degrees) – zenith angle of observation
geoTransform (tuple, size={(6,), (8,)}) – affine transformation coefficients of array Z
- Returns
dI, dJ – terrain displacements
- Return type
numpy.ndarray, size=(m,n), unit=pixels
Notes
It is important to know what type of coordinate systems exist, hence:
coordinate | coordinate ^ y system 'ij'| system 'xy' | | | | j | x --------+--------> --------+--------> | | | | | i | v |
The azimuth angle is declared in the following coordinate frame:
^ North & y | - <--|--> + | +----> East & x
The angles related to the satellite are as follows:
#*# #*# satellite ^ / ^ /| | / | / | nadir |-- zenith angle | / v | / | /| |/ |/ | elevation angle +----- surface +------
- dhdt.preprocessing.acquisition_geometry.get_template_acquisition_angles(Az, Zn, Det, i_samp, j_samp, t_size)¶
- Parameters
Az (numpy.ndarray, size=(m,n), float) – array with sensor azimuth angles
Zn (numpy.ndarray, size=(m,n), float) – array with sensor zenith angles
Det (numpy.ndarray, size=(m,n), bool) – array with sensor detector ids
i_samp (numpy.ndarray, size=(k,l), integer) – array with collumn coordinates of the template centers
j_samp (numpy.ndarray, size=(k,l), integer) – array with row coordinates of the template centers
t_size (integer, {x ∈ ℕ | x ≥ 1}, unit=pixels) – size of the template
- Returns
Azimuth, Zenith – mean observation angles, that is, azimuth and zenith
- Return type
np.array, size=(k,l), float
Notes
The azimuth angle declared in the following coordinate frame:
^ North & y | - <--|--> + | +----> East & x
- dhdt.preprocessing.acquisition_geometry.get_template_aspect_slope(Z, i_samp, j_samp, t_size, spac=10.0)¶
- Parameters
Z (numpy.ndarray, size=(m,n), float, unit=meters) – array with elevation values
i_samp (numpy.ndarray, size=(k,l), integer, unit=pixels) – array with collumn coordinates of the template centers
j_samp (numpy.ndarray, size=(k,l), integer, unit=pixels) – array with row coordinates of the template centers
t_size (integer, {x ∈ ℕ | x ≥ 1}, unit=pixels) – size of the template
spac (float, unit=meters) – spacing of the elevation grid
- Returns
Slope, Aspect – mean slope and aspect angle in the template
- Return type
numpy.ndarray, size=(k,l), float
See also
get_aspect_slope
Notes
It is important to know what type of coordinate systems exist, hence:
coordinate | coordinate ^ y system 'ij'| system 'xy' | | | | j | x --------+--------> --------+--------> | | | | | i | v |
- dhdt.preprocessing.acquisition_geometry.slope_along_perp(Z, Az, spac=10)¶
get the slope of the terrain along a specific direction
- Parameters
Z (numpy.ndarray, size=(m,n), units=meters) – elevation model
Az ({numpy.ndarray, float}, units=degrees) – argument or azimuth angle of interest
spac (float, units=meters) – spacing used for the elevation model
- Returns
Slp_para (numpy.ndarray, size=(m,n), units=degrees) – slope in the along direction of the azimuth angle
Slp_perp (numpy.ndarray, size=(m,n), units=degrees) – slope in the perpendicular direction of the azimuth angle
aquatic transforms¶
functions of use for enhancement of water or suspended material
- dhdt.preprocessing.aquatic_transforms.alternative_floating_algae_index(Red, Rededge, Near)¶
transform near and short wave infrared arrays to get a algae-index. Based upon [Wa18]
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Rededge (numpy.ndarray, size=(m,n)) – rededge band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
AFAI – array with floating algae transform
- Return type
numpy.ndarray, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[AFAI=(B_{06} - B_{04})/(B_{8A} - B_{04})\]References
- Wa18
Wang, et al. “Remote sensing of Sargassum biomass, nutrients, and pigments” Geophysical Research Letters, vol.45(22) pp.12-359, 2018.
- dhdt.preprocessing.aquatic_transforms.atmospherically_resistant_vegetation_index(Blue, Red, Near, gamma=1.0, L=0.5)¶
Based upon [KT92]
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
ARVI – array with atmospherically resistant vegetation index
- Return type
numpy.ndarray, size=(m,n)
References
- KT92
Kaufman & Tanre, “Atmospherically resistant vegetation index (ARVI)” IEEE transactions in geoscience and remote sensing, vol.30 pp.261–270, 1992.
- dhdt.preprocessing.aquatic_transforms.floating_algae_index(Red, Near, Swir, df)¶
Following [Hu09]
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
Swir (numpy.ndarray, size=(m,n)) – shortwave infrared band of satellite image
df (datafram) –
- Returns
FAI – array with floating algae index
- Return type
numpy.ndarray, size=(m,n)
References
- Hu09
Hu, “A novel ocean color index to detect floating algae in the global oceans” Remote sensing of the enviroment, vol.113(10) pp.2188-2129, 2009
- dhdt.preprocessing.aquatic_transforms.infrared_percentage_vegetation_index(Red, Near)¶
adjustment over the NDVI, see [Cr90].
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
NPVI – array with infrared percentage vegetation index
- Return type
numpy.ndarray, size=(m,n), range=0…1
References
- Cr90
Crippen “Calculating the vegetation index faster” Remote sensing of environment, vol.34 pp.71-73, 1990.
- dhdt.preprocessing.aquatic_transforms.modified_chlorophyll_absorption_index(Green, Red, Rededge)¶
Based upon [Da00]
- Parameters
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Rededgde (numpy.ndarray, size=(m,n)) – rededge band of satellite image
- Returns
MCARI – array with atmospherically resistant vegetation index
- Return type
numpy.ndarray, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[MCARI= ((B_{05} - B_{04}) - .2(B_{05} - B_{03})) * (B_{04})/(B_{03})\]References
- Da00
Daughtry et al. “Estimating corn leaf chlorophyll concentration from leaf and canopy reflectance” Remote sensing of environment vol.74 pp.229–239, 2000.
- dhdt.preprocessing.aquatic_transforms.modified_simple_ratio(Red, Near)¶
Based upon [Ch96]
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
SR – array with simple ratio
- Return type
numpy.ndarray, size=(m,n)
References
- Ch96
Chen, “Evaluation of vegetation indices and modified simple ratio for boreal applications” Canadian journal of remote sensing, vol.22 pp.229–242, 1996.
- dhdt.preprocessing.aquatic_transforms.modified_soil_adjusted_vegetation_index(Red, Near)¶
Based upon [Qi94]
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
MSAVI – array with modified soil adjusted vegetation index
- Return type
numpy.ndarray, size=(m,n)
References
- Qi94
Qi et al., “A modified soil vegetation adjusted index” Remote sensing of environment, vol.48, pp.119–126, 1994.
- dhdt.preprocessing.aquatic_transforms.normalized_difference_vegetation_index(Red, Near)¶
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
NDVI – array with normalized difference vegetation index
- Return type
numpy.ndarray, size=(m,n)
References
- Ro74
Rouse et al. “Monitoring the vernal advancements and retrogradation of natural vegetation” in NASA/GSFC final report, pp.1–137, 1974.
- dhdt.preprocessing.aquatic_transforms.renormalized_difference_vegetation_index(Red, Near)¶
Based upon [RB95]
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
NDVI – array with renormalized difference vegetation index
- Return type
numpy.ndarray, size=(m,n)
References
- RB95
Rougean & Breon, “Estimating PAR absorbed by vegetation from bidirectional reflectance measurements” Remote sensing of environment, vol.51 pp.375–384, 1995.
- dhdt.preprocessing.aquatic_transforms.simple_ratio(Red, Near)¶
Based upon [Jo69]
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
SR – array with simple ratio
- Return type
numpy.ndarray, size=(m,n)
References
- Jo69
Jordan, “Derivation of leaf area index from quality of light on the forest floor” Ecology vol.50 pp.663–666, 1969.
- dhdt.preprocessing.aquatic_transforms.soil_adjusted_vegetation_index(Red, Near, L=0.5)¶
Based upon [Hu88]
- Parameters
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
L (float, default=.5) – multiplication factor
- Returns
SAVI – array with soil adjusted vegetation index
- Return type
numpy.ndarray, size=(m,n)
References
- Hu88
Huete, “A soil vegetation adjusted index (SAVI)” Remote sensing of environment, vol.25 pp.295–309, 1988.
- dhdt.preprocessing.aquatic_transforms.triangular_vegetation_index(Green, Red, Rededge)¶
Based upon [BL00]
- Parameters
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Rededgde (numpy.ndarray, size=(m,n)) – rededge band of satellite image
- Returns
TVI – array with triangular vegetation index
- Return type
numpy.ndarray, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[TVI= .5 (120(B_{06} - B_{03}) - 200(B_{04} - B_{03}))\]References
- BL00
Broge & Leblanc, “Comparing prediction power and stability of broadband and hyperspectral vegetation indices for estimation of green leaf area index and canopy chlorophyll density” Remote sensing of environment, vol.76 pp.156–172, 2000.
snow transforms¶
functions of use for enhancement of snow and ice
- dhdt.preprocessing.snow_transforms.automated_glacier_extraction_index(Red, Near, Short, alpha=0.5)¶
transform red, near and shortwave infrared arrays to AGEI, based on [Zh19].
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Near (numpy.array, size=(m,n)) – near infrared band of satellite image
Short (numpy.array, size=(m,n)) – shortwave infrared band of satellite image
alpha (float) –
- Returns
AGEI – array with AGEI transform
- Return type
numpy.array, size=(m,n)
See also
Notes
Based on the bands of Sentinel-2:
\[AGEI = (B_{3}-B_{8}) / (B_{3}+B_{8})\]References
- Zh19
Zhang et al. “Automated glacier extraction index by optimization of red/SWIR and NIR/SWIR ratio index for glacier mapping using Landsat imagery” Water, vol.11 pp.1233, 2019.
- dhdt.preprocessing.snow_transforms.normalized_difference_snow_ice_index(Red, Short)¶
transform red and short wave infrared arrays NDSII-index, see also [Xi01]
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Short (numpy.array, size=(m,n)) – shortwave infrared band of satellite image
- Returns
NDSII – array with NDSII transform
- Return type
numpy.array, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[NDSII = (B_{4}-B_{11}) / (B_{4}+B_{11})\]The bands of Sentinel-2 correspond to the following spectral range:
\(B_{4}\) : red
\(B_{11}\) : shortwave infrared
References
- Xi01
Xiao et al. “Assessing the potential of Vegetation sensor data for mapping snow and ice cover: A normalized snow and ice index” International journal of remote sensing, vol.22(13) pp.2479-2487, 2001.
- dhdt.preprocessing.snow_transforms.normalized_difference_snow_index(Green, Short)¶
transform red and short wave infrared arrays NDSII-index, see also [Do89].
- Parameters
Green (numpy.array, size=(m,n)) – green band of satellite image
Short (numpy.array, size=(m,n)) – shortwave infrared band of satellite image
- Returns
NDSI – array with NDSI transform
- Return type
numpy.array, size=(m,n)
Notes
- Based on the bands of Sentinel-2:
- \[NDSI = (B_{3}-B_{8}) / (B_{3}+B_{8})\]
- The bands of Sentinel-2 correspond to the following spectral range:
\(B_{3}\) : green
\(B_{8}\) : near infrared
See also
References
- Do89
Dozier. “Spectral signature of alpine snow cover from the Landsat Thematic Mapper” Remote sensing of environment, vol.28 pp.9-22, 1989.
- dhdt.preprocessing.snow_transforms.s3(Red, Near, Short)¶
transform red, near and shortwave infrared arrays to s3, see [SY99].
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Near (numpy.array, size=(m,n)) – near infrared band of satellite image
Short (numpy.array, size=(m,n)) – shortwave infrared band of satellite image
- Returns
S3 – array with S3 transform
- Return type
numpy.array, size=(m,n)
See also
Notes
Based on the bands of Sentinel-2:
\[S3 = B_{8}(B_{4}-B_{11}) / ((B_{4} + B_{8})(B_{4} + B_{11}))\]The bands of Sentinel-2 correspond to the following spectral range:
\(B_{4}\) : red
\(B_{8}\) : near infrared
\(B_{11}\) : shortwave infrared
See also
References
- SY99
Saito & Yamazaki. “Characteristics of spectral reflectance for vegetation ground surfaces with snowcover; Vegetation indeces and snow indices” Journal of Japan society of hydrology and water resources, vol.12(1) pp.28-38, 1999.
- dhdt.preprocessing.snow_transforms.snow_cover_fraction(Near, Short)¶
transform near and short wave infrared arrays to get a snow-index, from [Pa02].
- Parameters
Near (numpy.array, size=(m,n)) – near infrared band of satellite image
Short (numpy.array, size=(m,n)) – shortwave infrared band of satellite image
- Returns
NDSII – array with NDSII transform
- Return type
numpy.array, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[SCF= B_{8} / B_{11}\]The bands of Sentinel-2 correspond to the following spectral range:
\(B_{8}\) : near infrared
\(B_{11}\) : shortwave infrared
References
- Pa02
Paul et al. “The new remote-sensing-derived Swiss glacier inventory: I methods” Annuals of glaciology, vol.34 pp.355-361, 2002.
- dhdt.preprocessing.snow_transforms.snow_fraction(Red, Short)¶
transform red and short wave infrared arrays to get a snow-index, see also [PK05].
- Parameters
Red (numpy.array, size=(m,n)) – red band of satellite image
Short (numpy.array, size=(m,n)) – shortwave infrared band of satellite image
- Returns
SF – array with snow fraction transform
- Return type
numpy.array, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[SF = B_{4} / B_{11}\]The bands of Sentinel-2 correspond to the following spectral range:
\(B_{4}\) : red
\(B_{11}\) : shortwave infrared
References
- PK05
Paul & Kääb “Perspectives on the production of a glacier inventory from multispectral satellite data in Arctic Canada: Cumberland peninsula, Baffin Island” Annuals of glaciology, vol.42 pp.59-66, 2005.
- dhdt.preprocessing.snow_transforms.snow_water_index(Green, Near, Short)¶
transform green, near and shortwave infrared arrays to snow water index, as described in [Di19].
- Parameters
Green (numpy.array, size=(m,n)) – green band of satellite image
Near (numpy.array, size=(m,n)) – near infrared band of satellite image
Short (numpy.array, size=(m,n)) – shortwave infrared band of satellite image
- Returns
SWI – array with SWI transform
- Return type
numpy.array, size=(m,n)
See also
Notes
Based on the bands of Sentinel-2:
\[SWI = B_{8}(B_{3}-B_{11}) / ((B_{3} + B_{8})(B_{3} + B_{11}))\]The bands of Sentinel-2 correspond to the following spectral range:
\(B_{3}\) : green
\(B_{8}\) : near infrared
\(B_{11}\) : shortwave infrared
References
- Di19
Dixit et al. “Development and evaluation of a new snow water index (SWI) for accurate snow cover delineation” Remote sensing, vol.11(23) pp.2774, 2019.
shadow transforms¶
functions of use for enhancement of intensities, making use of spectral properties
- dhdt.preprocessing.shadow_transforms.apply_shadow_transform(method, Blue, Green, Red, RedEdge, Near, Shw, **kwargs)¶
Given a specific method, employ shadow transform
- Parameters
method ({‘siz’,’ isi’,’ sei’,’ fcsdi’, ‘nsvi’, ‘nsvdi’, ‘sil’, ‘sr’, ‘sdi’, ‘sisdi’, ‘mixed’, ‘c3’, ‘entropy’, ‘shi’,’ sp’, 'nri'}) –
method name to be implemented,can be one of the following:
’siz’ : shadow index first version
’isi’ : improved shadow index
’sei’ : shadow enhancement index
’fcsdi’ : false color shadow difference index
’csi’ : combinational shadow index
’nsvdi’ : normalized saturation value difference index
’sil’ : shadow index second
’sr’ : specthem ratio
’sdi’ : shadow detector index
’sisdi’ : saturation intensity shadow detector index
’mixed’ : mixed property based shadow index
’msf’ : modified shadow fraction
’c3’ : color invariant
’entropy’ : entropy shade removal
’shi’ : shade index
’sp’ : shadow probabilities
’nri’ : normalized range shadow index
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
RedEdge (numpy.ndarray, size=(m,n)) – red edge band of satellite image
Near (numpy.ndarray, size=(m,n)) – near-infrared band of satellite image
Shw (numpy.ndarray, size=(m,n)) – synthetic shading from for example a digital elevation model
- Returns
M – shadow enhanced satellite image
- Return type
numpy.ndarray, size=(m,n)
See also
handler_multispec.get_shadow_bands
provides Blue, Green, Red bandnumbers for a specific satellite
handler_multispec.read_shadow_bands
reads the specific bands given
- dhdt.preprocessing.shadow_transforms.color_invariant(Blue, Green, Red)¶
transform red, green, blue arrays to color invariant (c3), see also [GS99] and [AA08].
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Red (numpy.array, size=(m,n)) – red band of satellite image
- Returns
c3 – array with shadow enhancement
- Return type
numpy.array, size=(m,n)
References
- GS99
Gevers & Smeulders. “Color-based object recognition” Pattern recognition, vol.32 pp.453–464, 1999
- AA08
Arévalo González & Ambrosio. “Shadow detection in colour high‐resolution satellite images” International journal of remote sensing, vol.29(7) pp.1945–1963, 2008
- dhdt.preprocessing.shadow_transforms.combinational_shadow_index(Blue, Green, Red, Near)¶
transform red, green, blue arrays to combined shadow index. See also [Su19].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
CSI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
See also
References
- Su19
Sun et al. “Combinational shadow index for building shadow extraction in urban areas from Sentinel-2A MSI imagery” International journal of applied earth observation and geoinformation, vol.78 pp.53–65, 2019.
- dhdt.preprocessing.shadow_transforms.entropy_shade_removal(Ia, Ib, Ic, a=None)¶
Make use of Wiens’ law to get illumination component, see also [Fi09].
- Parameters
Ia (numpy.array, size=(m,n)) – highest frequency band of satellite image
Ib (numpy.array, size=(m,n)) – lower frequency band of satellite image
Ic (numpy.array, size=(m,n)) – lowest frequency band of satellite image
a (float, default=-1, unit=degrees) – angle of the coordinate rotation
- Returns
S – shadow enhanced band
- Return type
numpy.array, size=(m,n)
References
- Fi09
Finlayson et al. “Entropy minimization for shadow removal” International journal of computer vision vol.85(1) pp.35-57 2009.
- dhdt.preprocessing.shadow_transforms.false_color_shadow_difference_index(Green, Red, Near)¶
transform red and near infrared arrays to shadow index. See also [Te11].
- Parameters
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
FCSI – array with shadow enhancement
- Return type
numpy.ndarray, size=(m,n)
See also
Notes
\[FCSI = (I_{sat} - I_{int}) / (I_{sat} + I_{int})\]References
- Te11
Teke et al. “Multi-spectral false color shadow detection”, Proceedings of the ISPRS conference on photogrammetric image analysis, pp.109-119, 2011.
- dhdt.preprocessing.shadow_transforms.fractional_range_shadow_index(*args)¶
transform multi-spectral arrays to shadow index
- Parameters
*args (numpy.array, size=(m,n)) – different bands of a satellite image, all of the same extent
- Returns
SI – shadow enhanced band
- Return type
numpy.array, size=(m,n)
Notes
Schematically this looks like:
I1 ┌-----┐ -┬---┤ | I2|┌--┤ | Imin -┼┼┬-┤ min ├------┐ |||┌┤ | |┌-------┐ ||||└-----┘ └┤ a / b | SI I3||||┌-----┐ ┌┤ ├------ -┼┴┼┼┤ | Imax |└-------┘ I4└-┼┼┤ ├------┘ ----┼┴┤ max | └-┤ | └-----┘
- dhdt.preprocessing.shadow_transforms.improved_shadow_index(Blue, Green, Red, Near)¶
transform red, green, blue arrays to improved shadow index. Based upon [Zh21].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
ISI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
See also
Notes
\[ISI = (SI + (1-I_{near}))/(SI + (1+I_{near}))\]References
- Zh21
Zhou et al. “Shadow detection and compensation from remote sensing images under complex urban conditions”, 2021.
- dhdt.preprocessing.shadow_transforms.mixed_property_based_shadow_index(Blue, Green, Red)¶
transform red, green, blue arrays to Mixed property-based shadow index, see also [Ha18].
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Red (numpy.array, size=(m,n)) – red band of satellite image
- Returns
MPSI – array with shadow enhancement
- Return type
numpy.array, size=(m,n)
References
- Ha18
Han et al. “A mixed property-based automatic shadow detection approach for VHR multispectral remote sensing images” Applied sciences vol.8(10) pp.1883, 2018.
- dhdt.preprocessing.shadow_transforms.modified_color_invariant(Blue, Green, Red, Near)¶
transform red, green, blue arrays to color invariant (c3), see also [GS99] and [BA15].
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Red (numpy.array, size=(m,n)) – red band of satellite image
Near (numpy.array, size=(m,n)) – near infrared band of satellite image
- Returns
c3 – array with shadow enhancement
- Return type
numpy.array, size=(m,n)
References
- GS99
Gevers & Smeulders. “Color-based object recognition” Pattern recognition, vol.32 pp.453–464, 1999
- BA15
Besheer & Abdelhafiz. “Modified invariant colour model for shadow detection” International journal of remote sensing, vol.36(24) pp.6214–6223, 2015.
- dhdt.preprocessing.shadow_transforms.modified_normalized_difference_water_index(Green, Shortwave)¶
transform green and shortwave infrared arrays to modified NDW-index. See also [Xu06].
- Parameters
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Shortwave (numpy.ndarray, size=(m,n)) – shortwave infrared band of satellite image, for Sentinel-2 this is B11, originally the Landsat TM band 5 was used in [Xu06].
- Returns
NDWI – array with NDWI transform
- Return type
numpy.ndarray, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[NDWI = [B_{3}-B_{11}]/[B_{3}+B_{11}]\]References
- dhdt.preprocessing.shadow_transforms.modified_shadow_fraction(Blue, Green, Red, P_S=0.95)¶
transform red, green, blue arrays to shadow index. Based upon [Ch09].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
P_S (float, range=0...1) – threshold value for classification, based on the quantile
- Returns
SI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
See also
References
- Ch09
Chung, et al. “Efficient shadow detection of color aerial images based on successive thresholding scheme.” IEEE transactions on geoscience and remote sensing vol.47, pp.671–682, 2009.
- dhdt.preprocessing.shadow_transforms.normalized_difference_blue_water_index(Blue, Near)¶
transform green and near infrared arrays NDW-index. See also [Qu11].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
NDWI – array with NDWI transform
- Return type
numpy.ndarray, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[NDWI = [B_{3}-B_{8}]/[B_{3}+B_{8}]\]References
- Qu11
Qu et al, “Research on automatic extraction of water bodies and wetlands on HJU satellite CCD images” Remote sensing information, vol.4 pp.28–33, 2011
- dhdt.preprocessing.shadow_transforms.normalized_difference_moisture_index(Near, Shortwave)¶
transform near and shortwave infrared arrays to NDM-index. See also [Wi02].
- Parameters
Near (numpy.ndarray, size=(m,n)) – green band of satellite image
Shortwave (numpy.ndarray, size=(m,n)) – shortwave infrared band of satellite image
- Returns
NDMI – array with NDMI transform
- Return type
numpy.ndarray, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[NDMI = [B_{8}-B_{11}]/[B_{8}+B_{11}]\]References
- Wi02
Wilson, “Detection of forest harvest type using multiple dates of Landsat TM imagery” Remote sensing of environment, vol.80 pp.385–396, 2002.
- dhdt.preprocessing.shadow_transforms.normalized_difference_water_index(Green, Near)¶
transform green and near infrared arrays NDW-index. See also [Ga96].
- Parameters
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
NDWI – array with NDWI transform
- Return type
numpy.ndarray, size=(m,n)
Notes
Based on the bands of Sentinel-2:
\[NDWI = [B_{3}-B_{8}]/[B_{3}+B_{8}]\]References
- Ga96
Gao, “NDWI - a normalized difference water index for remote sensing of vegetation liquid water from space” Remote sensing of environonment, vol.58 pp.257–266, 1996
- dhdt.preprocessing.shadow_transforms.normalized_range_shadow_index(*args)¶
transform multi-spectral arrays to shadow index
- Parameters
*args (numpy.array, size=(m,n)) – different bands of a satellite image, all of the same extent
- Returns
SI – shadow enhanced band
- Return type
numpy.array, size=(m,n)
Notes
Schematically this looks like:
I1 ┌-----┐ -┬---┤ | I2|┌--┤ | Imin -┼┼┬-┤ min ├------┐ |||┌┤ | |┌-------┐ ||||└-----┘ └┤(a-b)/ | SI I3||||┌-----┐ ┌┤ (a+b)├------ -┼┴┼┼┤ | Imax |└-------┘ I4└-┼┼┤ ├------┘ ----┼┴┤ max | └-┤ | └-----┘
- dhdt.preprocessing.shadow_transforms.normalized_sat_value_difference_index(Blue, Green, Red)¶
transform red, green, blue arrays to normalized sat-value difference. See also [Ma08].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
NSVDI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
See also
false_color_shadow_difference_index
,rgb2hsi
Notes
\[NSVDI = (I_{sat} - I_{int}) / (I_{sat} + I_{int})\]References
- Ma08
Ma et al. “Shadow segmentation and compensation in high resolution satellite images”, Proceedings of IEEE IGARSS, pp.II-1036–II-1039, 2008.
- dhdt.preprocessing.shadow_transforms.reinhard(Blue, Green, Red)¶
transform red, green, blue arrays to luminance, see also [Re01].
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Red (numpy.array, size=(m,n)) – red band of satellite image
- Returns
l – shadow enhanced band
- Return type
numpy.array, size=(m,n)
References
- Re01
Reinhard et al. “Color transfer between images” IEEE Computer graphics and applications vol.21(5) pp.34-41, 2001.
- dhdt.preprocessing.shadow_transforms.sat_int_shadow_detector_index(Blue, RedEdge, Near)¶
transform red, green, blue arrays to sat-int shadow detector index, see also [MA18].
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
RedEdge (numpy.array, size=(m,n)) – rededge band of satellite image
Near (numpy.array, size=(m,n)) – near infrared band of satellite image
- Returns
SISDI – array with shadow enhancement
- Return type
numpy.array, size=(m,n)
References
- MA18
Mustafa & Abedelwahab. “Corresponding regions for shadow restoration in satellite high-resolution images” International journal of remote sensing, vol.39(20) pp.7014–7028, 2018.
- dhdt.preprocessing.shadow_transforms.shade_index(*args)¶
transform multi-spectral arrays to shade index, see also [HS10] and [Al12].
- Parameters
*args (numpy.array, size=(m,n)) – different bands of a satellite image, all of the same extent
- Returns
SI – shadow enhanced band
- Return type
numpy.array, size=(m,n)
Notes
Amplify dark regions by simple multiplication:
\[SI = Pi_{i}^{N}{[1 - B_{i}]}\]References
- HS10
Hogan & Smith, “Refinement of digital elvation models from shadowing cues” IEEE computer society conference on computer vision and pattern recognition, 2010.
- Al12
Altena, “Filling the white gap on the map: Photoclinometry for glacier elevation modelling” MSc thesis TU Delft, 2012.
- dhdt.preprocessing.shadow_transforms.shadow_detector_index(Blue, Green, Red)¶
transform red, green, blue arrays to shadow detector index, see [MA17].
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Red (numpy.array, size=(m,n)) – red band of satellite image
- Returns
SDI – array with shadow enhancement
- Return type
numpy.array, size=(m,n)
References
- MA17
Mustafa & Abedehafez. “Accurate shadow detection from high-resolution satellite images” IEEE geoscience and remote sensing letters, vol.14(4) pp.494–498, 2017.
- dhdt.preprocessing.shadow_transforms.shadow_enhancement_index(Blue, Green, Red, Near)¶
transform red, green, blue arrays to SEI index. See also [Su19].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
Near (numpy.ndarray, size=(m,n)) – near infrared band of satellite image
- Returns
SEI – array with sei transform
- Return type
numpy.ndarray, size=(m,n)
See also
Notes
Based on the bands of Sentinel-2:
\[SEI = [(B_{1}+B_{9})-(B_{3}+B_{8})]/[(B_{1}+B_{9})+(B_{3}+B_{8})]\]References
- Su19
Sun et al. “Combinational shadow index for building shadow extraction in urban areas from Sentinel-2A MSI imagery” International journal of applied earth observation and geoinformation, vol.78 pp.53–65, 2019.
- dhdt.preprocessing.shadow_transforms.shadow_hcv_fraction(Blue, Green, Red)¶
transform red, green, blue arrays to shadow index. Based upon [Ts06].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
SI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
References
- Ts06
Tsai. “A comparative study on shadow compensation of color aerial images in invariant color models.” IEEE Transactions on geoscience and remote sensing vol.44 pp.1661–1671, 2006.
- dhdt.preprocessing.shadow_transforms.shadow_hsv_fraction(Blue, Green, Red)¶
transform red, green, blue arrays to shadow index. Based upon [Ts06].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
SI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
References
- Ts06
Tsai. “A comparative study on shadow compensation of color aerial images in invariant color models.” IEEE Transactions on geoscience and remote sensing vol.44 pp.1661–1671, 2006.
- dhdt.preprocessing.shadow_transforms.shadow_identification(Blue, Green, Red)¶
transform red, green, blue arrays to sat-value difference, following [Po03].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
SI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
See also
false_color_shadow_difference_index
,rgb2hsi
Notes
\[SI = I_{sat} - I_{int}\]References
- Po03
Polidorio et al. “Automatic shadow segmentation in aerial color images”, Proceedings of the 16th Brazilian symposium on computer graphics and image processing, pp.270-277, 2003.
- dhdt.preprocessing.shadow_transforms.shadow_index_liu(Blue, Green, Red)¶
transform red, green, blue arrays to shadow index, see [LX13] and [Su19].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
SI – array with shadow enhancement
- Return type
numpy.ndarray, size=(m,n)
References
- LX13
Liu & Xie. “Study on shadow detection in high resolution remote sensing image of PCA and HIS model”, Remote sensing technology and application, vol.28(1) pp. 78–84, 2013.
- Su19
Sun et al. “Combinational shadow index for building shadow extraction in urban areas from Sentinel-2A MSI imagery” International journal of applied earth observation and geoinformation, vol.78 pp.53–65, 2019.
- dhdt.preprocessing.shadow_transforms.shadow_index_zhou(Blue, Green, Red)¶
transform red, green, blue arrays to shadow index, based upon [Zh21].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
SI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
See also
References
- Zh21
Zhou et al. “Shadow detection and compensation from remote sensing images under complex urban conditions”, 2021.
- dhdt.preprocessing.shadow_transforms.shadow_probabilities(Blue, Green, Red, Near, **kwargs)¶
transform blue, green, red and near infrared to shade probabilities, see also [FS10] and [Rü13].
- Parameters
Blue (numpy.array, size=(m,n)) – blue band of satellite image
Green (numpy.array, size=(m,n)) – green band of satellite image
Red (numpy.array, size=(m,n)) – red band of satellite image
Near (numpy.array, size=(m,n)) – near infrared band of satellite image
- Returns
M – shadow enhanced band
- Return type
numpy.array, size=(m,n)
References
- FS10
Fredembach & Süsstrunk. “Automatic and accurate shadow detection from (potentially) a single image using near-infrared information” EPFL Tech Report 165527, 2010.
- Rü13
Rüfenacht et al. “Automatic and accurate shadow detection using near-infrared information” IEEE transactions on pattern analysis and machine intelligence, vol.36(8) pp.1672-1678, 2013.
- dhdt.preprocessing.shadow_transforms.shadow_ycbcr_fraction(Blue, Green, Red)¶
transform red, green, blue arrays to shadow index. Based upon [Ts06].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
SI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
References
- Ts06
Tsai. “A comparative study on shadow compensation of color aerial images in invariant color models.” IEEE Transactions on geoscience and remote sensing vol.44 pp.1661–1671, 2006.
- dhdt.preprocessing.shadow_transforms.shadow_yiq_fraction(Blue, Green, Red)¶
transform red, green, blue arrays to shadow index. Based upon [Ts06].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
SI – array with shadow transform
- Return type
numpy.ndarray, size=(m,n)
See also
References
- Ts06
Tsai. “A comparative study on shadow compensation of color aerial images in invariant color models.” IEEE Transactions on geoscience and remote sensing vol.44 pp.1661–1671, 2006.
- dhdt.preprocessing.shadow_transforms.shannon_entropy(X, band_width=100)¶
calculate shanon entropy from data sample
- Parameters
X (numpy.ndarray, size=(m,1)) – data array
band_width (integer, optional) – sample bandwith
- Returns
shannon – entropy value of the total data sample
- Return type
float
- dhdt.preprocessing.shadow_transforms.specthem_ratio(Blue, Green, Red)¶
transform red, green, blue arrays to specthem ratio, see [Si18].
- Parameters
Blue (numpy.ndarray, size=(m,n)) – blue band of satellite image
Green (numpy.ndarray, size=(m,n)) – green band of satellite image
Red (numpy.ndarray, size=(m,n)) – red band of satellite image
- Returns
Sr – array with shadow enhancement
- Return type
numpy.ndarray, size=(m,n)
Notes
In the paper the logarithmic is also used further to enhance the shadows,
\[Sr = log{[Sr+1]}\]References
- Si18
Silva et al. “Near real-time shadow detection and removal in aerial motion imagery application” ISPRS journal of photogrammetry and remote sensing, vol.140 pp.104–121, 2018.
shadow filters¶
functions of use for enhancement of intensities, making use of spatial properties
- dhdt.preprocessing.shadow_filters.L0_smoothing(Z, lamb=0.02, kappa=2.0, beta_max=100000.0)¶
- Parameters
Z (numpy.ndarray, size={(m,n), (m,n,b)}, dtype=float) – grid with intensity values
lamb (float, default=.002) –
kappa (float, default=2.) –
beta_max (float, default=10000) –
- Returns
Z – modified intensity image
- Return type
numpy.ndarray, size={(m,n), (m,n,b)}, dtype=float
References
- Xu11
Xu et al. “Image smoothing via L0 gradient minimization” ACM transactions on graphics, 2011.
- dhdt.preprocessing.shadow_filters.anistropic_diffusion_scalar(Z, iter=10, K=0.15, s=0.25, n=4)¶
non-linear anistropic diffusion filter of a scalar field
- Parameters
Z (np.array, size=(m,n), dtype={float,complex}, ndim={2,3}) – intensity array
iter (integer, {x ∈ ℕ | x ≥ 0}) – amount of iterations
K (float, default=.15) – parameter based upon the noise level
s (float, default=.25) – parameter for the integration constant
n ({4,8}) – amount of neighbors used
- Returns
I
- Return type
np.array, size=(m,n)
References
- PM87
Perona & Malik “Scale space and edge detection using anisotropic diffusion” Proceedings of the IEEE workshop on computer vision, pp.16-22, 1987
- Ge92
Gerig et al. “Nonlinear anisotropic filtering of MRI data” IEEE transactions on medical imaging, vol.11(2), pp.221-232, 1992
- dhdt.preprocessing.shadow_filters.diffusion_strength_1(Z, K)¶
first diffusion function, proposed by [PM87], when a complex array is provided the absolute magnitude is used, following [Ge92]
- Parameters
Z (np.array, size=(m,n), dtype={float,complex}, ndim={2,3}) – array with intensities
K (float) – parameter based upon the noise level
- Returns
g – array with diffusion parameters
- Return type
np.array, size=(m,n), dtype=float
References
- PM87
Perona & Malik “Scale space and edge detection using anisotropic diffusion” Proceedings of the IEEE workshop on computer vision, pp.16-22, 1987
- Ge92
Gerig et al. “Nonlinear anisotropic filtering of MRI data” IEEE transactions on medical imaging, vol.11(2), pp.221-232, 1992
- dhdt.preprocessing.shadow_filters.diffusion_strength_2(Z, K)¶
second diffusion function, proposed by [PM87], when a complex array is provided the absolute magnitude is used, following [Ge82].
- Parameters
Z (np.array, size=(m,n), dtype={float,complex}) – array with intensities
K (float) – parameter based upon the noise level
- Returns
g – array with diffusion parameters
- Return type
np.array, size=(m,n), dtype=float
References
- PM87
Perona & Malik “Scale space and edge detection using anisotropic diffusion” Proceedings of the IEEE workshop on computer vision, pp.16-22, 1987.
- Ge92
Gerig et al. “Nonlinear anisotropic filtering of MRI data” IEEE transactions on medical imaging, vol.11(2), pp.221-232, 1992.
- dhdt.preprocessing.shadow_filters.enhance_shadows(Shw, method, **kwargs)¶
Given a specific method, employ shadow transform
- Parameters
Shw (np.array, size=(m,n), dtype={float,integer}) – array with intensities of shading and shadowing
method ({‘mean’,’kuwahara’,’median’,’otsu’,'anistropic'}) –
method name to be implemented,can be one of the following:
’mean’ : mean shift filter
’kuwahara’ : kuwahara filter
’median’ : iterative median filter
’anistropic’ : anistropic diffusion filter
- Returns
M – shadow enhanced image, done through the given method
- Return type
np.array, size=(m,n), dtype={float,integer}
- dhdt.preprocessing.shadow_filters.fade_shadow_cast(Shw, az, t_size=9)¶
- Parameters
Shw (np.array, size=(m,n)) – image with mask of shadows
az (float, unit=degrees) – illumination orientation
t_size (integer, {x ∈ ℕ | x ≥ 1}) – buffer size of the Gaussian blur
- Returns
Shw – image with faded shadows on the occluding end
- Return type
np.array, size=(m,n), dtype=float, range=0…1
- dhdt.preprocessing.shadow_filters.iterative_median_filter(Z, tsize=5, loop=50)¶
Transform intensity to more clustered intensity, through iterative filtering with a median operation
- Parameters
Z (np.array, size=(m,n), dtype=float) – array with intensity values
tsize (integer, {x ∈ ℕ | x ≥ 1}) – dimension of the kernel
loop (integer, {x ∈ ℕ | x ≥ 0}) – amount of iterations
- Returns
I_new – array with stark edges
- Return type
np.array, size=(m,n), dtype=float
See also
- dhdt.preprocessing.shadow_filters.kuwahara_filter(Z, tsize=5)¶
Transform intensity to more clustered intensity
- Parameters
Z (np.array, size=(m,n), dtype=float) – array with intensity values
tsize (integer, {x ∈ ℕ | x ≥ 1}) – dimension of the kernel
- Returns
Z_new
- Return type
np.array, size=(m,n), dtype=float
Notes
The template is subdivided into four blocks, where the sampled mean and standard deviation are calculated. Then mean of the block with the lowest variance is assigned to the central pixel. The configuration of the blocks look like:
┌-----┬-┬-----┐ ^ | A | | B | | ├-----┼-┼-----┤ | tsize ├-----┼-┼-----┤ | | C | | D | | └-----┴-┴-----┘ v
See also
- dhdt.preprocessing.shadow_filters.mean_shift_filter(Z, quantile=0.1)¶
Transform intensity to more clustered intensity, through mean-shift
- Parameters
Z (np.array, size=(m,n), dtype=float) – array with intensity values
quantile (float, range=0...1) –
- Returns
labels – array with numbered labels
- Return type
np.array, size=(m,n), dtype=integer
composition and geometry of the atmosphere¶
functions of use for describing the refractive properties of the atmosphere
- dhdt.preprocessing.atmospheric_geometry.ciddor_eq5(ρ_a, ρ_w, n_axs, n_ws, ρ_axs, ρ_ws)¶
estimating the refractive index through summation
- Parameters
ρ_a ({numpy.array, float}, unit=kg m-3) – density of dry air
ρ_w ((numpy.array, float}, unit=kg m-3) – density of moist air
n_axs (numpy.array, size=(k,)) – refractive index for dry air at sea level at a specific wavelength
n_ws (numpy.array, size=(k,)) – refractive index for standard moist air at sea level for a specific wavelength
ρ_axs (float, unit=kg m-3) – density of standard dry air at sea level
ρ_ws (float, unit=kg m-3) – density of standard moist air at sea level
- Returns
n_prop – refractive index
- Return type
float
References
- Ow67
Owens, “Optical refractive index of air: Dependence on pressure, temperature and composition”, Applied optics, vol.6(1) pp.51-59, 1967.
- dhdt.preprocessing.atmospheric_geometry.ciddor_eq6(ρ_a, ρ_w, n_axs, n_ws, ρ_axs, ρ_ws)¶
estimating the refractive index through the Lorenz-Lorentz relation
- Parameters
ρ_a (float, unit=kg m-3) – density of dry air
ρ_w (float, unit=kg m-3) –
n_axs (float) – refractive index for dry air at sea level
n_ws (float) – refractive index for standard moist air at sea level
ρ_axs (float, unit=kg m-3) – density of standard dry air at sea level
ρ_ws (float, unit=kg m-3) – density of standard moist air at sea level
- Returns
n_LL – refractive index
- Return type
float
References
- Ow67
Owens, “Optical refractive index of air: Dependence on pressure, temperature and composition”, Applied optics, vol.6(1) pp.51-59, 1967.
- Ci96
Ciddor, “Refractive index of air: new equations for the visible and near infrared”, Applied optics, vol.35(9) pp.1566-1573, 1996.
- dhdt.preprocessing.atmospheric_geometry.compressability_moist_air(p, T, x_w)¶
- Parameters
p (float, unit=Pascal) – total pressure
T (float, unit=Kelvin) – temperature
x_w (float) – molar fraction of water vapor in moist air
- Returns
Z – compressability factor
- Return type
float
See also
References
- Ci96
Ciddor, “Refractive index of air: new equations for the visible and near infrared”, Applied optics, vol.35(9) pp.1566-1573, 1996.
- DaXX
Davis, “Equation for the determination of the density of moist air (1981/91)” Metrologica, vol.29 pp.67-70.
- dhdt.preprocessing.atmospheric_geometry.get_CO2_level(lat, date, m=3, nb=3)¶
Rough estimation of CO₂, based upon measurements of global ‘Station data’__
The annual trend is based upon the ‘Global trend’__ in CO₂ while the seasonal trend is based upon relations as shown in Figure 1 in [Ba16]. Where a latitudal component is present in the ampltitude, as well as, an asymmetric triangular wave.
- Parameters
lat (float, unit=degrees, range=-90...+90) – latitude of the location of interest
date (np.datetime64, range=1950...2030) – time of interest
m (float, default=3) – assymetry of the wave - m<=2: crest is leftwards leaning - m==2: symmetric wave - m>=2: the ramp is leftwards leaning
nb (integer, default=3) – number of base functions to use, higher numbers will give a more pointy wave function
- Returns
total_co2 – estimate of spatial temporal carbon dioxide (CO₂) value
- Return type
float, unit=ppm
Notes
Data where this fitting is based on, can be found here:
The following acronyms are used:
CO2 : carbon dioxide
ppm : parts per million
Examples
>>> import numpy as np >>> import matplotlib.pyplot as plt
Plot the annual trend at 60ºN
>>> t = np.arange("1950-01-01","2022-01-01",dtype="M8[Y]") >>> co2 = get_CO2_level(60., t) >>> plt.plot(t, co2)
Plot seasonal signal
>>> t = np.arange("2018-01-01","2022-01-01",dtype="M8[D]") >>> co2 = get_CO2_level(60., t) >>> plt.plot(t, co2)
References
- Ba16
Barnes et al. “Isentropic transport and the seasonal cycle amplitude of CO₂” Journal of geophysical research: atmosphere, vol.121 pp.8106-8124, 2016.
- dhdt.preprocessing.atmospheric_geometry.get_T_in_troposphere(h, T_0)¶
Using parameters of the first layer of the eight layer atmospheric model of ISO 2533:1970 to estimate the temperature at different altitudes
- Parameters
h ({float, numpy.array}, unit=m) – altitude above sea-level
T_0 (float, unit=Kelvin) – temperature at the surface (reference)
- Returns
T_h – temperature in the troposphere
- Return type
{float, numpy.array}, unit=Kelvin
References
- Ya16
Yan et al, “Correction of atmospheric refraction geolocation error of high resolution optical satellite pushbroom images” Photogrammetric engineering & remote sensing, vol.82(6) pp.427-435, 2016.
- dhdt.preprocessing.atmospheric_geometry.get_T_sealevel(lat, method='noerdlinger')¶
The temperature at sea level can be related to its latitude [Al73] and a least squares fit is given in [No99] and [Ya16].
- Parameters
lat ({float, numpy.ndarray}, unit=degrees, range=-90...+90) – latitude of interest
- Returns
T_sl ({float, numpy.array}, unit=Kelvin) – estimate of temperature at sea level
method ({‘noerdlinger’, ‘yan’}) –
‘noerdlinger’, from [2] has no dependency between hemispheres
’yan’, from [3] is hemisphere dependent
References
- Al73
Allen, “Astrophysical quantities”, pp.121-124, 1973.
- No99
Noerdlinger, “Atmospheric refraction effects in Earth remote sensing” ISPRS journal of photogrammetry and remote sensing, vol.54, pp.360-373, 1999.
- Ya16
Yan et al, “Correction of atmospheric refraction geolocation error of high resolution optical satellite pushbroom images” Photogrammetric engineering & remote sensing, vol.82(6) pp.427-435, 2016.
- dhdt.preprocessing.atmospheric_geometry.get_density_air(T, P, fH, moist=True, CO2=450, M_w=0.018015, R=8.31451, saturation=False)¶
- Parameters
T (float, unit=Kelvin) – temperature at location
P (float, unit=Pascal) – the total pressure in Pascals
fH (float, range=0...1) – fractional humidity
moist (bool) – estimate moist (True) or dry air (False)
CO2 (float, unit=ppm of CO₂) – parts per million of CO₂ in the atmosphere
M_w (float, unit=kg mol-1) – molar mass of water vapor
R (float, unit=J mol-1 K-1) – universal gas constant
- Returns
ρ – density of moist air
- Return type
float, unit=kg m3
Notes
The following parameters are used in this function:
M_a : unit=kg mol-1, molar mass of dry air
x_w : molecular fraction of water vapor
Z : compressibility of moist air
CO2 : unit=ppm, carbon dioxide
ppm : parts per million volume
References
- Ci96
Ciddor, “Refractive index of air: new equations for the visible and near infrared”, Applied optics, vol.35(9) pp.1566-1573, 1996.
- Da92
Davis, “Equation for the determination of the density of moist air (1981/91)” Metrologica, vol.29 pp.67-70, 1992.
- dhdt.preprocessing.atmospheric_geometry.get_density_fraction(φ, h_0, alpha=0.0065, R=8314.36, M_t=28.825)¶
- Parameters
ϕ (unit=degrees) – latitude of point of interest
h_0 (float, unit=meters) – initial altitude above geoid
alpha (float, unit=K km-1, default=0.0065) – lapse rate in the Troposphere
R (float, unit=J kmol-1 K-1, default=8314.36) – universal gas constant
M_t (float, unit=kg kmol-1, default=28.825) – mean moleculear mass in the Troposphere
- Returns
ρ_frac – fraction of density in relation to density as sea-level
- Return type
float
Notes
The following parameters are used
T_sl : unit=Kelvin, temperature at sea level T_t : unit=Kelvin, temperature at the Tropopause h_t : unit=meters, altitude of the Tropopause g : unit=m s-2, acceleration due to gravity
References
- No99
Noerdlinger, “Atmospheric refraction effects in Earth remote sensing” ISPRS journal of photogrammetry and remote sensing, vol.54, pp.360-373, 1999.
- dhdt.preprocessing.atmospheric_geometry.get_height_tropopause(φ)¶
The tropopause can be related to its latitude [Al73] and a least squares fit is given in [No99].
- Parameters
ϕ ({float, numpy.ndarray}, unit=degrees, range=-90...+90) – latitude of interest
- Returns
h_t – estimate of tropopause, that is the layer between the troposphere and the stratosphere
- Return type
{float, numpy.ndarray}, unit=meters
See also
get_mean_lat
Notes
The Strato- and Troposphere have different refractive behaviour. A simplified view of this can be like:
z_s Mesosphere / ┌-----------------/-┐ ↑ 80 km | * T_s |* | / | | | * |* | Stratosphere / | | | * |* | / | | | * |* ├--Tropopause-/-----┤ | ↑ 11 km |_*_________T_t |*________ | | | | | | ** | * | Tropo | | | | | ** | ** | sphere | z_h | | | | ** T_h | **** └------------#------┘ ↓ ↓ └-----------T_0 └--------- Earth surface Temperature Pressure
References
- Al73
Allen, “Astrophysical quantities”, pp.121-124, 1973.
- No99
Noerdlinger, “Atmospheric refraction effects in Earth remote sensing” ISPRS journal of photogrammetry and remote sensing, vol.54, pp.360-373, 1999.
- dhdt.preprocessing.atmospheric_geometry.get_sat_vapor_press(T, method='IAPWS')¶
uses either [1] or the international standard of IAPWS [2].
- Parameters
T ({float, numpy.array}, unit=Kelvin) – temperature
- Returns
svp – saturation vapor pressure
- Return type
{float, numpy.array}
References
- Gi82
Giacomo, “Equation for the determination of the density of moist air” Metrologica, vol.18, pp.33-40, 1982.
- WP02
Wagner and Pruß, “The IAPWS formulation 1995 for the thermodynamic properties of ordinary water substance for general and scientific use.” Journal of physical and chemical reference data, vol.31(2), pp.387-535, 2002.
- dhdt.preprocessing.atmospheric_geometry.get_sat_vapor_press_ice(T)¶
following the international standard of IAPWS [2].
- Parameters
T ({float, numpy.array}, unit=Kelvin) – temperature
- Returns
svp – saturation vapor pressure
- Return type
{float, numpy.array}
References
- WP02
Wagner and Pruß, “The IAPWS formulation 1995 for the thermodynamic properties of ordinary water substance for general and scientific use.” Journal of physical and chemical reference data, vol.31(2), pp.387-535, 2002.
- dhdt.preprocessing.atmospheric_geometry.get_simple_density_fraction(φ)¶
The density fraction is related to latitude through a least squares fit as is given in [No99].
- Parameters
ϕ ({float, numpy.array}, unit=degrees, range=-90...+90) – latitude of interest
- Returns
densFrac – estimate of the density fraction
- Return type
{float, numpy.array}
See also
get_mean_lat
References
- No99
Noerdlinger, “Atmospheric refraction effects in Earth remote sensing” ISPRS journal of photogrammetry and remote sensing, vol.54, pp.360-373, 1999.
- dhdt.preprocessing.atmospheric_geometry.get_water_vapor_enhancement(t, p)¶
estimate the enhancement factor for the humidity calculation
- Parameters
t ({float, numpy.array}, unit=Celsius) – temperature
p ({float, numpy.array}, unit=Pascal) – atmospheric pressure
- Returns
f – enhancement factor
- Return type
{float, numpy.array}
References
- Gi82
Giacomo, “Equation for the determination of the density of moist air” Metrologica, vol.18, pp.33-40, 1982.
- dhdt.preprocessing.atmospheric_geometry.get_water_vapor_press(T)¶
- Parameters
t ({float, numpy.array}, unit=Kelvin) – temperature
- Returns
P_w – water vapor pressure
- Return type
{float, numpy.array}, unit=Pascal
See also
References
- Ya16
Yan et al, “Correction of atmospheric refraction geolocation error of high resolution optical satellite pushbroom images” Photogrammetric engineering & remote sensing, vol.82(6) pp.427-435, 2016.
- dhdt.preprocessing.atmospheric_geometry.get_water_vapor_press_ice(T)¶
- Parameters
T ({float, numpy.array}, unit=Kelvin) – temperature
- Returns
P_w – water vapor pressure
- Return type
{float, numpy.array}, unit=Pascal
See also
References
- MM93
Marti & Mauersberger, “A survey and new measurements of ice vapor pressure at temperatures between 170 and 250K” Geophysical research letters, vol.20(5) pp.363-366, 1993.
- dhdt.preprocessing.atmospheric_geometry.refraction_angle_analytical(zn_0, n_0)¶
- Parameters
zn_0 ({float, numpy.array}, unit=degrees) – observation angle, when atmosphere is not taken into account
n_0 ({float, numpy.array}) – refractive index at the surface
- Returns
zn – corrected observation angle, taking refraction into account
- Return type
{float, numpy.array}, unit=degrees
References
- No99
Noerdlinger, “Atmospheric refraction effects in Earth remote sensing” ISPRS journal of photogrammetry and remote sensing, vol.54, pp.360-373, 1999.
- dhdt.preprocessing.atmospheric_geometry.refraction_spherical_symmetric(zn_0, lat=None, h_0=None)¶
- Parameters
zn_0 – observation angle, when atmosphere is not taken into account
float – observation angle, when atmosphere is not taken into account
lat (float, unit=degrees, range=-90...+90) – latitude of point of interest
h_0 (float, unit=meters) – initial altitude above geoid
- Returns
zn_hat – analytical refraction angle
- Return type
float
See also
get_mean_lat
Notes
The angles related to the sun are as follows:
surface normal * sun ↑ ↑ / | | / |-- zenith angle | / | / | /| |/ |/ | elevation angle +---- +------
References
- No99
Noerdlinger, “Atmospheric refraction effects in Earth remote sensing” ISPRS journal of photogrammetry and remote sensing, vol.54, pp.360-373, 1999.
- BJ94
Birch and Jones, “Correction to the updated Edlen equation for the refractive index of air”, Metrologica, vol.31(4) pp.315-316, 1994.
- dhdt.preprocessing.atmospheric_geometry.refractive_index_CO2(n_s, CO2)¶
correct standard refractivity for diffferent carbon dioxide concentrations, see [2] for more information.
- Parameters
n_s ({float, numpy.array}) – refractivity of standard air
CO2 ({float, numpy.array}, unit=ppm) – concentration of carbon dioxide (CO₂) in the atmosphere
- Returns
n_x – refractivity of standard air at a given CO₂ concentration
- Return type
{float, numpy.array}
Notes
The following acronyms are used:
CO2 : carbon dioxide
ppm : parts per million
References
- BD93
Birch & Downs, “An updated Edlen equation for the refractive index of air” Metrologica, vol.30(155) pp.155-162, 1993.
- Ed66
Edlén, “The refractive index of air”, Metrologica, vol.2(2) pp.71-80, 1966.
- dhdt.preprocessing.atmospheric_geometry.refractive_index_broadband(df, T_0, P_0, fH_0, CO2=450.0, T_a=15.0, P_a=101325.0, T_w=20.0, P_w=1333.0, LorentzLorenz=False)¶
Estimate the refractive index in for wavelengths in the range of 300…1700 nm, see also [1].
- Parameters
df ({dataframe, float}, unit=µm) – central wavelengths of the spectral bands
CO2 (float, unit=ppm of CO₂, default=450.) – parts per million of CO₂
T_0 (float, unit=Celsius) – temperature at location
P_0 (float, unit=Pascal) – the total pressure in Pascal
fH_0 (float, range=0...1) – fractional humidity
T_a (float, unit=Celcius) – temperature of standard dry air, see Notes
P_a (float, unit=Pascal) – pressure of standard dry air, see Notes
T_w (float, unit=Celsius) – temperature of standard moist air, see [1]
P_w (float, unit=Pascal) – pressure of standard moist air, see [1]
LorentzLorenz (bool) –
Two ways are possible to calculate the proportional refraction index, when
True : the Lorentz-Lorenz relationship is used.
False : component wise addition is used
- Returns
n – index of refraction
- Return type
{float, numpy.array}
Notes
Standard air [2] is defined to be dry air at a temperature of 15 degrees Celsius and a total pressure of 1013.25 mb, having the following composition by molar percentage:
78.09% nitrogen (N₂)
20.95% oxygen (O₂)
0.93% argon (Ar)
0.03% carbon dioxide (C0₂)
See also
ciddor_eq5
calculation used for the Lorentz-Lorenz relationship
ciddor_eq6
calculation used for component wise addition
References
- Ci96
Ciddor, “Refractive index of air: new equations for the visible and near infrared”, Applied optics, vol.35(9) pp.1566-1573, 1996.
- Ow67
Owens, “Optical refractive index of air: Dependence on pressure, temperature and composition”, Applied optics, vol.6(1) pp.51-59, 1967.
- dhdt.preprocessing.atmospheric_geometry.refractive_index_broadband_vapour(σ)¶
- Parameters
σ ({float, numpy.array}, unit=µm-1) – vacuum wavenumber, that is the reciprocal of the vacuum wavelength
- Returns
n_ws – refractive index for standard moist air at sea level
- Return type
{float, numpy.array}
References
- Ow67
Owens, “Optical refractive index of air: Dependence on pressure, temperature and composition”, Applied optics, vol.6(1) pp.51-59, 1967.
- Ci96
Ciddor, “Refractive index of air: new equations for the visible and near infrared”, Applied optics, vol.35(9) pp.1566-1573, 1996.
- dhdt.preprocessing.atmospheric_geometry.refractive_index_simple(p, RH, t)¶
simple shop-floor formula fro refraction, see appendix B in [1].
- Parameters
p ({numpy.array, float}, unit=kPa) – pressure
RH ({numpy.array, float}, unit=percent) – relative humidity
t ({numpy.array, float}, unit=Celsius) – temperature
- Returns
n – refractive index
- Return type
{numpy.array, float}
References
- dhdt.preprocessing.atmospheric_geometry.refractive_index_visible(df, T, P, p_w=None, CO2=None)¶
calculate refractive index, based on measurement around 633nm, see [1].
It applies to ambient atmospheric conditions over the range of wavelengths 350 nm to 650 nm.
- Parameters
df ({dataframe, float}, unit=µm) – central wavelengths of the spectral bands
T ({float, numpy.array}, unit=Celcius) – temperature at location
P ({float, numpy.array}, unit=Pascal) – atmospheric pressure at location
CO2 ({float, numpy.array}, unit=ppm) – concentration of carbon dioxide (CO₂) in the atmosphere
- Returns
n_0 – index of refraction
- Return type
numpy.array
See also
References
- BJ94
Birch and Jones, “Correction to the updated Edlén equation for the refractive index of air”, Metrologica, vol.31(4) pp.315-316, 1994.
- Ed66
Edlén, “The refractive index of air”, Metrologica, vol.2(2) pp.71-80, 1966.
- BJ94
Birch and Jones, “An updated Edlén equation for the refractive index of air”, Metrologica, vol.30 pp.155-162, 1993.
- dhdt.preprocessing.atmospheric_geometry.update_list_with_corr_zenith_pd(dh, file_dir, file_name='conn.txt')¶
- Parameters
dh (pandas.DataFrame) – array with photohypsometric information
file_dir (string) – directory where the file should be situated
file_name (string, default='conn.txt') – name of the file to export the data columns to
Notes
The nomenclature is as follows:
* sun \ <-- sun trace |\ <-- caster | | \ ^ Z | \ | ----┴----+ <-- casted └-> {X,Y}
The angles related to the sun are as follows:
surface normal * sun ^ ^ / | | / |-- zenith angle | / | / | /| |/ |/ | elevation angle └---- └------ {X,Y}
The azimuth angle declared in the following coordinate frame:
^ North & Y | - <--|--> + | +----> East & X
The text file has at least the following columns:
* conn.txt ├ caster_X <- map location of the start of the shadow trace ├ caster_Y ├ casted_X ├ casted_Y <- map location of the end of the shadow trace ├ azimuth <- argument of the illumination direction ├ zenith <- off-normal angle of the sun without atmosphere └ zenith_refrac <- off-normal angle of the sun with atmosphere
See also
get_refraction_angle
- dhdt.preprocessing.atmospheric_geometry.water_vapor_frac(p, T, fH)¶
- Parameters
p (float, unit=Pascal) – total pressure
T (float, unit=Kelvin) – temperature
fH (float, range=0...1) – fractional humidity
- Returns
x_w – molar fraction of water vapor in moist air
- Return type
float
Notes
The following parameters are used in this function:
T : unit=Kelvin, temperature svp : unit=Pascal, saturation vapor pressure t : unit=Celcius, temperature
References
- Ci96
Ciddor, “Refractive index of air: new equations for the visible and near infrared”, Applied optics, vol.35(9) pp.1566-1573, 1996.
- MM93
Marti & Mauersberger, “A survey and new measurements of ice vapor pressure at temperatures between 170 and 250K” Geophysical research letters, vol.20(5) pp.363-366, 1993.
- Gi82
Giacomo, “Equation for the determination of the density of moist air” Metrologica, vol.18, pp.33-40, 1982.