@@ -500,6 +500,40 @@ def geostrophic_wind(heights, f, dx, dy):
500
500
return - norm_factor * dhdy , norm_factor * dhdx
501
501
502
502
503
+ @exporter .export
504
+ @ensure_yx_order
505
+ def ageostrophic_wind (heights , f , dx , dy , u , v , dim_order = 'yx' ):
506
+ r"""Calculate the ageostrophic wind given from the heights or geopotential.
507
+
508
+ Parameters
509
+ ----------
510
+ heights : (M, N) ndarray
511
+ The height field, with either leading dimensions of (x, y) or trailing dimensions
512
+ of (y, x), depending on the value of ``dim_order``.
513
+ f : array_like
514
+ The coriolis parameter. This can be a scalar to be applied
515
+ everywhere or an array of values.
516
+ dx : scalar
517
+ The grid spacing in the x-direction
518
+ dy : scalar
519
+ The grid spacing in the y-direction
520
+ u : (M, N) ndarray
521
+ The u wind field, with either leading dimensions of (x, y) or trailing dimensions
522
+ of (y, x), depending on the value of ``dim_order``.
523
+ v : (M, N) ndarray
524
+ The u wind field, with either leading dimensions of (x, y) or trailing dimensions
525
+ of (y, x), depending on the value of ``dim_order``.
526
+
527
+ Returns
528
+ -------
529
+ A 2-item tuple of arrays
530
+ A tuple of the u-component and v-component of the ageostrophic wind.
531
+
532
+ """
533
+ u_geostrophic , v_geostrophic = geostrophic_wind (heights , f , dx , dy , dim_order = dim_order )
534
+ return u - u_geostrophic , v - v_geostrophic
535
+
536
+
503
537
@exporter .export
504
538
@check_units ('[length]' , '[temperature]' )
505
539
def montgomery_streamfunction (height , temperature ):
0 commit comments