FV3 Bundle
gsw_melting_ice_sa_ct_ratio.f90
Go to the documentation of this file.
1 !==========================================================================
2 elemental function gsw_melting_ice_sa_ct_ratio (sa, ct, p, t_ih)
3 !==========================================================================
4 !
5 ! Calculates the ratio of SA to CT changes when ice melts into seawater.
6 ! It is assumed that a small mass of ice melts into an infinite mass of
7 ! seawater. Because of the infinite mass of seawater, the ice will always
8 ! melt.
9 !
10 ! The output, melting_seaice_SA_CT_ratio, is dSA/dCT rather than dCT/dSA.
11 ! This is done so that when SA = 0, the output, dSA/dCT is zero whereas
12 ! dCT/dSA would be infinite.
13 !
14 ! SA = Absolute Salinity of seawater [ g/kg ]
15 ! CT = Conservative Temperature of seawater (ITS-90) [ deg C ]
16 ! p = sea pressure at which the melting occurs [ dbar ]
17 ! ( i.e. absolute pressure - 10.1325d0 dbar )
18 ! t_Ih = the in-situ temperature of the ice (ITS-90) [ deg C ]
19 !
20 ! melting_ice_SA_CT_ratio = the ratio of SA to CT changes when ice melts
21 ! into a large mass of seawater
22 ! [ g kg^-1 K^-1 ]
23 !--------------------------------------------------------------------------
24 
28 
30 
31 use gsw_mod_kinds
32 
33 implicit none
34 
35 real (r8), intent(in) :: sa, ct, p, t_ih
36 
37 real (r8) :: gsw_melting_ice_sa_ct_ratio
38 
39 real (r8) :: ctf, h, h_ih, tf, h_hat_sa, h_hat_ct
40 
41 real (r8), parameter :: saturation_fraction = 0.0_r8
42 
43 character (*), parameter :: func_name = "gsw_melting_ice_sa_ct_ratio"
44 
45 ctf = gsw_ct_freezing_exact(sa,p,saturation_fraction)
46 if (ct .lt. ctf) then
47  ! the seawater ct input is below the freezing temperature
49  return
50 end if
51 
52 tf = gsw_t_freezing_exact(0.0_r8,p,saturation_fraction)
53 if (t_ih .gt. tf) then
54  ! t_ih exceeds the freezing temperature at sa = 0
56  return
57 end if
58 
59 h = gsw_enthalpy_ct_exact(sa,ct,p)
60 h_ih = gsw_enthalpy_ice(t_ih,p)
61 call gsw_enthalpy_first_derivatives_ct_exact(sa,ct,p,h_hat_sa,h_hat_ct)
62  ! Note that h_hat_CT is equal to cp0*(273.15 + t)/(273.15 + pt0)
63 
64 gsw_melting_ice_sa_ct_ratio = sa*h_hat_ct/(h - h_ih - sa*h_hat_sa)
65 
66 return
67 end function
68 
69 !--------------------------------------------------------------------------
elemental real(r8) function gsw_melting_ice_sa_ct_ratio(sa, ct, p, t_ih)
real(r8), parameter, public gsw_error_limit
elemental real(r8) function, public gsw_error_code(err_num, func_name, error_code)