FV3 Bundle
gsw_specvol_second_derivatives_wrt_enthalpy.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental subroutine gsw_specvol_second_derivatives_wrt_enthalpy (sa, ct, &
3  p, v_sa_sa, v_sa_h, v_h_h, iflag)
4 ! =========================================================================
5 !
6 ! Calculates three first-order derivatives of specific volume (v) with
7 ! respect to enthalpy. Note that this function uses the using the
8 ! computationally-efficient expression for specific volume
9 ! (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 ! v_SA_SA = The second-order derivative of specific volume with respect to
17 ! Absolute Salinity at constant h & p. [ J/(kg (g/kg)^2) ]
18 ! v_SA_h = The second-order derivative of specific volume with respect to
19 ! SA and h at constant p. [ J/(kg K(g/kg)) ]
20 ! v_h_h = The second-order derivative with respect to h at
21 ! constant SA & p.
22 !--------------------------------------------------------------------------
23 
28 
29 use gsw_mod_kinds
30 
31 implicit none
32 
33 real (r8), intent(in) :: sa, ct, p
34 integer, intent(in), optional :: iflag
35 real (r8), intent(out), optional :: v_sa_sa, v_sa_h, v_h_h
36 
37 logical :: flags(3)
38 real (r8) :: h_ct, h_ct_ct, h_sa, h_sa_ct, h_sa_sa, rec_h_ct, v_h_h_part
39 real (r8) :: rec_h_ct2, v_ct, vct_ct_ct, vct_sa_ct, vct_sa_sa, v_sa_h_part
40 
41 if (present(iflag)) then
42  flags(1) = present(v_sa_sa) .and. btest(iflag,1)
43  flags(2) = present(v_sa_h) .and. btest(iflag,2)
44  flags(3) = present(v_h_h) .and. btest(iflag,3)
45 else
46  flags(1) = present(v_sa_sa)
47  flags(2) = present(v_sa_h)
48  flags(3) = present(v_h_h)
49 end if
50 
51 call gsw_specvol_first_derivatives(sa,ct,p,v_ct=v_ct)
52 
53 if (flags(1) .or. flags(2)) then
54  call gsw_enthalpy_first_derivatives(sa,ct,p,h_sa,h_ct)
55 else
56  call gsw_enthalpy_first_derivatives(sa,ct,p,h_ct=h_ct)
57 end if
58 
59 if (flags(1)) then
60  call gsw_specvol_second_derivatives(sa,ct,p,vct_sa_sa,vct_sa_ct,vct_ct_ct)
61 else if (flags(2)) then
62  call gsw_specvol_second_derivatives(sa,ct,p,v_sa_ct=vct_sa_ct, &
63  v_ct_ct=vct_ct_ct)
64 else
65  call gsw_specvol_second_derivatives(sa,ct,p,v_ct_ct=vct_ct_ct)
66 end if
67 
68 if (flags(1)) then
69  call gsw_enthalpy_second_derivatives(sa,ct,p,h_sa_sa,h_sa_ct,h_ct_ct)
70 else if (flags(2)) then
71  call gsw_enthalpy_second_derivatives(sa,ct,p,h_sa_ct=h_sa_ct,h_ct_ct=h_ct_ct)
72 else
73  call gsw_enthalpy_second_derivatives(sa,ct,p,h_ct_ct=h_ct_ct)
74 end if
75 
76 rec_h_ct = 1.0_r8/h_ct
77 rec_h_ct2 = rec_h_ct**2.0_r8
78 
79 v_h_h_part = (vct_ct_ct*h_ct - h_ct_ct*v_ct)*(rec_h_ct2*rec_h_ct)
80 
81 if (flags(3)) v_h_h = v_h_h_part
82 
83 if (flags(1) .or. flags(2)) then
84 
85  v_sa_h_part = (vct_sa_ct*h_ct - v_ct*h_sa_ct)*rec_h_ct2 - h_sa*v_h_h_part
86 
87  if (flags(2)) v_sa_h = v_sa_h_part
88 
89  if (flags(1)) v_sa_sa = vct_sa_sa - (h_ct*(vct_sa_ct*h_sa &
90  - v_ct*h_sa_sa) + v_ct*h_sa*h_sa_ct)*rec_h_ct2 - h_sa*v_sa_h_part
91 end if
92 
93 return
94 end subroutine
95 
96 !--------------------------------------------------------------------------
elemental subroutine gsw_specvol_second_derivatives_wrt_enthalpy(sa, ct, p, v_sa_sa, v_sa_h, v_h_h, iflag)