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