FV3 Bundle
gsw_pt_first_derivatives.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental subroutine gsw_pt_first_derivatives (sa, ct, pt_sa, pt_ct)
3 ! =========================================================================
4 !
5 ! Calculates the following two partial derivatives of potential temperature
6 ! (the regular potential temperature whose reference sea pressure is 0 dbar)
7 ! (1) pt_SA, the derivative with respect to Absolute Salinity at
8 ! constant Conservative Temperature, and
9 ! (2) pt_CT, the derivative with respect to Conservative Temperature at
10 ! constant Absolute Salinity.
11 !
12 ! SA = Absolute Salinity [ g/kg ]
13 ! CT = Conservative Temperature (ITS-90) [ deg C ]
14 !
15 ! pt_SA = The derivative of potential temperature with respect to
16 ! Absolute Salinity at constant Conservative Temperature.
17 ! [ K/(g/kg)]
18 ! pt_CT = The derivative of potential temperature with respect to
19 ! Conservative Temperature at constant Absolute Salinity.
20 ! pt_CT is dimensionless. [ unitless ]
21 !--------------------------------------------------------------------------
22 
24 
26 
27 use gsw_mod_kinds
28 
29 implicit none
30 
31 real (r8), intent(in) :: sa, ct
32 real (r8), intent(out), optional :: pt_sa, pt_ct
33 
34 real (r8) :: abs_pt, ct_pt, ct_sa, pt
35 
36 integer, parameter :: n0=0, n1=1, n2=2
37 real (r8), parameter :: pr0 = 0.0_r8
38 
39 pt = gsw_pt_from_ct(sa,ct)
40 abs_pt = (gsw_t0 + pt)
41 
42 ct_pt = -(abs_pt*gsw_gibbs(n0,n2,n0,sa,pt,pr0))/gsw_cp0
43 
44 if (present(pt_sa)) then
45 
46  ct_sa = (gsw_gibbs(n1,n0,n0,sa,pt,pr0) -&
47  abs_pt*gsw_gibbs(n1,n1,n0,sa,pt,pr0))/gsw_cp0
48 
49  pt_sa = -ct_sa/ct_pt
50 
51 end if
52 
53 if (present(pt_ct)) pt_ct = 1.0_r8/ct_pt
54 
55 return
56 end subroutine
57 
58 !--------------------------------------------------------------------------
elemental subroutine gsw_pt_first_derivatives(sa, ct, pt_sa, pt_ct)