FV3 Bundle
gsw_rho_first_derivatives_wrt_enthalpy.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental subroutine gsw_rho_first_derivatives_wrt_enthalpy (sa, ct, p, &
3  rho_sa, rho_h)
4 ! =========================================================================
5 !
6 ! Calculates two first-order derivatives of specific volume (v).
7 ! Note that this function uses the using the computationally-efficient
8 ! expression for specific volume (Roquet et al., 2014).
9 !
10 ! SA = Absolute Salinity [ g/kg ]
11 ! CT = Conservative Temperature (ITS-90) [ deg C ]
12 ! p = sea pressure [ dbar ]
13 ! ( i.e. absolute pressure - 10.1325 dbar )
14 !
15 ! rho_SA = The first derivative of rho with respect to
16 ! Absolute Salinity at constant CT & p. [ J/(kg (g/kg)^2) ]
17 ! rho_h = The first derivative of rho with respect to
18 ! SA and CT at constant p. [ J/(kg K(g/kg)) ]
19 !--------------------------------------------------------------------------
20 
21 use gsw_mod_toolbox, only : gsw_specvol
23 
24 use gsw_mod_kinds
25 
26 implicit none
27 
28 real (r8), intent(in) :: sa, ct, p
29 real (r8), intent(out), optional :: rho_sa, rho_h
30 
31 real (r8) :: rec_v2, v_h, v_sa
32 
33 if (present(rho_sa) .and. present(rho_h)) then
34 
35  call gsw_specvol_first_derivatives_wrt_enthalpy(sa,ct,p,v_sa,v_h)
36 
37 else if (present(rho_sa)) then
38 
39  call gsw_specvol_first_derivatives_wrt_enthalpy(sa,ct,p,v_sa=v_sa)
40 
41 else if (present(rho_h)) then
42 
43  call gsw_specvol_first_derivatives_wrt_enthalpy(sa,ct,p,v_h=v_h)
44 
45 end if
46 
47 rec_v2 = (1.0_r8/gsw_specvol(sa,ct,p))**2
48 
49 if (present(rho_sa)) rho_sa = -v_sa*rec_v2
50 
51 if (present(rho_h)) rho_h = -v_h*rec_v2
52 
53 return
54 end subroutine
55 
56 !--------------------------------------------------------------------------
elemental subroutine gsw_rho_first_derivatives_wrt_enthalpy(sa, ct, p, rho_sa, rho_h)