FV3 Bundle
gsw_nsquared_min_const_t.f90
Go to the documentation of this file.
1 !==========================================================================
2 pure subroutine gsw_nsquared_min_const_t (sa, t, p, lat, n2, n2_p, &
3  n2_specvol, n2_alpha, n2_beta, dsa, dct, dp, n2_beta_ratio)
4 !==========================================================================
5 !
6 ! Calculates the minimum buoyancy frequency squared (N^2) (i.e. the
7 ! Brunt-Vaisala frequency squared) between two bottles from the equation,
8 !
9 ! 2 2 beta.dSA - alpha.dCT
10 ! N = g . -------------------------
11 ! specvol_local.dP
12 !
13 ! The pressure increment, dP, in the above formula is in Pa, so that it is
14 ! 10^4 times the pressure increment dp in dbar.
15 !
16 ! SA = Absolute Salinity [ g/kg ]
17 ! t = In-situ temperature (ITS-90) [ deg C ]
18 ! p = Sea pressure [ dbar ]
19 ! ( i.e. absolute pressure - 10.1325 dbar )
20 ! lat = Latitude in decimal degrees north [ -90 ... +90 ]
21 ! Note: If lat is outside this range, a default gravitational
22 ! acceleration of 9.7963 m/s^2 (Griffies, 2004) will be applied.
23 !
24 ! N2 = minimum Brunt-Vaisala Frequency squared [ 1/s^2 ]
25 ! N2_p = pressure of minimum N2 [ dbar ]
26 ! N2_specvol = specific volume at the minimum N2 [ kg/m3 ]
27 ! N2_alpha = thermal expansion coefficient with respect [ 1/K ]
28 ! to Conservative Temperature at the minimum N2
29 ! N2_beta = saline contraction coefficient at constant [ kg/g ]
30 ! Conservative Temperature at the minimum N2
31 ! dSA = difference in salinity between bottles [ g/kg ]
32 ! dCT = difference in Conservative Temperature between [ deg C ]
33 ! bottles
34 ! dp = difference in pressure between bottles [ dbar ]
35 ! N2_beta_ratio = ratio of the saline contraction [ unitless ]
36 ! coefficient at constant Conservative Temperature to
37 ! the saline contraction coefficient at constant in-situ
38 ! temperature at the minimum N2
39 !
40 !==========================================================================
41 
44 
45 use gsw_mod_kinds
46 
48 
49 implicit none
50 
51 real (r8), intent(in) :: sa(:), t(:), p(:), lat
52 real (r8), intent(out) :: n2(:), n2_p(:), n2_specvol(:), n2_alpha(:)
53 real (r8), intent(out) :: n2_beta(:), dsa(:), dct(:), dp(:)
54 real (r8), intent(out) :: n2_beta_ratio(:)
55 
56 integer :: i, ideep, ishallow, mp
57 
58 real (r8) :: n2_deep, n2_shallow
59 real (r8), allocatable :: alpha(:), beta(:), specvol(:), grav2(:)
60 real (r8), allocatable :: ct(:), beta_ratio(:)
61 
62 mp = size(sa)
63 allocate(grav2(mp),specvol(mp),alpha(mp),beta(mp),ct(mp),beta_ratio(mp))
64 
65 if (lat .lt. -90.0_r8 .or. lat .gt. +90.0_r8) then
66  grav2 = (/ (9.7963_r8**2, i=1,mp) /)
67 else
68  grav2 = gsw_grav(lat,p)**2
69 end if
70 
71 ct = gsw_ct_from_t(sa,t,p)
72 
73 dp = p(2:mp) - p(1:mp-1)
74 dsa = sa(2:mp) - sa(1:mp-1)
75 dct = ct(2:mp) - ct(1:mp-1)
76 
77 call gsw_specvol_alpha_beta(sa,ct,p,specvol,alpha,beta)
78 
79 beta_ratio = beta/gsw_beta_const_t_exact(sa,t,p)
80 
81 ishallow = 1
82 ideep = 2
83 do i = 1, mp-1
84  n2_shallow = grav2(ishallow)/(specvol(ishallow)*db2pa*dp(i))* &
85  (beta(ishallow)*dsa(i) - alpha(ishallow)*dct(i))
86  n2_deep = grav2(ideep)/(specvol(ideep)*db2pa*dp(i))* &
87  (beta(ideep)*dsa(i) - alpha(ideep)*dct(i))
88  if (n2_shallow .lt. n2_deep) then
89  n2(i) = n2_shallow
90  n2_p(i) = p(ishallow)
91  n2_specvol(i) = specvol(ishallow)
92  n2_alpha(i) = alpha(ishallow)
93  n2_beta(i) = beta(ishallow)
94  n2_beta_ratio(i) = beta_ratio(ishallow)
95  else
96  n2(i) = n2_deep
97  n2_p(i) = p(ideep)
98  n2_specvol(i) = specvol(ideep)
99  n2_alpha(i) = alpha(ideep)
100  n2_beta(i) = beta(ideep)
101  n2_beta_ratio(i) = beta_ratio(ideep)
102  end if
103  ishallow = ishallow + 1
104  ideep = ideep + 1
105 end do
106 
107 return
108 end subroutine
109 
110 !--------------------------------------------------------------------------
pure subroutine gsw_nsquared_min_const_t(sa, t, p, lat, n2, n2_p, n2_specvol, n2_alpha, n2_beta, dsa, dct, dp, n2_beta_ratio)