FV3 Bundle
gsw_pt_second_derivatives.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental subroutine gsw_pt_second_derivatives (sa, ct, pt_sa_sa, &
3  pt_sa_ct, pt_ct_ct)
4 ! =========================================================================
5 !
6 ! Calculates the following three second-order derivatives of potential
7 ! temperature (the regular potential temperature which has a reference
8 ! sea pressure of 0 dbar),
9 ! (1) pt_SA_SA, the second derivative with respect to Absolute Salinity
10 ! at constant Conservative Temperature,
11 ! (2) pt_SA_CT, the derivative with respect to Conservative Temperature
12 ! and Absolute Salinity, and
13 ! (3) pt_CT_CT, the second derivative with respect to Conservative
14 ! Temperature at constant Absolute Salinity.
15 !
16 ! SA = Absolute Salinity [ g/kg ]
17 ! CT = Conservative Temperature (ITS-90) [ deg C ]
18 !
19 ! pt_SA_SA = The second derivative of potential temperature (the
20 ! regular potential temperature which has reference sea
21 ! pressure of 0 dbar) with respect to Absolute Salinity
22 ! at constant Conservative Temperature.
23 ! pt_SA_SA has units of: [ K/((g/kg)^2) ]
24 ! pt_SA_CT = The derivative of potential temperature with respect
25 ! to Absolute Salinity and Conservative Temperature.
26 ! pt_SA_CT has units of: [ 1/(g/kg) ]
27 ! pt_CT_CT = The second derivative of potential temperature (the
28 ! regular one with p_ref = 0 dbar) with respect to
29 ! Conservative Temperature at constant SA.
30 ! pt_CT_CT has units of: [ 1/K ]
31 !--------------------------------------------------------------------------
32 
34 
35 use gsw_mod_kinds
36 
37 implicit none
38 
39 real (r8), intent(in) :: sa, ct
40 real (r8), intent(out), optional :: pt_sa_sa, pt_sa_ct, pt_ct_ct
41 
42 real (r8) :: ct_l, ct_u, pt_ct_l, pt_ct_u, pt_sa_l, pt_sa_u, sa_l, sa_u
43 
44 real (r8), parameter :: dct = 1e-2_r8, dsa = 1e-3_r8
45 
46 if (present(pt_sa_sa)) then
47 
48  sa_l = max(sa - dsa, 0.0_r8)
49  sa_u = sa + dsa
50 
51  call gsw_pt_first_derivatives(sa_l,ct,pt_sa=pt_sa_l)
52  call gsw_pt_first_derivatives(sa_u,ct,pt_sa=pt_sa_u)
53 
54  pt_sa_sa = (pt_sa_u - pt_sa_l)/(sa_u - sa_l)
55 
56 end if
57 
58 if (present(pt_sa_ct) .or. present(pt_ct_ct)) then
59 
60  ct_l = ct - dct
61  ct_u = ct + dct
62 
63  if (present(pt_sa_ct) .and. present(pt_ct_ct)) then
64 
65  call gsw_pt_first_derivatives(sa,ct_l,pt_sa_l,pt_ct_l)
66  call gsw_pt_first_derivatives(sa,ct_u,pt_sa_u,pt_ct_u)
67 
68  pt_sa_ct = (pt_sa_u - pt_sa_l)/(ct_u - ct_l)
69  pt_ct_ct = (pt_ct_u - pt_ct_l)/(ct_u - ct_l)
70 
71  else if (present(pt_sa_ct) .and. .not. present(pt_ct_ct)) then
72 
73  call gsw_pt_first_derivatives(sa,ct_l,pt_sa=pt_sa_l)
74  call gsw_pt_first_derivatives(sa,ct_u,pt_sa=pt_sa_u)
75 
76  pt_sa_ct = (pt_sa_u - pt_sa_l)/(ct_u - ct_l)
77 
78  else if (.not. present(pt_sa_ct) .and. present(pt_ct_ct)) then
79 
80  call gsw_pt_first_derivatives(sa,ct_l,pt_ct=pt_ct_l)
81  call gsw_pt_first_derivatives(sa,ct_u,pt_ct=pt_ct_u)
82 
83  pt_ct_ct = (pt_ct_u - pt_ct_l)/(ct_u - ct_l)
84 
85  end if
86 
87 end if
88 
89 return
90 end subroutine
91 
92 !--------------------------------------------------------------------------
elemental subroutine gsw_pt_second_derivatives(sa, ct, pt_sa_sa, pt_sa_ct, pt_ct_ct)
#define max(a, b)
Definition: mosaic_util.h:33