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