FV3 Bundle
gsw_nsquared_min.f90
Go to the documentation of this file.
1 !==========================================================================
2 pure subroutine gsw_nsquared_min (sa, ct, p, lat, n2, n2_p, &
3  n2_specvol, n2_alpha, n2_beta, dsa, dct, dp)
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 ! CT = Conservative 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 !
22 ! N2 = minimum Brunt-Vaisala Frequency squared [ 1/s^2 ]
23 ! N2_p = pressure of minimum N2 [ dbar ]
24 ! N2_specvol = specific volume at the minimum N2 [ kg/m3 ]
25 ! N2_alpha = thermal expansion coefficient with respect [ 1/K ]
26 ! to Conservative Temperature at the minimum N2
27 ! N2_beta = saline contraction coefficient at constant [ kg/g ]
28 ! Conservative Temperature at the minimum N2
29 ! dSA = difference in salinity between bottles [ g/kg ]
30 ! dCT = difference in Conservative Temperature between [ deg C ]
31 ! bottles
32 ! dp = difference in pressure between bottles [ dbar ]
33 !
34 !==========================================================================
35 
37 
39 
40 use gsw_mod_kinds
41 
43 
44 implicit none
45 
46 real (r8), intent(in) :: sa(:), ct(:), p(:), lat(:)
47 real (r8), intent(out) :: n2(:), n2_p(:), n2_specvol(:), n2_alpha(:)
48 real (r8), intent(out) :: n2_beta(:), dsa(:), dct(:), dp(:)
49 
50 integer :: i, ideep, ishallow, mp
51 
52 real (r8) :: n2_deep, n2_shallow
53 real (r8), allocatable :: alpha(:), beta(:), specvol(:), grav2(:)
54 
55 character (*), parameter :: func_name = "gsw_nsquared_min"
56 
57 mp = size(sa)
58 if (size(n2).lt.mp-1 .or. size(n2_p).lt.mp-1 .or. &
59  size(n2_specvol).lt.mp-1 .or. size(n2_alpha).lt.mp-1 .or. &
60  size(n2_beta).lt.mp-1 .or. size(dsa).lt.mp-1 .or. &
61  size(dct).lt.mp-1 .or. size(dp).lt.mp-1) then
62  n2 = gsw_error_code(1,func_name)
63  n2_p = n2(1)
64  n2_specvol = n2(1)
65  n2_alpha = n2(1)
66  n2_beta = n2(1)
67  dsa = n2(1)
68  dct = n2(1)
69  dp = n2(1)
70  return
71 end if
72 
73 allocate(grav2(mp),specvol(mp),alpha(mp),beta(mp))
74 
75 grav2 = gsw_grav(lat,p)**2
76 
77 dp(1:mp-1) = p(2:mp) - p(1:mp-1)
78 dsa(1:mp-1) = sa(2:mp) - sa(1:mp-1)
79 dct(1:mp-1) = ct(2:mp) - ct(1:mp-1)
80 
81 call gsw_specvol_alpha_beta(sa,ct,p,specvol,alpha,beta)
82 
83 ishallow = 1
84 ideep = 2
85 do i = 1, mp-1
86  n2_shallow = grav2(ishallow)/(specvol(ishallow)*db2pa*dp(i))* &
87  (beta(ishallow)*dsa(i) - alpha(ishallow)*dct(i))
88  n2_deep = grav2(ideep)/(specvol(ideep)*db2pa*dp(i))* &
89  (beta(ideep)*dsa(i) - alpha(ideep)*dct(i))
90  if (n2_shallow .lt. n2_deep) then
91  n2(i) = n2_shallow
92  n2_p(i) = p(ishallow)
93  n2_specvol(i) = specvol(ishallow)
94  n2_alpha(i) = alpha(ishallow)
95  n2_beta(i) = beta(ishallow)
96  else
97  n2(i) = n2_deep
98  n2_p(i) = p(ideep)
99  n2_specvol(i) = specvol(ideep)
100  n2_alpha(i) = alpha(ideep)
101  n2_beta(i) = beta(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(sa, ct, p, lat, n2, n2_p, n2_specvol, n2_alpha, n2_beta, dsa, dct, dp)
elemental real(r8) function, public gsw_error_code(err_num, func_name, error_code)