FV3 Bundle
gsw_alpha_on_beta.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental function gsw_alpha_on_beta (sa, ct, p)
3 !==========================================================================
4 !
5 ! Calculates alpha divided by beta, where alpha is the thermal expansion
6 ! coefficient and beta is the saline contraction coefficient of seawater
7 ! from Absolute Salinity and Conservative Temperature. This function uses
8 ! the computationally-efficient expression for specific volume in terms of
9 ! SA, CT and p (Roquet et al., 2014).
10 !
11 ! SA = Absolute Salinity [ g/kg ]
12 ! CT = Conservative Temperature (ITS-90) [ deg C ]
13 ! p = sea pressure [ dbar ]
14 ! ( i.e. absolute pressure - 10.1325 dbar )
15 !
16 ! alpha_on_beta = thermal expansion coefficient with respect to
17 ! Conservative Temperature divided by the saline
18 ! contraction coefficient at constant Conservative
19 ! Temperature [ kg g^-1 K^-1 ]
20 !--------------------------------------------------------------------------
21 
23 
25 
26 use gsw_mod_kinds
27 
28 implicit none
29 
30 real (r8), intent(in) :: sa, ct, p
31 
32 real (r8) :: gsw_alpha_on_beta
33 
34 real (r8) :: xs, ys, z, v_ct_part, v_sa_part
35 
36 xs = sqrt(gsw_sfac*sa + offset)
37 ys = ct*0.025_r8
38 z = p*1e-4_r8
39 
40 v_ct_part = a000 + xs*(a100 + xs*(a200 + xs*(a300 + xs*(a400 + a500*xs)))) &
41  + ys*(a010 + xs*(a110 + xs*(a210 + xs*(a310 + a410*xs))) &
42  + ys*(a020 + xs*(a120 + xs*(a220 + a320*xs)) + ys*(a030 &
43  + xs*(a130 + a230*xs) + ys*(a040 + a140*xs + a050*ys )))) &
44  + z*(a001 + xs*(a101 + xs*(a201 + xs*(a301 + a401*xs))) &
45  + ys*(a011 + xs*(a111 + xs*(a211 + a311*xs)) + ys*(a021 &
46  + xs*(a121 + a221*xs) + ys*(a031 + a131*xs + a041*ys))) &
47  + z*(a002 + xs*(a102 + xs*(a202 + a302*xs)) + ys*(a012 &
48  + xs*(a112 + a212*xs) + ys*(a022 + a122*xs + a032*ys)) &
49  + z*(a003 + a103*xs + a013*ys + a004*z)))
50 
51 v_sa_part = b000 + xs*(b100 + xs*(b200 + xs*(b300 + xs*(b400 + b500*xs)))) &
52  + ys*(b010 + xs*(b110 + xs*(b210 + xs*(b310 + b410*xs))) &
53  + ys*(b020 + xs*(b120 + xs*(b220 + b320*xs)) + ys*(b030 &
54  + xs*(b130 + b230*xs) + ys*(b040 + b140*xs + b050*ys)))) &
55  + z*(b001 + xs*(b101 + xs*(b201 + xs*(b301 + b401*xs))) &
56  + ys*(b011 + xs*(b111 + xs*(b211 + b311*xs)) + ys*(b021 &
57  + xs*(b121 + b221*xs) + ys*(b031 + b131*xs + b041*ys))) &
58  + z*(b002 + xs*(b102 + xs*(b202 + b302*xs))+ ys*(b012 &
59  + xs*(b112 + b212*xs) + ys*(b022 + b122*xs + b032*ys)) &
60  + z*(b003 + b103*xs + b013*ys + b004*z)))
61 
62 gsw_alpha_on_beta = -(v_ct_part*xs)/(20.0_r8*gsw_sfac*v_sa_part)
63 
64 return
65 end function
66 
67 !--------------------------------------------------------------------------
elemental real(r8) function gsw_alpha_on_beta(sa, ct, p)