FV3 Bundle
gsw_entropy_second_derivatives.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental subroutine gsw_entropy_second_derivatives (sa, ct, eta_sa_sa, &
3  eta_sa_ct, eta_ct_ct)
4 ! =========================================================================
5 !
6 ! Calculates the following three second-order partial derivatives of
7 ! specific entropy (eta)
8 ! (1) eta_SA_SA, the second derivative with respect to Absolute
9 ! Salinity at constant Conservative Temperature, and
10 ! (2) eta_SA_CT, the derivative with respect to Absolute Salinity and
11 ! Conservative Temperature.
12 ! (3) eta_CT_CT, the second derivative with respect to Conservative
13 ! Temperature at constant Absolute Salinity.
14 !
15 ! SA = Absolute Salinity [ g/kg ]
16 ! CT = Conservative Temperature (ITS-90) [ deg C ]
17 !
18 ! eta_SA_SA = The second derivative of specific entropy with respect
19 ! to Absolute Salinity (in units of g kg^-1) at constant
20 ! Conservative Temperature.
21 ! eta_SA_SA has units of: [ J/(kg K(g/kg)^2)]
22 ! eta_SA_CT = The second derivative of specific entropy with respect
23 ! to Conservative Temperature at constant Absolute
24 ! Salinity. eta_SA_CT has units of: [ J/(kg (g/kg) K^2) ]
25 ! eta_CT_CT = The second derivative of specific entropy with respect
26 ! to Conservative Temperature at constant Absolute
27 ! Salinity. eta_CT_CT has units of: [ J/(kg K^3) ]
28 !--------------------------------------------------------------------------
29 
31 
33 
34 use gsw_mod_kinds
35 
36 implicit none
37 
38 real (r8), intent(in) :: sa, ct
39 real (r8), intent(out), optional :: eta_sa_sa, eta_sa_ct, eta_ct_ct
40 
41 real (r8) :: abs_pt, ct_pt, ct_sa, pt, ct_ct
42 
43 integer, parameter :: n0=0, n1=1, n2=2
44 real (r8), parameter :: pr0 = 0.0_r8
45 
46 pt = gsw_pt_from_ct(sa,ct)
47 abs_pt = gsw_t0 + pt
48 
49 ct_pt = -(abs_pt*gsw_gibbs(n0,n2,n0,sa,pt,pr0))/gsw_cp0
50 
51 ct_ct = -gsw_cp0/(ct_pt*abs_pt*abs_pt)
52 
53 if (present(eta_sa_ct) .or. present(eta_sa_sa)) then
54 
55  ct_sa = (gsw_gibbs(n1,n0,n0,sa,pt,pr0) - &
56  (abs_pt*gsw_gibbs(n1,n1,n0,sa,pt,pr0)))/gsw_cp0
57 
58  if (present(eta_sa_ct)) eta_sa_ct = -ct_sa*ct_ct
59 
60  if (present(eta_sa_sa)) eta_sa_sa = -gsw_gibbs(n2,n0,n0,sa,pt,pr0)/abs_pt +&
61  ct_sa*ct_sa*ct_ct
62 end if
63 
64 if (present(eta_ct_ct)) eta_ct_ct = ct_ct
65 
66 return
67 end subroutine
68 
69 !--------------------------------------------------------------------------
elemental subroutine gsw_entropy_second_derivatives(sa, ct, eta_sa_sa, eta_sa_ct, eta_ct_ct)