FV3 Bundle
gsw_rho_first_derivatives.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental subroutine gsw_rho_first_derivatives (sa, ct, p, drho_dsa, &
3  drho_dct, drho_dp)
4 !==========================================================================
5 !
6 ! Calculates the three (3) partial derivatives of in-situ density with
7 ! respect to Absolute Salinity, Conservative Temperature and pressure.
8 ! Note that the pressure derivative is done with respect to pressure in
9 ! Pa, not dbar. This function uses the computationally-efficient expression
10 ! for specific volume in terms of 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 ! drho_dSA = partial derivatives of density [ kg^2/(g m^3) ]
18 ! with respect to Absolute Salinity
19 ! drho_dCT = partial derivatives of density [ kg/(K m^3) ]
20 ! with respect to Conservative Temperature
21 ! drho_dP = partial derivatives of density [ kg/(Pa m^3) ]
22 ! with respect to pressure in Pa
23 !--------------------------------------------------------------------------
24 
26 
28 
29 use gsw_mod_kinds
30 
31 implicit none
32 
33 real (r8), intent(in) :: sa, ct, p
34 real (r8), intent(out), optional :: drho_dsa, drho_dct, drho_dp
35 
36 real (r8) :: rho2, v_ct, v_p, v_sa, xs, ys, z, v
37 
38 xs = sqrt(gsw_sfac*sa + offset)
39 ys = ct*0.025_r8
40 z = p*1e-4_r8
41 
42 v = v000 + xs*(v010 + xs*(v020 + xs*(v030 + xs*(v040 + xs*(v050 &
43  + v060*xs))))) + ys*(v100 + xs*(v110 + xs*(v120 + xs*(v130 + xs*(v140 &
44  + v150*xs)))) + ys*(v200 + xs*(v210 + xs*(v220 + xs*(v230 + v240*xs))) &
45  + ys*(v300 + xs*(v310 + xs*(v320 + v330*xs)) + ys*(v400 + xs*(v410 &
46  + v420*xs) + ys*(v500 + v510*xs + v600*ys))))) + z*(v001 + xs*(v011 &
47  + xs*(v021 + xs*(v031 + xs*(v041 + v051*xs)))) + ys*(v101 + xs*(v111 &
48  + xs*(v121 + xs*(v131 + v141*xs))) + ys*(v201 + xs*(v211 + xs*(v221 &
49  + v231*xs)) + ys*(v301 + xs*(v311 + v321*xs) + ys*(v401 + v411*xs &
50  + v501*ys)))) + z*(v002 + xs*(v012 + xs*(v022 + xs*(v032 + v042*xs))) &
51  + ys*(v102 + xs*(v112 + xs*(v122 + v132*xs)) + ys*(v202 + xs*(v212 &
52  + v222*xs) + ys*(v302 + v312*xs + v402*ys))) + z*(v003 + xs*(v013 &
53  + v023*xs) + ys*(v103 + v113*xs + v203*ys) + z*(v004 + v014*xs + v104*ys &
54  + z*(v005 + v006*z)))))
55 
56 rho2 = (1.0_r8/v)**2
57 
58 if (present(drho_dsa)) then
59 
60  v_sa = b000 + xs*(b100 + xs*(b200 + xs*(b300 + xs*(b400 + b500*xs)))) &
61  + ys*(b010 + xs*(b110 + xs*(b210 + xs*(b310 + b410*xs))) &
62  + ys*(b020 + xs*(b120 + xs*(b220 + b320*xs)) + ys*(b030 &
63  + xs*(b130 + b230*xs) + ys*(b040 + b140*xs + b050*ys)))) &
64  + z*(b001 + xs*(b101 + xs*(b201 + xs*(b301 + b401*xs))) &
65  + ys*(b011 + xs*(b111 + xs*(b211 + b311*xs)) + ys*(b021 &
66  + xs*(b121 + b221*xs) + ys*(b031 + b131*xs + b041*ys))) &
67  + z*(b002 + xs*(b102 + xs*(b202 + b302*xs))+ ys*(b012 &
68  + xs*(b112 + b212*xs) + ys*(b022 + b122*xs + b032*ys)) &
69  + z*(b003 + b103*xs + b013*ys + b004*z)))
70 
71  drho_dsa = -rho2*0.5_r8*gsw_sfac*v_sa/xs
72 
73 end if
74 
75 if (present(drho_dct)) then
76 
77  v_ct = a000 + xs*(a100 + xs*(a200 + xs*(a300 + xs*(a400 + a500*xs)))) &
78  + ys*(a010 + xs*(a110 + xs*(a210 + xs*(a310 + a410*xs))) &
79  + ys*(a020 + xs*(a120 + xs*(a220 + a320*xs)) + ys*(a030 &
80  + xs*(a130 + a230*xs) + ys*(a040 + a140*xs + a050*ys )))) &
81  + z*(a001 + xs*(a101 + xs*(a201 + xs*(a301 + a401*xs))) &
82  + ys*(a011 + xs*(a111 + xs*(a211 + a311*xs)) + ys*(a021 &
83  + xs*(a121 + a221*xs) + ys*(a031 + a131*xs + a041*ys))) &
84  + z*(a002 + xs*(a102 + xs*(a202 + a302*xs)) + ys*(a012 &
85  + xs*(a112 + a212*xs) + ys*(a022 + a122*xs + a032*ys)) &
86  + z*(a003 + a103*xs + a013*ys + a004*z)))
87 
88  drho_dct = -rho2*0.025_r8*v_ct
89 
90 end if
91 
92 if (present(drho_dp)) then
93 
94  v_p = c000 + xs*(c100 + xs*(c200 + xs*(c300 + xs*(c400 + c500*xs)))) &
95  + ys*(c010 + xs*(c110 + xs*(c210 + xs*(c310 + c410*xs))) + ys*(c020 &
96  + xs*(c120 + xs*(c220 + c320*xs)) + ys*(c030 + xs*(c130 + c230*xs) &
97  + ys*(c040 + c140*xs + c050*ys)))) + z*(c001 + xs*(c101 + xs*(c201 &
98  + xs*(c301 + c401*xs))) + ys*(c011 + xs*(c111 + xs*(c211 + c311*xs)) &
99  + ys*(c021 + xs*(c121 + c221*xs) + ys*(c031 + c131*xs + c041*ys))) &
100  + z*(c002 + xs*(c102 + c202*xs) + ys*(c012 + c112*xs + c022*ys) &
101  + z*(c003 + c103*xs + c013*ys + z*(c004 + c005*z))))
102 
103  drho_dp = -rho2*1e-4_r8*pa2db*v_p
104 
105 end if
106 
107 return
108 end subroutine
109 
110 !--------------------------------------------------------------------------
elemental subroutine gsw_rho_first_derivatives(sa, ct, p, drho_dsa, drho_dct, drho_dp)