FV3 Bundle
gsw_specvol_alpha_beta.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental subroutine gsw_specvol_alpha_beta (sa, ct, p, specvol, alpha, &
3  beta)
4 !==========================================================================
5 !
6 ! Calculates specific volume, the appropiate thermal expansion coefficient
7 ! and the appropriate saline contraction coefficient of seawater from
8 ! Absolute Salinity and Conservative Temperature. This function uses the
9 ! computationally-efficient expression for specific volume in terms of
10 ! SA, CT and p (Roquet et al., 2014).
11 !
12 ! SA = Absolute Salinity [ g/kg ]
13 ! CT = Conservative Temperature (ITS-90) [ deg C ]
14 ! p = sea pressure [ dbar ]
15 ! ( i.e. absolute pressure - 10.1325 dbar )
16 !
17 ! specvol = specific volume [ m/kg ]
18 ! alpha = thermal expansion coefficient [ 1/K ]
19 ! with respect to Conservative Temperature
20 ! beta = saline (i.e. haline) contraction [ kg/g ]
21 ! coefficient at constant Conservative Temperature
22 !--------------------------------------------------------------------------
23 
25 
27 
28 use gsw_mod_kinds
29 
30 implicit none
31 
32 real (r8), intent(in) :: sa, ct, p
33 real (r8), intent(out), optional :: specvol, alpha, beta
34 
35 real (r8) :: v, v_ct, v_sa_part, xs, ys, z
36 
37 xs = sqrt(gsw_sfac*sa + offset)
38 ys = ct*0.025_r8
39 z = p*1e-4_r8
40 
41 v = v000 + xs*(v010 + xs*(v020 + xs*(v030 + xs*(v040 + xs*(v050 &
42  + v060*xs))))) + ys*(v100 + xs*(v110 + xs*(v120 + xs*(v130 + xs*(v140 &
43  + v150*xs)))) + ys*(v200 + xs*(v210 + xs*(v220 + xs*(v230 + v240*xs))) &
44  + ys*(v300 + xs*(v310 + xs*(v320 + v330*xs)) + ys*(v400 + xs*(v410 &
45  + v420*xs) + ys*(v500 + v510*xs + v600*ys))))) + z*(v001 + xs*(v011 &
46  + xs*(v021 + xs*(v031 + xs*(v041 + v051*xs)))) + ys*(v101 + xs*(v111 &
47  + xs*(v121 + xs*(v131 + v141*xs))) + ys*(v201 + xs*(v211 + xs*(v221 &
48  + v231*xs)) + ys*(v301 + xs*(v311 + v321*xs) + ys*(v401 + v411*xs &
49  + v501*ys)))) + z*(v002 + xs*(v012 + xs*(v022 + xs*(v032 + v042*xs))) &
50  + ys*(v102 + xs*(v112 + xs*(v122 + v132*xs)) + ys*(v202 + xs*(v212 &
51  + v222*xs) + ys*(v302 + v312*xs + v402*ys))) + z*(v003 + xs*(v013 &
52  + v023*xs) + ys*(v103 + v113*xs + v203*ys) + z*(v004 + v014*xs + v104*ys &
53  + z*(v005 + v006*z)))))
54 
55 if (present(specvol)) specvol = v
56 
57 if (present(alpha)) then
58 
59  v_ct = a000 + xs*(a100 + xs*(a200 + xs*(a300 + xs*(a400 + a500*xs)))) &
60  + ys*(a010 + xs*(a110 + xs*(a210 + xs*(a310 + a410*xs))) &
61  + ys*(a020 + xs*(a120 + xs*(a220 + a320*xs)) + ys*(a030 &
62  + xs*(a130 + a230*xs) + ys*(a040 + a140*xs + a050*ys )))) &
63  + z*(a001 + xs*(a101 + xs*(a201 + xs*(a301 + a401*xs))) &
64  + ys*(a011 + xs*(a111 + xs*(a211 + a311*xs)) + ys*(a021 &
65  + xs*(a121 + a221*xs) + ys*(a031 + a131*xs + a041*ys))) &
66  + z*(a002 + xs*(a102 + xs*(a202 + a302*xs)) + ys*(a012 &
67  + xs*(a112 + a212*xs) + ys*(a022 + a122*xs + a032*ys)) &
68  + z*(a003 + a103*xs + a013*ys + a004*z)))
69 
70  alpha = 0.025_r8*v_ct/v
71 
72 end if
73 
74 if (present(beta)) then
75 
76  v_sa_part = b000 + xs*(b100 + xs*(b200 + xs*(b300 + xs*(b400 + b500*xs)))) &
77  + ys*(b010 + xs*(b110 + xs*(b210 + xs*(b310 + b410*xs))) &
78  + ys*(b020 + xs*(b120 + xs*(b220 + b320*xs)) + ys*(b030 &
79  + xs*(b130 + b230*xs) + ys*(b040 + b140*xs + b050*ys)))) &
80  + z*(b001 + xs*(b101 + xs*(b201 + xs*(b301 + b401*xs))) &
81  + ys*(b011 + xs*(b111 + xs*(b211 + b311*xs)) + ys*(b021 &
82  + xs*(b121 + b221*xs) + ys*(b031 + b131*xs + b041*ys))) &
83  + z*(b002 + xs*(b102 + xs*(b202 + b302*xs))+ ys*(b012 &
84  + xs*(b112 + b212*xs) + ys*(b022 + b122*xs + b032*ys)) &
85  + z*(b003 + b103*xs + b013*ys + b004*z)))
86 
87  beta = -v_sa_part*0.5_r8*gsw_sfac/(v*xs)
88 
89 end if
90 
91 return
92 end subroutine
93 
94 !--------------------------------------------------------------------------
elemental subroutine gsw_specvol_alpha_beta(sa, ct, p, specvol, alpha, beta)