1
|
|
|
# Copyright (c) 2008-2016 MetPy Developers. |
2
|
|
|
# Distributed under the terms of the BSD 3-Clause License. |
3
|
|
|
# SPDX-License-Identifier: BSD-3-Clause |
4
|
|
|
""" |
5
|
|
|
Inverse Distance Verification: Cressman and Barnes |
6
|
|
|
================================================== |
7
|
|
|
|
8
|
|
|
Compare inverse distance interpolation methods |
9
|
|
|
|
10
|
|
|
Two popular interpolation schemes that use inverse distance weighting of observations are the |
11
|
|
|
Barnes and Cressman analyses. The Cressman analysis is relatively straightforward and uses |
12
|
|
|
the ratio between distance of an observation from a grid cell and the maximum allowable |
13
|
|
|
distance to calculate the relative importance of an observation for calculating an |
14
|
|
|
interpolation value. Barnes uses the inverse exponential ratio of each distance between |
15
|
|
|
an observation and a grid cell and the average spacing of the observations over the domain. |
16
|
|
|
|
17
|
|
|
Algorithmically: |
18
|
|
|
|
19
|
|
|
1. A KDTree data structure is built using the locations of each observation. |
20
|
|
|
2. All observations within a maximum allowable distance of a particular grid cell are found in |
21
|
|
|
O(log n) time. |
22
|
|
|
3. Using the weighting rules for Cressman or Barnes analyses, the observations are given a |
23
|
|
|
proportional value, primarily based on their distance from the grid cell. |
24
|
|
|
4. The sum of these proportional values is calculated and this value is used as the |
25
|
|
|
interpolated value. |
26
|
|
|
5. Steps 2 through 4 are repeated for each grid cell. |
27
|
|
|
""" |
28
|
|
|
import matplotlib.pyplot as plt |
29
|
|
|
import numpy as np |
30
|
|
|
from scipy.spatial import cKDTree |
31
|
|
|
from scipy.spatial.distance import cdist |
32
|
|
|
|
33
|
|
|
from metpy.gridding.gridding_functions import calc_kappa |
34
|
|
|
from metpy.gridding.interpolation import barnes_point, cressman_point |
35
|
|
|
from metpy.gridding.triangles import dist_2 |
36
|
|
|
|
37
|
|
|
plt.rcParams['figure.figsize'] = (15, 10) |
38
|
|
|
|
39
|
|
|
|
40
|
|
|
def draw_circle(x, y, r, m, label): |
41
|
|
|
nx = x + r * np.cos(np.deg2rad(list(range(360)))) |
42
|
|
|
ny = y + r * np.sin(np.deg2rad(list(range(360)))) |
43
|
|
|
plt.plot(nx, ny, m, label=label) |
44
|
|
|
|
45
|
|
|
|
46
|
|
|
########################################### |
47
|
|
|
# Generate random x and y coordinates, and observation values proportional to x * y. |
48
|
|
|
# |
49
|
|
|
# Set up two test grid locations at (30, 30) and (60, 60). |
50
|
|
|
np.random.seed(100) |
51
|
|
|
|
52
|
|
|
pts = np.random.randint(0, 100, (10, 2)) |
53
|
|
|
xp = pts[:, 0] |
54
|
|
|
yp = pts[:, 1] |
55
|
|
|
zp = xp * xp / 1000 |
56
|
|
|
|
57
|
|
|
sim_gridx = [30, 60] |
58
|
|
|
sim_gridy = [30, 60] |
59
|
|
|
|
60
|
|
|
########################################### |
61
|
|
|
# Set up a cKDTree object and query all of the observations within "radius" of each grid point. |
62
|
|
|
# |
63
|
|
|
# The variable ``indices`` represents the index of each matched coordinate within the |
64
|
|
|
# cKDTree's ``data`` list. |
65
|
|
|
grid_points = np.array(list(zip(sim_gridx, sim_gridy))) |
66
|
|
|
|
67
|
|
|
radius = 40 |
68
|
|
|
obs_tree = cKDTree(list(zip(xp, yp))) |
69
|
|
|
indices = obs_tree.query_ball_point(grid_points, r=radius) |
70
|
|
|
|
71
|
|
|
########################################### |
72
|
|
|
# For grid 0, we will use Cressman to interpolate its value. |
73
|
|
|
x1, y1 = obs_tree.data[indices[0]].T |
74
|
|
|
cress_dist = dist_2(sim_gridx[0], sim_gridy[0], x1, y1) |
75
|
|
|
cress_obs = zp[indices[0]] |
76
|
|
|
|
77
|
|
|
cress_val = cressman_point(cress_dist, cress_obs, radius) |
78
|
|
|
|
79
|
|
|
########################################### |
80
|
|
|
# For grid 1, we will use barnes to interpolate its value. |
81
|
|
|
# |
82
|
|
|
# We need to calculate kappa--the average distance between observations over the domain. |
83
|
|
|
x2, y2 = obs_tree.data[indices[1]].T |
84
|
|
|
barnes_dist = dist_2(sim_gridx[1], sim_gridy[1], x2, y2) |
85
|
|
|
barnes_obs = zp[indices[1]] |
86
|
|
|
|
87
|
|
|
ave_spacing = np.mean((cdist(list(zip(xp, yp)), list(zip(xp, yp))))) |
88
|
|
|
kappa = calc_kappa(ave_spacing) |
89
|
|
|
|
90
|
|
|
barnes_val = barnes_point(barnes_dist, barnes_obs, kappa) |
91
|
|
|
|
92
|
|
|
########################################### |
93
|
|
|
# Plot all of the affiliated information and interpolation values. |
94
|
|
|
for i, zval in enumerate(zp): |
95
|
|
|
plt.plot(pts[i, 0], pts[i, 1], '.') |
96
|
|
|
plt.annotate(str(zval) + ' F', xy=(pts[i, 0] + 2, pts[i, 1])) |
97
|
|
|
|
98
|
|
|
plt.plot(sim_gridx, sim_gridy, '+', markersize=10) |
99
|
|
|
|
100
|
|
|
plt.plot(x1, y1, 'ko', fillstyle='none', markersize=10, label='grid 0 matches') |
101
|
|
|
plt.plot(x2, y2, 'ks', fillstyle='none', markersize=10, label='grid 1 matches') |
102
|
|
|
|
103
|
|
|
draw_circle(sim_gridx[0], sim_gridy[0], m='k-', r=radius, label='grid 0 radius') |
104
|
|
|
draw_circle(sim_gridx[1], sim_gridy[1], m='b-', r=radius, label='grid 1 radius') |
105
|
|
|
|
106
|
|
|
plt.annotate('grid 0: cressman {:.3f}'.format(cress_val), xy=(sim_gridx[0] + 2, sim_gridy[0])) |
107
|
|
|
plt.annotate('grid 1: barnes {:.3f}'.format(barnes_val), xy=(sim_gridx[1] + 2, sim_gridy[1])) |
108
|
|
|
|
109
|
|
|
plt.axes().set_aspect('equal', 'datalim') |
110
|
|
|
plt.legend() |
111
|
|
|
|
112
|
|
|
########################################### |
113
|
|
|
# For each point, we will do a manual check of the interpolation values by doing a step by |
114
|
|
|
# step and visual breakdown. |
115
|
|
|
# |
116
|
|
|
# Plot the grid point, observations within radius of the grid point, their locations, and |
117
|
|
|
# their distances from the grid point. |
118
|
|
|
plt.annotate('grid 0: ({}, {})'.format(sim_gridx[0], sim_gridy[0]), |
119
|
|
|
xy=(sim_gridx[0] + 2, sim_gridy[0])) |
120
|
|
|
plt.plot(sim_gridx[0], sim_gridy[0], '+', markersize=10) |
121
|
|
|
|
122
|
|
|
mx, my = obs_tree.data[indices[0]].T |
123
|
|
|
mz = zp[indices[0]] |
124
|
|
|
|
125
|
|
|
for x, y, z in zip(mx, my, mz): |
126
|
|
|
d = np.sqrt((sim_gridx[0] - x)**2 + (y - sim_gridy[0])**2) |
127
|
|
|
plt.plot([sim_gridx[0], x], [sim_gridy[0], y], '--') |
128
|
|
|
|
129
|
|
|
xave = np.mean([sim_gridx[0], x]) |
130
|
|
|
yave = np.mean([sim_gridy[0], y]) |
131
|
|
|
|
132
|
|
|
plt.annotate('distance: {}'.format(d), xy=(xave, yave)) |
133
|
|
|
plt.annotate('({}, {}) : {} F'.format(x, y, z), xy=(x, y)) |
134
|
|
|
|
135
|
|
|
plt.xlim(0, 80) |
136
|
|
|
plt.ylim(0, 80) |
137
|
|
|
plt.axes().set_aspect('equal', 'datalim') |
138
|
|
|
|
139
|
|
|
########################################### |
140
|
|
|
# Step through the cressman calculations. |
141
|
|
|
dists = np.array([22.803508502, 7.21110255093, 31.304951685, 33.5410196625]) |
142
|
|
|
values = np.array([0.064, 1.156, 3.364, 0.225]) |
143
|
|
|
|
144
|
|
|
cres_weights = (radius * radius - dists * dists) / (radius * radius + dists * dists) |
145
|
|
|
total_weights = np.sum(cres_weights) |
146
|
|
|
proportion = cres_weights / total_weights |
147
|
|
|
value = values * proportion |
148
|
|
|
|
149
|
|
|
val = cressman_point(cress_dist, cress_obs, radius) |
150
|
|
|
|
151
|
|
|
print('Manual cressman value for grid 1:\t', np.sum(value)) |
152
|
|
|
print('Metpy cressman value for grid 1:\t', val) |
153
|
|
|
|
154
|
|
|
########################################### |
155
|
|
|
# Now repeat for grid 1, except use barnes interpolation. |
156
|
|
|
|
157
|
|
|
plt.annotate('grid 1: ({}, {})'.format(sim_gridx[1], sim_gridy[1]), |
158
|
|
|
xy=(sim_gridx[1] + 2, sim_gridy[1])) |
159
|
|
|
plt.plot(sim_gridx[1], sim_gridy[1], '+', markersize=10) |
160
|
|
|
|
161
|
|
|
mx, my = obs_tree.data[indices[1]].T |
162
|
|
|
mz = zp[indices[1]] |
163
|
|
|
|
164
|
|
|
for x, y, z in zip(mx, my, mz): |
165
|
|
|
d = np.sqrt((sim_gridx[1] - x)**2 + (y - sim_gridy[1])**2) |
166
|
|
|
plt.plot([sim_gridx[1], x], [sim_gridy[1], y], '--') |
167
|
|
|
|
168
|
|
|
xave = np.mean([sim_gridx[1], x]) |
169
|
|
|
yave = np.mean([sim_gridy[1], y]) |
170
|
|
|
|
171
|
|
|
plt.annotate('distance: {}'.format(d), xy=(xave, yave)) |
172
|
|
|
plt.annotate('({}, {}) : {} F'.format(x, y, z), xy=(x, y)) |
173
|
|
|
|
174
|
|
|
plt.xlim(40, 80) |
175
|
|
|
plt.ylim(40, 100) |
176
|
|
|
plt.axes().set_aspect('equal', 'datalim') |
177
|
|
|
|
178
|
|
|
########################################### |
179
|
|
|
# Step through barnes calculations. |
180
|
|
|
dists = np.array([9.21954445729, 22.4722050542, 27.892651362, 38.8329756779]) |
181
|
|
|
values = np.array([2.809, 6.241, 4.489, 2.704]) |
182
|
|
|
|
183
|
|
|
weights = np.exp(-dists**2 / kappa) |
184
|
|
|
total_weights = np.sum(weights) |
185
|
|
|
value = np.sum(values * (weights / total_weights)) |
186
|
|
|
|
187
|
|
|
print('Manual barnes value:\t', value) |
188
|
|
|
print('Metpy barnes value:\t', barnes_point(barnes_dist, barnes_obs, kappa)) |
189
|
|
|
|