1
|
|
|
# Copyright (c) 2008-2015 MetPy Developers. |
2
|
|
|
# Distributed under the terms of the BSD 3-Clause License. |
3
|
|
|
# SPDX-License-Identifier: BSD-3-Clause |
4
|
|
|
"""Contains calculation of kinematic parameters (e.g. divergence or vorticity).""" |
5
|
|
|
from __future__ import division |
6
|
|
|
|
7
|
|
|
import functools |
8
|
|
|
import warnings |
9
|
|
|
|
10
|
|
|
import numpy as np |
11
|
|
|
|
12
|
|
|
from ..cbook import is_string_like, iterable |
13
|
|
|
from ..constants import g |
14
|
|
|
from ..package_tools import Exporter |
15
|
|
|
from ..units import atleast_2d, concatenate, units |
16
|
|
|
|
17
|
|
|
exporter = Exporter(globals()) |
18
|
|
|
|
19
|
|
|
|
20
|
|
|
def _gradient(f, *args, **kwargs): |
21
|
|
|
"""Wrapper for :func:`numpy.gradient` that handles units.""" |
22
|
|
|
if len(args) < f.ndim: |
23
|
|
|
args = list(args) |
24
|
|
|
args.extend([units.Quantity(1., 'dimensionless')] * (f.ndim - len(args))) |
25
|
|
|
grad = np.gradient(f, *(a.magnitude for a in args), **kwargs) |
26
|
|
|
if f.ndim == 1: |
27
|
|
|
return units.Quantity(grad, f.units / args[0].units) |
28
|
|
|
return [units.Quantity(g, f.units / dx.units) for dx, g in zip(args, grad)] |
29
|
|
|
|
30
|
|
|
|
31
|
|
|
def _stack(arrs): |
32
|
|
|
return concatenate([a[np.newaxis] for a in arrs], axis=0) |
33
|
|
|
|
34
|
|
|
|
35
|
|
|
def _get_gradients(u, v, dx, dy): |
36
|
|
|
"""Helper function for getting convergence and vorticity from 2D arrays.""" |
37
|
|
|
dudy, dudx = _gradient(u, dy, dx) |
38
|
|
|
dvdy, dvdx = _gradient(v, dy, dx) |
39
|
|
|
return dudx, dudy, dvdx, dvdy |
40
|
|
|
|
41
|
|
|
|
42
|
|
|
def _is_x_first_dim(dim_order): |
43
|
|
|
"""Determine whether x is the first dimension based on the value of dim_order.""" |
44
|
|
|
if dim_order is None: |
45
|
|
|
warnings.warn('dim_order is using the default setting (currently "xy"). This will ' |
46
|
|
|
'change to "yx" in the next version. It is recommended that you ' |
47
|
|
|
'specify the appropriate ordering ("xy", "yx") for your data by ' |
48
|
|
|
'passing the `dim_order` argument to the calculation.', FutureWarning) |
49
|
|
|
dim_order = 'xy' |
50
|
|
|
return dim_order == 'xy' |
51
|
|
|
|
52
|
|
|
|
53
|
|
|
def _check_and_flip(arr): |
54
|
|
|
"""Transpose array or list of arrays if they are 2D.""" |
55
|
|
|
if hasattr(arr, 'ndim'): |
56
|
|
|
if arr.ndim >= 2: |
57
|
|
|
return arr.T |
58
|
|
|
else: |
59
|
|
|
return arr |
60
|
|
|
elif not is_string_like(arr) and iterable(arr): |
61
|
|
|
return tuple(_check_and_flip(a) for a in arr) |
62
|
|
|
else: |
63
|
|
|
return arr |
64
|
|
|
|
65
|
|
|
|
66
|
|
|
def ensure_yx_order(func): |
67
|
|
|
"""Wrap a function to ensure all array arguments are y, x ordered, based on kwarg.""" |
68
|
|
|
@functools.wraps(func) |
69
|
|
|
def wrapper(*args, **kwargs): |
70
|
|
|
# Check what order we're given |
71
|
|
|
dim_order = kwargs.pop('dim_order', None) |
72
|
|
|
x_first = _is_x_first_dim(dim_order) |
73
|
|
|
|
74
|
|
|
# If x is the first dimension, flip (transpose) every array within the function args. |
75
|
|
|
if x_first: |
76
|
|
|
args = tuple(_check_and_flip(arr) for arr in args) |
77
|
|
|
for k, v in kwargs: |
78
|
|
|
kwargs[k] = _check_and_flip(v) |
79
|
|
|
|
80
|
|
|
ret = func(*args, **kwargs) |
81
|
|
|
|
82
|
|
|
# If we flipped on the way in, need to flip on the way out so that output array(s) |
83
|
|
|
# match the dimension order of the original input. |
84
|
|
|
if x_first: |
85
|
|
|
return _check_and_flip(ret) |
86
|
|
|
else: |
87
|
|
|
return ret |
88
|
|
|
|
89
|
|
|
# Inject a docstring for the dim_order argument into the function's docstring. |
90
|
|
|
dim_order_doc = ''' |
91
|
|
|
dim_order : str or ``None``, optional |
92
|
|
|
The ordering of dimensions in passed in arrays. Can be one of ``None``, ``'xy'``, |
93
|
|
|
or ``'yx'``. ``'xy'`` indicates that the dimension corresponding to x is the leading |
94
|
|
|
dimension, followed by y. ``'yx'`` indicates that x is the last dimension, preceded |
95
|
|
|
by y. ``None`` indicates that the default ordering should be assumed, |
96
|
|
|
which will change in version 0.6 from 'xy' to 'yx'. Can only be passed as a keyword |
97
|
|
|
argument, i.e. func(..., dim_order='xy').''' |
98
|
|
|
|
99
|
|
|
# Find the first blank line after the start of the parameters section |
100
|
|
|
params = wrapper.__doc__.find('Parameters') |
101
|
|
|
blank = wrapper.__doc__.find('\n\n', params) |
102
|
|
|
wrapper.__doc__ = wrapper.__doc__[:blank] + dim_order_doc + wrapper.__doc__[blank:] |
103
|
|
|
|
104
|
|
|
return wrapper |
105
|
|
|
|
106
|
|
|
|
107
|
|
|
@exporter.export |
108
|
|
|
@ensure_yx_order |
109
|
|
|
def v_vorticity(u, v, dx, dy): |
110
|
|
|
r"""Calculate the vertical vorticity of the horizontal wind. |
111
|
|
|
|
112
|
|
|
The grid must have a constant spacing in each direction. |
113
|
|
|
|
114
|
|
|
Parameters |
115
|
|
|
---------- |
116
|
|
|
u : (M, N) ndarray |
117
|
|
|
x component of the wind |
118
|
|
|
v : (M, N) ndarray |
119
|
|
|
y component of the wind |
120
|
|
|
dx : float |
121
|
|
|
The grid spacing in the x-direction |
122
|
|
|
dy : float |
123
|
|
|
The grid spacing in the y-direction |
124
|
|
|
|
125
|
|
|
Returns |
126
|
|
|
------- |
127
|
|
|
(M, N) ndarray |
128
|
|
|
vertical vorticity |
129
|
|
|
|
130
|
|
|
See Also |
131
|
|
|
-------- |
132
|
|
|
h_convergence, convergence_vorticity |
133
|
|
|
""" |
134
|
|
|
_, dudy, dvdx, _ = _get_gradients(u, v, dx, dy) |
135
|
|
|
return dvdx - dudy |
136
|
|
|
|
137
|
|
|
|
138
|
|
|
@exporter.export |
139
|
|
|
@ensure_yx_order |
140
|
|
|
def h_convergence(u, v, dx, dy): |
141
|
|
|
r"""Calculate the horizontal convergence of the horizontal wind. |
142
|
|
|
|
143
|
|
|
The grid must have a constant spacing in each direction. |
144
|
|
|
|
145
|
|
|
Parameters |
146
|
|
|
---------- |
147
|
|
|
u : (M, N) ndarray |
148
|
|
|
x component of the wind |
149
|
|
|
v : (M, N) ndarray |
150
|
|
|
y component of the wind |
151
|
|
|
dx : float |
152
|
|
|
The grid spacing in the x-direction |
153
|
|
|
dy : float |
154
|
|
|
The grid spacing in the y-direction |
155
|
|
|
|
156
|
|
|
Returns |
157
|
|
|
------- |
158
|
|
|
(M, N) ndarray |
159
|
|
|
The horizontal convergence |
160
|
|
|
|
161
|
|
|
See Also |
162
|
|
|
-------- |
163
|
|
|
v_vorticity, convergence_vorticity |
164
|
|
|
""" |
165
|
|
|
dudx, _, _, dvdy = _get_gradients(u, v, dx, dy) |
166
|
|
|
return dudx + dvdy |
167
|
|
|
|
168
|
|
|
|
169
|
|
|
@exporter.export |
170
|
|
|
@ensure_yx_order |
171
|
|
|
def convergence_vorticity(u, v, dx, dy): |
172
|
|
|
r"""Calculate the horizontal convergence and vertical vorticity of the horizontal wind. |
173
|
|
|
|
174
|
|
|
The grid must have a constant spacing in each direction. |
175
|
|
|
|
176
|
|
|
Parameters |
177
|
|
|
---------- |
178
|
|
|
u : (M, N) ndarray |
179
|
|
|
x component of the wind |
180
|
|
|
v : (M, N) ndarray |
181
|
|
|
y component of the wind |
182
|
|
|
dx : float |
183
|
|
|
The grid spacing in the x-direction |
184
|
|
|
dy : float |
185
|
|
|
The grid spacing in the y-direction |
186
|
|
|
|
187
|
|
|
Returns |
188
|
|
|
------- |
189
|
|
|
convergence, vorticity : tuple of (M, N) ndarrays |
190
|
|
|
The horizontal convergence and vertical vorticity, respectively |
191
|
|
|
|
192
|
|
|
See Also |
193
|
|
|
-------- |
194
|
|
|
v_vorticity, h_convergence |
195
|
|
|
|
196
|
|
|
Notes |
197
|
|
|
----- |
198
|
|
|
This is a convenience function that will do less work than calculating |
199
|
|
|
the horizontal convergence and vertical vorticity separately. |
200
|
|
|
""" |
201
|
|
|
dudx, dudy, dvdx, dvdy = _get_gradients(u, v, dx, dy) |
202
|
|
|
return dudx + dvdy, dvdx - dudy |
203
|
|
|
|
204
|
|
|
|
205
|
|
|
@exporter.export |
206
|
|
|
@ensure_yx_order |
207
|
|
|
def advection(scalar, wind, deltas): |
208
|
|
|
r"""Calculate the advection of a scalar field by the wind. |
209
|
|
|
|
210
|
|
|
The order of the dimensions of the arrays must match the order in which |
211
|
|
|
the wind components are given. For example, if the winds are given [u, v], |
212
|
|
|
then the scalar and wind arrays must be indexed as x,y (which puts x as the |
213
|
|
|
rows, not columns). |
214
|
|
|
|
215
|
|
|
Parameters |
216
|
|
|
---------- |
217
|
|
|
scalar : N-dimensional array |
218
|
|
|
Array (with N-dimensions) with the quantity to be advected. |
219
|
|
|
wind : sequence of arrays |
220
|
|
|
Length N sequence of N-dimensional arrays. Represents the flow, |
221
|
|
|
with a component of the wind in each dimension. For example, for |
222
|
|
|
horizontal advection, this could be a list: [u, v], where u and v |
223
|
|
|
are each a 2-dimensional array. |
224
|
|
|
deltas : sequence |
225
|
|
|
A (length N) sequence containing the grid spacing in each dimension. |
226
|
|
|
|
227
|
|
|
Returns |
228
|
|
|
------- |
229
|
|
|
N-dimensional array |
230
|
|
|
An N-dimensional array containing the advection at all grid points. |
231
|
|
|
""" |
232
|
|
|
# This allows passing in a list of wind components or an array. |
233
|
|
|
wind = _stack(wind) |
234
|
|
|
|
235
|
|
|
# If we have more than one component, we need to reverse the order along the first |
236
|
|
|
# dimension so that the wind components line up with the |
237
|
|
|
# order of the gradients from the ..., y, x ordered array. |
238
|
|
|
if wind.ndim > scalar.ndim: |
239
|
|
|
wind = wind[::-1] |
240
|
|
|
|
241
|
|
|
# Gradient returns a list of derivatives along each dimension. We convert |
242
|
|
|
# this to an array with dimension as the first index. Reverse the deltas to line up |
243
|
|
|
# with the order of the dimensions. |
244
|
|
|
grad = _stack(_gradient(scalar, *deltas[::-1])) |
245
|
|
|
|
246
|
|
|
# Make them be at least 2D (handling the 1D case) so that we can do the |
247
|
|
|
# multiply and sum below |
248
|
|
|
grad, wind = atleast_2d(grad, wind) |
249
|
|
|
|
250
|
|
|
return (-grad * wind).sum(axis=0) |
251
|
|
|
|
252
|
|
|
|
253
|
|
|
@exporter.export |
254
|
|
|
@ensure_yx_order |
255
|
|
|
def geostrophic_wind(heights, f, dx, dy): |
256
|
|
|
r"""Calculate the geostrophic wind given from the heights or geopotential. |
257
|
|
|
|
258
|
|
|
Parameters |
259
|
|
|
---------- |
260
|
|
|
heights : (M, N) ndarray |
261
|
|
|
The height field, with either leading dimensions of (x, y) or trailing dimensions |
262
|
|
|
of (y, x), depending on the value of ``dim_order``. |
263
|
|
|
f : array_like |
264
|
|
|
The coriolis parameter. This can be a scalar to be applied |
265
|
|
|
everywhere or an array of values. |
266
|
|
|
dx : scalar |
267
|
|
|
The grid spacing in the x-direction |
268
|
|
|
dy : scalar |
269
|
|
|
The grid spacing in the y-direction |
270
|
|
|
|
271
|
|
|
Returns |
272
|
|
|
------- |
273
|
|
|
A 2-item tuple of arrays |
274
|
|
|
A tuple of the u-component and v-component of the geostrophic wind. |
275
|
|
|
""" |
276
|
|
|
if heights.dimensionality['[length]'] == 2.0: |
277
|
|
|
norm_factor = 1. / f |
278
|
|
|
else: |
279
|
|
|
norm_factor = g / f |
280
|
|
|
|
281
|
|
|
# If heights has more than 2 dimensions, we need to pass in some dummy |
282
|
|
|
# grid deltas so that we can still use np.gradient. It may be better to |
283
|
|
|
# to loop in this case, but that remains to be done. |
284
|
|
|
deltas = [dy, dx] |
285
|
|
|
if heights.ndim > 2: |
286
|
|
|
deltas = [units.Quantity(1., units.m)] * (heights.ndim - 2) + deltas |
287
|
|
|
|
288
|
|
|
grad = _gradient(heights, *deltas) |
289
|
|
|
dy, dx = grad[-2:] # Throw away unused gradient components |
290
|
|
|
return -norm_factor * dy, norm_factor * dx |
291
|
|
|
|