FV3 Bundle
gsw_turner_rsubrho.f90
Go to the documentation of this file.
1 !==========================================================================
2 pure subroutine gsw_turner_rsubrho (sa, ct, p, tu, rsubrho, p_mid)
3 !==========================================================================
4 !
5 ! Calculates the Turner angle and the Rsubrho as a function of pressure
6 ! down a vertical water column. These quantities express the relative
7 ! contributions of the vertical gradients of Conservative Temperature
8 ! and Absolute Salinity to the vertical stability (the square of the
9 ! Brunt-Vaisala Frequency squared, N^2). Tu and Rsubrho are evaluated at
10 ! the mid pressure between the individual data points in the vertical.
11 ! density in terms of SA, CT and p (IOC et al., 2010).
12 !
13 ! Note that in the double-diffusive literature, papers concerned with
14 ! the "diffusive" form of double-diffusive convection often define the
15 ! stability ratio as the reciprocal of what is defined here as the
16 ! stability ratio.
17 !
18 ! sa : Absolute Salinity (a profile (length nz)) [g/kg]
19 ! ct : Conservative Temperature (a profile (length nz)) [deg C]
20 ! p : sea pressure (a profile (length nz)) [dbar]
21 ! tu : Turner angle, on the same (nz-1) grid as p_mid.
22 ! Turner angle has units of: [ degrees of rotation ]
23 ! rsubrho : Stability Ratio, on the same (nz-1) grid as p_mid.
24 ! rsubrho is dimensionless. [ unitless ]
25 ! p_mid : Mid pressure between p grid (length nz-1) [dbar]
26 !--------------------------------------------------------------------------
27 
29 
31 
33 
34 use gsw_mod_kinds
35 
36 implicit none
37 
38 real (r8), intent(in) :: sa(:), ct(:), p(:)
39 real (r8), intent(out) :: tu(:), rsubrho(:), p_mid(:)
40 
41 integer :: nz, k
42 real (r8), dimension(:), allocatable :: dsa, sa_mid, dct, ct_mid, dp
43 real (r8), dimension(:), allocatable :: alpha_mid, beta_mid
44 
45 character (*), parameter :: func_name = "gsw_turner_rsubrho"
46 
47 nz = size(sa)
48 if (size(tu).lt.nz-1 .or. size(rsubrho).lt.nz-1 .or. size(p_mid).lt.nz-1) then
49  tu = gsw_error_code(1,func_name)
50  rsubrho = tu(1)
51  p_mid = tu(1)
52  return
53 end if
54 
55 allocate (dsa(nz-1), sa_mid(nz-1), dct(nz-1), ct_mid(nz-1), dp(nz-1))
56 allocate (alpha_mid(nz-1), beta_mid(nz-1))
57 
58 forall (k = 1: nz-1)
59  dsa(k) = (sa(k) - sa(k+1))
60  sa_mid(k) = 0.5_r8*(sa(k) + sa(k+1))
61  dct(k) = (ct(k) - ct(k+1))
62  ct_mid(k) = 0.5_r8*(ct(k) + ct(k+1))
63  dp(k) = (p(k) - p(k+1))
64  p_mid(k) = 0.5_r8*(p(k) + p(k+1))
65 end forall
66 
67 call gsw_specvol_alpha_beta(sa_mid,ct_mid,p_mid(1:nz-1), &
68  alpha=alpha_mid,beta=beta_mid)
69 
70 tu(1:nz-1) = rad2deg*atan2((alpha_mid*dct + beta_mid*dsa), &
71  (alpha_mid*dct - beta_mid*dsa))
72 
73 where (dsa .ne. 0.0_r8)
74  rsubrho = (alpha_mid*dct)/(beta_mid*dsa)
75 elsewhere
76  rsubrho = gsw_error_code(2,func_name)
77 end where
78 
79 deallocate (dsa, sa_mid, dct, ct_mid, dp)
80 deallocate (alpha_mid, beta_mid)
81 
82 return
83 end subroutine
84 
85 !--------------------------------------------------------------------------
pure subroutine gsw_turner_rsubrho(sa, ct, p, tu, rsubrho, p_mid)
elemental real(r8) function, public gsw_error_code(err_num, func_name, error_code)