MOM6
MOM_wave_structure.F90
1 !> Vertical structure functions for first baroclinic mode wave speed
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 ! By Benjamin Mater & Robert Hallberg, 2015
7 
8 ! The subroutine in this module calculates the vertical structure
9 ! functions of the first baroclinic mode internal wave speed.
10 ! Calculation of interface values is the same as done in
11 ! MOM_wave_speed by Hallberg, 2008.
12 
13 use mom_debugging, only : isnan => is_nan
14 use mom_diag_mediator, only : post_data, query_averaging_enabled, diag_ctrl
15 use mom_diag_mediator, only : register_diag_field, safe_alloc_ptr, time_type
17 use mom_error_handler, only : mom_error, fatal, warning
19 use mom_grid, only : ocean_grid_type
23 
24 implicit none ; private
25 
26 #include <MOM_memory.h>
27 
28 public wave_structure, wave_structure_init
29 
30 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
31 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
32 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
33 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
34 
35 !> The control structure for the MOM_wave_structure module
36 type, public :: wave_structure_cs ; !private
37  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
38  !! regulate the timing of diagnostic output.
39  real, allocatable, dimension(:,:,:) :: w_strct
40  !< Vertical structure of vertical velocity (normalized) [m s-1].
41  real, allocatable, dimension(:,:,:) :: u_strct
42  !< Vertical structure of horizontal velocity (normalized) [m s-1].
43  real, allocatable, dimension(:,:,:) :: w_profile
44  !< Vertical profile of w_hat(z), where
45  !! w(x,y,z,t) = w_hat(z)*exp(i(kx+ly-freq*t)) is the full time-
46  !! varying vertical velocity with w_hat(z) = W0*w_strct(z) [Z T-1 ~> m s-1].
47  real, allocatable, dimension(:,:,:) :: uavg_profile
48  !< Vertical profile of the magnitude of horizontal velocity,
49  !! (u^2+v^2)^0.5, averaged over a period [L T-1 ~> m s-1].
50  real, allocatable, dimension(:,:,:) :: z_depths
51  !< Depths of layer interfaces [m].
52  real, allocatable, dimension(:,:,:) :: n2
53  !< Squared buoyancy frequency at each interface [s-2].
54  integer, allocatable, dimension(:,:):: num_intfaces
55  !< Number of layer interfaces (including surface and bottom)
56  real :: int_tide_source_x !< X Location of generation site
57  !! for internal tide for testing (BDM)
58  real :: int_tide_source_y !< Y Location of generation site
59  !! for internal tide for testing (BDM)
60 
61 end type wave_structure_cs
62 
63 contains
64 
65 !> This subroutine determines the internal wave velocity structure for any mode.
66 !!
67 !! This subroutine solves for the eigen vector [vertical structure, e(k)] associated with
68 !! the first baroclinic mode speed [i.e., smallest eigen value (lam = 1/c^2)] of the
69 !! system d2e/dz2 = -(N2/cn2)e, or (A-lam*I)e = 0, where A = -(1/N2)(d2/dz2), lam = 1/c^2,
70 !! and I is the identity matrix. 2nd order discretization in the vertical lets this system
71 !! be represented as
72 !!
73 !! -Igu(k)*e(k-1) + (Igu(k)+Igl(k)-lam)*e(k) - Igl(k)*e(k+1) = 0.0
74 !!
75 !! with rigid lid boundary conditions e(1) = e(nz+1) = 0.0 giving
76 !!
77 !! (Igu(2)+Igl(2)-lam)*e(2) - Igl(2)*e(3) = 0.0
78 !! -Igu(nz)*e(nz-1) + (Igu(nz)+Igl(nz)-lam)*e(nz) = 0.0
79 !!
80 !! where, upon noting N2 = reduced gravity/layer thickness, we get
81 !! Igl(k) = 1.0/(gprime(k)*H(k)) ; Igu(k) = 1.0/(gprime(k)*H(k-1))
82 !!
83 !! The eigen value for this system is approximated using "wave_speed." This subroutine uses
84 !! these eigen values (mode speeds) to estimate the corresponding eigen vectors (velocity
85 !! structure) using the "inverse iteration with shift" method. The algorithm is
86 !!
87 !! Pick a starting vector reasonably close to mode structure and with unit magnitude, b_guess
88 !! For n=1,2,3,...
89 !! Solve (A-lam*I)e = e_guess for e
90 !! Set e_guess=e/|e| and repeat, with each iteration refining the estimate of e
91 subroutine wave_structure(h, tv, G, GV, US, cn, ModeNum, freq, CS, En, full_halos)
92  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
93  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
94  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
95  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
96  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
97  !! thermodynamic variables.
98  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: cn !< The (non-rotational) mode internal
99  !! gravity wave speed [L T-1 ~> m s-1].
100  integer, intent(in) :: modenum !< Mode number
101  real, intent(in) :: freq !< Intrinsic wave frequency [T-1 ~> s-1].
102  type(wave_structure_cs), pointer :: cs !< The control structure returned by a
103  !! previous call to wave_structure_init.
104  real, dimension(SZI_(G),SZJ_(G)), &
105  optional, intent(in) :: en !< Internal wave energy density [R Z3 T-2 ~> J m-2]
106  logical, optional, intent(in) :: full_halos !< If true, do the calculation
107  !! over the entire computational domain.
108  ! Local variables
109  real, dimension(SZK_(G)+1) :: &
110  drho_dt, & ! Partial derivative of density with temperature [R degC-1 ~> kg m-3 degC-1]
111  drho_ds, & ! Partial derivative of density with salinity [R ppt-1 ~> kg m-3 ppt-1]
112  pres, & ! Interface pressure [R L2 T-2 ~> Pa]
113  t_int, & ! Temperature interpolated to interfaces [degC]
114  s_int, & ! Salinity interpolated to interfaces [ppt]
115  gprime ! The reduced gravity across each interface [m2 Z-1 s-2 ~> m s-2].
116  real, dimension(SZK_(G)) :: &
117  igl, igu ! The inverse of the reduced gravity across an interface times
118  ! the thickness of the layer below (Igl) or above (Igu) it [s2 m-2].
119  real, dimension(SZK_(G),SZI_(G)) :: &
120  hf, & ! Layer thicknesses after very thin layers are combined [Z ~> m]
121  tf, & ! Layer temperatures after very thin layers are combined [degC]
122  sf, & ! Layer salinities after very thin layers are combined [ppt]
123  rf ! Layer densities after very thin layers are combined [R ~> kg m-3]
124  real, dimension(SZK_(G)) :: &
125  hc, & ! A column of layer thicknesses after convective istabilities are removed [Z ~> m]
126  tc, & ! A column of layer temperatures after convective istabilities are removed [degC]
127  sc, & ! A column of layer salinites after convective istabilities are removed [ppt]
128  rc, & ! A column of layer densities after convective istabilities are removed [R ~> kg m-3]
129  det, ddet
130  real, dimension(SZI_(G),SZJ_(G)) :: &
131  htot ! The vertical sum of the thicknesses [Z ~> m]
132  real :: lam
133  real :: min_h_frac
134  real :: z_to_pres ! A conversion factor from thicknesses to pressure [R L2 T-2 Z-1 ~> Pa m-1]
135  real, dimension(SZI_(G)) :: &
136  hmin, & ! Thicknesses [Z ~> m]
137  h_here, & ! A thickness [Z ~> m]
138  hxt_here, & ! A layer integrated temperature [degC Z ~> degC m]
139  hxs_here, & ! A layer integrated salinity [ppt Z ~> ppt m]
140  hxr_here ! A layer integrated density [R Z ~> kg m-2]
141  real :: speed2_tot
142  real :: i_hnew ! The inverse of a new layer thickness [Z-1 ~> m-1]
143  real :: drxh_sum ! The sum of density diffrences across interfaces times thicknesses [R Z ~> kg m-2]
144  real, parameter :: tol1 = 0.0001, tol2 = 0.001
145  real, pointer, dimension(:,:,:) :: t => null(), s => null()
146  real :: g_rho0 ! G_Earth/Rho0 in [m2 s-2 Z-1 R-1 ~> m4 s-2 kg-1].
147  ! real :: rescale, I_rescale
148  integer :: kf(szi_(g))
149  integer, parameter :: max_itt = 1 ! number of times to iterate in solving for eigenvector
150  real :: cg_subro ! A tiny wave speed to prevent division by zero [L T-1 ~> m s-1]
151  real, parameter :: a_int = 0.5 ! value of normalized integral: \int(w_strct^2)dz = a_int
152  real :: i_a_int ! inverse of a_int
153  real :: f2 ! squared Coriolis frequency [T-2 ~> s-2]
154  real :: kmag2 ! magnitude of horizontal wave number squared
155  logical :: use_eos ! If true, density is calculated from T & S using an
156  ! equation of state.
157  real, dimension(SZK_(G)+1) :: w_strct, u_strct, w_profile, uavg_profile, z_int, n2
158  ! local representations of variables in CS; note,
159  ! not all rows will be filled if layers get merged!
160  real, dimension(SZK_(G)+1) :: w_strct2, u_strct2
161  ! squared values
162  real, dimension(SZK_(G)) :: dz ! thicknesses of merged layers (same as Hc I hope)
163  ! real, dimension(SZK_(G)+1) :: dWdz_profile ! profile of dW/dz
164  real :: w2avg ! average of squared vertical velocity structure funtion
165  real :: int_dwdz2
166  real :: int_w2
167  real :: int_n2w2
168  real :: ke_term ! terms in vertically averaged energy equation
169  real :: pe_term ! terms in vertically averaged energy equation
170  real :: w0 ! A vertical velocity magnitude [Z T-1 ~> m s-1]
171  real :: gp_unscaled ! A version of gprime rescaled to [m s-2].
172  real, dimension(SZK_(G)-1) :: lam_z ! product of eigen value and gprime(k); one value for each
173  ! interface (excluding surface and bottom)
174  real, dimension(SZK_(G)-1) :: a_diag, b_diag, c_diag
175  ! diagonals of tridiagonal matrix; one value for each
176  ! interface (excluding surface and bottom)
177  real, dimension(SZK_(G)-1) :: e_guess ! guess at eigen vector with unit amplitde (for TDMA)
178  real, dimension(SZK_(G)-1) :: e_itt ! improved guess at eigen vector (from TDMA)
179  real :: pi
180  integer :: kc
181  integer :: i, j, k, k2, itt, is, ie, js, je, nz, nzm, row, ig, jg, ig_stop, jg_stop
182 
183  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
184  i_a_int = 1/a_int
185 
186  !if (present(CS)) then
187  if (.not. associated(cs)) call mom_error(fatal, "MOM_wave_structure: "// &
188  "Module must be initialized before it is used.")
189  !endif
190 
191  if (present(full_halos)) then ; if (full_halos) then
192  is = g%isd ; ie = g%ied ; js = g%jsd ; je = g%jed
193  endif ; endif
194 
195  pi = (4.0*atan(1.0))
196 
197  s => tv%S ; t => tv%T
198  g_rho0 = us%L_T_to_m_s**2 * gv%g_Earth / gv%Rho0
199  cg_subro = 1e-100*us%m_s_to_L_T ! The hard-coded value here might need to increase.
200  use_eos = associated(tv%eqn_of_state)
201 
202  ! Simplifying the following could change answers at roundoff.
203  z_to_pres = gv%Z_to_H * (gv%H_to_RZ * gv%g_Earth)
204  ! rescale = 1024.0**4 ; I_rescale = 1.0/rescale
205 
206  min_h_frac = tol1 / real(nz)
207 
208  do j=js,je
209  ! First merge very thin layers with the one above (or below if they are
210  ! at the top). This also transposes the row order so that columns can
211  ! be worked upon one at a time.
212  do i=is,ie ; htot(i,j) = 0.0 ; enddo
213  do k=1,nz ; do i=is,ie ; htot(i,j) = htot(i,j) + h(i,j,k)*gv%H_to_Z ; enddo ; enddo
214 
215  do i=is,ie
216  hmin(i) = htot(i,j)*min_h_frac ; kf(i) = 1 ; h_here(i) = 0.0
217  hxt_here(i) = 0.0 ; hxs_here(i) = 0.0 ; hxr_here(i) = 0.0
218  enddo
219  if (use_eos) then
220  do k=1,nz ; do i=is,ie
221  if ((h_here(i) > hmin(i)) .and. (h(i,j,k)*gv%H_to_Z > hmin(i))) then
222  hf(kf(i),i) = h_here(i)
223  tf(kf(i),i) = hxt_here(i) / h_here(i)
224  sf(kf(i),i) = hxs_here(i) / h_here(i)
225  kf(i) = kf(i) + 1
226 
227  ! Start a new layer
228  h_here(i) = h(i,j,k)*gv%H_to_Z
229  hxt_here(i) = (h(i,j,k)*gv%H_to_Z)*t(i,j,k)
230  hxs_here(i) = (h(i,j,k)*gv%H_to_Z)*s(i,j,k)
231  else
232  h_here(i) = h_here(i) + h(i,j,k)*gv%H_to_Z
233  hxt_here(i) = hxt_here(i) + (h(i,j,k)*gv%H_to_Z)*t(i,j,k)
234  hxs_here(i) = hxs_here(i) + (h(i,j,k)*gv%H_to_Z)*s(i,j,k)
235  endif
236  enddo ; enddo
237  do i=is,ie ; if (h_here(i) > 0.0) then
238  hf(kf(i),i) = h_here(i)
239  tf(kf(i),i) = hxt_here(i) / h_here(i)
240  sf(kf(i),i) = hxs_here(i) / h_here(i)
241  endif ; enddo
242  else
243  do k=1,nz ; do i=is,ie
244  if ((h_here(i) > hmin(i)) .and. (h(i,j,k)*gv%H_to_Z > hmin(i))) then
245  hf(kf(i),i) = h_here(i) ; rf(kf(i),i) = hxr_here(i) / h_here(i)
246  kf(i) = kf(i) + 1
247 
248  ! Start a new layer
249  h_here(i) = h(i,j,k)*gv%H_to_Z
250  hxr_here(i) = (h(i,j,k)*gv%H_to_Z)*gv%Rlay(k)
251  else
252  h_here(i) = h_here(i) + h(i,j,k)*gv%H_to_Z
253  hxr_here(i) = hxr_here(i) + (h(i,j,k)*gv%H_to_Z)*gv%Rlay(k)
254  endif
255  enddo ; enddo
256  do i=is,ie ; if (h_here(i) > 0.0) then
257  hf(kf(i),i) = h_here(i) ; rf(kf(i),i) = hxr_here(i) / h_here(i)
258  endif ; enddo
259  endif ! use_EOS?
260 
261  ! From this point, we can work on individual columns without causing memory
262  ! to have page faults.
263  do i=is,ie ; if (cn(i,j)>0.0)then
264  !----for debugging, remove later----
265  ig = i + g%idg_offset ; jg = j + g%jdg_offset
266  !if (ig == CS%int_tide_source_x .and. jg == CS%int_tide_source_y) then
267  !-----------------------------------
268  if (g%mask2dT(i,j) > 0.5) then
269 
270  lam = 1/(us%L_T_to_m_s**2 * cn(i,j)**2)
271 
272  ! Calculate drxh_sum
273  if (use_eos) then
274  pres(1) = 0.0
275  do k=2,kf(i)
276  pres(k) = pres(k-1) + z_to_pres*hf(k-1,i)
277  t_int(k) = 0.5*(tf(k,i)+tf(k-1,i))
278  s_int(k) = 0.5*(sf(k,i)+sf(k-1,i))
279  enddo
280  call calculate_density_derivs(t_int, s_int, pres, drho_dt, drho_ds, &
281  tv%eqn_of_state, (/2,kf(i)/) )
282 
283  ! Sum the reduced gravities to find out how small a density difference
284  ! is negligibly small.
285  drxh_sum = 0.0
286  do k=2,kf(i)
287  drxh_sum = drxh_sum + 0.5*(hf(k-1,i)+hf(k,i)) * &
288  max(0.0,drho_dt(k)*(tf(k,i)-tf(k-1,i)) + &
289  drho_ds(k)*(sf(k,i)-sf(k-1,i)))
290  enddo
291  else
292  drxh_sum = 0.0
293  do k=2,kf(i)
294  drxh_sum = drxh_sum + 0.5*(hf(k-1,i)+hf(k,i)) * &
295  max(0.0,rf(k,i)-rf(k-1,i))
296  enddo
297  endif ! use_EOS?
298 
299  ! Find gprime across each internal interface, taking care of convective
300  ! instabilities by merging layers.
301  if (drxh_sum >= 0.0) then
302  ! Merge layers to eliminate convective instabilities or exceedingly
303  ! small reduced gravities.
304  if (use_eos) then
305  kc = 1
306  hc(1) = hf(1,i) ; tc(1) = tf(1,i) ; sc(1) = sf(1,i)
307  do k=2,kf(i)
308  if ((drho_dt(k)*(tf(k,i)-tc(kc)) + drho_ds(k)*(sf(k,i)-sc(kc))) * &
309  (hc(kc) + hf(k,i)) < 2.0 * tol2*drxh_sum) then
310  ! Merge this layer with the one above and backtrack.
311  i_hnew = 1.0 / (hc(kc) + hf(k,i))
312  tc(kc) = (hc(kc)*tc(kc) + hf(k,i)*tf(k,i)) * i_hnew
313  sc(kc) = (hc(kc)*sc(kc) + hf(k,i)*sf(k,i)) * i_hnew
314  hc(kc) = (hc(kc) + hf(k,i))
315  ! Backtrack to remove any convective instabilities above... Note
316  ! that the tolerance is a factor of two larger, to avoid limit how
317  ! far back we go.
318  do k2=kc,2,-1
319  if ((drho_dt(k2)*(tc(k2)-tc(k2-1)) + drho_ds(k2)*(sc(k2)-sc(k2-1))) * &
320  (hc(k2) + hc(k2-1)) < tol2*drxh_sum) then
321  ! Merge the two bottommost layers. At this point kc = k2.
322  i_hnew = 1.0 / (hc(kc) + hc(kc-1))
323  tc(kc-1) = (hc(kc)*tc(kc) + hc(kc-1)*tc(kc-1)) * i_hnew
324  sc(kc-1) = (hc(kc)*sc(kc) + hc(kc-1)*sc(kc-1)) * i_hnew
325  hc(kc-1) = (hc(kc) + hc(kc-1))
326  kc = kc - 1
327  else ; exit ; endif
328  enddo
329  else
330  ! Add a new layer to the column.
331  kc = kc + 1
332  drho_ds(kc) = drho_ds(k) ; drho_dt(kc) = drho_dt(k)
333  tc(kc) = tf(k,i) ; sc(kc) = sf(k,i) ; hc(kc) = hf(k,i)
334  endif
335  enddo
336  ! At this point there are kc layers and the gprimes should be positive.
337  do k=2,kc ! Revisit this if non-Boussinesq.
338  gprime(k) = g_rho0 * (drho_dt(k)*(tc(k)-tc(k-1)) + &
339  drho_ds(k)*(sc(k)-sc(k-1)))
340  enddo
341  else ! .not.use_EOS
342  ! Do the same with density directly...
343  kc = 1
344  hc(1) = hf(1,i) ; rc(1) = rf(1,i)
345  do k=2,kf(i)
346  if ((rf(k,i) - rc(kc)) * (hc(kc) + hf(k,i)) < 2.0*tol2*drxh_sum) then
347  ! Merge this layer with the one above and backtrack.
348  rc(kc) = (hc(kc)*rc(kc) + hf(k,i)*rf(k,i)) / (hc(kc) + hf(k,i))
349  hc(kc) = (hc(kc) + hf(k,i))
350  ! Backtrack to remove any convective instabilities above... Note
351  ! that the tolerance is a factor of two larger, to avoid limit how
352  ! far back we go.
353  do k2=kc,2,-1
354  if ((rc(k2)-rc(k2-1)) * (hc(k2)+hc(k2-1)) < tol2*drxh_sum) then
355  ! Merge the two bottommost layers. At this point kc = k2.
356  rc(kc-1) = (hc(kc)*rc(kc) + hc(kc-1)*rc(kc-1)) / (hc(kc) + hc(kc-1))
357  hc(kc-1) = (hc(kc) + hc(kc-1))
358  kc = kc - 1
359  else ; exit ; endif
360  enddo
361  else
362  ! Add a new layer to the column.
363  kc = kc + 1
364  rc(kc) = rf(k,i) ; hc(kc) = hf(k,i)
365  endif
366  enddo
367  ! At this point there are kc layers and the gprimes should be positive.
368  do k=2,kc ! Revisit this if non-Boussinesq.
369  gprime(k) = g_rho0 * (rc(k)-rc(k-1))
370  enddo
371  endif ! use_EOS?
372 
373  !-----------------NOW FIND WAVE STRUCTURE-------------------------------------
374  ! Construct and solve tridiagonal system for the interior interfaces
375  ! Note that kc = number of layers,
376  ! kc+1 = nzm = number of interfaces,
377  ! kc-1 = number of interior interfaces (excluding surface and bottom)
378  ! Also, note that "K" refers to an interface, while "k" refers to the layer below.
379  ! Need at least 3 layers (2 internal interfaces) to generate a matrix, also
380  ! need number of layers to be greater than the mode number
381  if (kc >= max(3, modenum + 1)) then
382  ! Set depth at surface
383  z_int(1) = 0.0
384  ! Calculate Igu, Igl, depth, and N2 at each interior interface
385  ! [excludes surface (K=1) and bottom (K=kc+1)]
386  do k=2,kc
387  igl(k) = 1.0/(gprime(k)*hc(k)) ; igu(k) = 1.0/(gprime(k)*hc(k-1))
388  z_int(k) = z_int(k-1) + hc(k-1)
389  n2(k) = us%m_to_Z**2*gprime(k)/(0.5*(hc(k)+hc(k-1)))
390  enddo
391  ! Set stratification for surface and bottom (setting equal to nearest interface for now)
392  n2(1) = n2(2) ; n2(kc+1) = n2(kc)
393  ! Calcualte depth at bottom
394  z_int(kc+1) = z_int(kc)+hc(kc)
395  ! check that thicknesses sum to total depth
396  if (abs(z_int(kc+1)-htot(i,j)) > 1.e-14*htot(i,j)) then
397  call mom_error(fatal, "wave_structure: mismatch in total depths")
398  endif
399 
400  ! Note that many of the calcluation from here on revert to using vertical
401  ! distances in m, not Z.
402 
403  ! Populate interior rows of tridiagonal matrix; must multiply through by
404  ! gprime to get tridiagonal matrix to the symmetrical form:
405  ! [-1/H(k-1)]e(k-1) + [1/H(k-1)+1/H(k)-lam_z]e(k) + [-1/H(k)]e(k+1) = 0,
406  ! where lam_z = lam*gprime is now a function of depth.
407  ! Frist, populate interior rows
408  do k=3,kc-1
409  row = k-1 ! indexing for TD matrix rows
410  gp_unscaled = us%m_to_Z*gprime(k)
411  lam_z(row) = lam*gp_unscaled
412  a_diag(row) = gp_unscaled*(-igu(k))
413  b_diag(row) = gp_unscaled*(igu(k)+igl(k)) - lam_z(row)
414  c_diag(row) = gp_unscaled*(-igl(k))
415  if (isnan(lam_z(row)))then ; print *, "Wave_structure: lam_z(row) is NAN" ; endif
416  if (isnan(a_diag(row)))then ; print *, "Wave_structure: a(k) is NAN" ; endif
417  if (isnan(b_diag(row)))then ; print *, "Wave_structure: b(k) is NAN" ; endif
418  if (isnan(c_diag(row)))then ; print *, "Wave_structure: c(k) is NAN" ; endif
419  enddo
420  ! Populate top row of tridiagonal matrix
421  k=2 ; row = k-1 ;
422  gp_unscaled = us%m_to_Z*gprime(k)
423  lam_z(row) = lam*gp_unscaled
424  a_diag(row) = 0.0
425  b_diag(row) = gp_unscaled*(igu(k)+igl(k)) - lam_z(row)
426  c_diag(row) = gp_unscaled*(-igl(k))
427  ! Populate bottom row of tridiagonal matrix
428  k=kc ; row = k-1
429  gp_unscaled = us%m_to_Z*gprime(k)
430  lam_z(row) = lam*gp_unscaled
431  a_diag(row) = gp_unscaled*(-igu(k))
432  b_diag(row) = gp_unscaled*(igu(k)+igl(k)) - lam_z(row)
433  c_diag(row) = 0.0
434 
435  ! Guess a vector shape to start with (excludes surface and bottom)
436  e_guess(1:kc-1) = sin((z_int(2:kc)/htot(i,j)) *pi)
437  e_guess(1:kc-1) = e_guess(1:kc-1)/sqrt(sum(e_guess(1:kc-1)**2))
438 
439  ! Perform inverse iteration with tri-diag solver
440  do itt=1,max_itt
441  call tridiag_solver(a_diag(1:kc-1),b_diag(1:kc-1),c_diag(1:kc-1), &
442  -lam_z(1:kc-1),e_guess(1:kc-1),"TDMA_H",e_itt)
443  e_guess(1:kc-1) = e_itt(1:kc-1) / sqrt(sum(e_itt(1:kc-1)**2))
444  enddo ! itt-loop
445  w_strct(2:kc) = e_guess(1:kc-1)
446  w_strct(1) = 0.0 ! rigid lid at surface
447  w_strct(kc+1) = 0.0 ! zero-flux at bottom
448 
449  ! Check to see if solver worked
450  ig_stop = 0 ; jg_stop = 0
451  if (isnan(sum(w_strct(1:kc+1))))then
452  print *, "Wave_structure: w_strct has a NAN at ig=", ig, ", jg=", jg
453  if (i<g%isc .or. i>g%iec .or. j<g%jsc .or. j>g%jec)then
454  print *, "This is occuring at a halo point."
455  endif
456  ig_stop = ig ; jg_stop = jg
457  endif
458 
459  ! Normalize vertical structure function of w such that
460  ! \int(w_strct)^2dz = a_int (a_int could be any value, e.g., 0.5)
461  nzm = kc+1 ! number of layer interfaces after merging
462  !(including surface and bottom)
463  w2avg = 0.0
464  do k=1,nzm-1
465  dz(k) = us%Z_to_m*hc(k)
466  w2avg = w2avg + 0.5*(w_strct(k)**2+w_strct(k+1)**2)*dz(k)
467  enddo
468  !### Some mathematical cancellations could occur in the next two lines.
469  w2avg = w2avg / htot(i,j)
470  w_strct(:) = w_strct(:) / sqrt(htot(i,j)*w2avg*i_a_int)
471 
472  ! Calculate vertical structure function of u (i.e. dw/dz)
473  do k=2,nzm-1
474  u_strct(k) = 0.5*((w_strct(k-1) - w_strct(k) )/dz(k-1) + &
475  (w_strct(k) - w_strct(k+1))/dz(k))
476  enddo
477  u_strct(1) = (w_strct(1) - w_strct(2) )/dz(1)
478  u_strct(nzm) = (w_strct(nzm-1)- w_strct(nzm))/dz(nzm-1)
479 
480  ! Calculate wavenumber magnitude
481  f2 = g%CoriolisBu(i,j)**2
482  !f2 = 0.25*US%s_to_T**2 *((G%CoriolisBu(I,J)**2 + G%CoriolisBu(I-1,J-1)**2) + &
483  ! (G%CoriolisBu(I,J-1)**2 + G%CoriolisBu(I-1,J)**2))
484  kmag2 = us%m_to_L**2 * (freq**2 - f2) / (cn(i,j)**2 + cg_subro**2)
485 
486  ! Calculate terms in vertically integrated energy equation
487  int_dwdz2 = 0.0 ; int_w2 = 0.0 ; int_n2w2 = 0.0
488  u_strct2(1:nzm) = u_strct(1:nzm)**2
489  w_strct2(1:nzm) = w_strct(1:nzm)**2
490  ! vertical integration with Trapezoidal rule
491  do k=1,nzm-1
492  int_dwdz2 = int_dwdz2 + 0.5*(u_strct2(k)+u_strct2(k+1)) * us%m_to_Z*dz(k)
493  int_w2 = int_w2 + 0.5*(w_strct2(k)+w_strct2(k+1)) * us%m_to_Z*dz(k)
494  int_n2w2 = int_n2w2 + 0.5*(w_strct2(k)*n2(k)+w_strct2(k+1)*n2(k+1)) * us%m_to_Z*dz(k)
495  enddo
496 
497  ! Back-calculate amplitude from energy equation
498  if (present(en) .and. (freq**2*kmag2 > 0.0)) then
499  ! Units here are [R
500  ke_term = 0.25*gv%Rho0*( ((freq**2 + f2) / (freq**2*kmag2))*int_dwdz2 + int_w2 )
501  pe_term = 0.25*gv%Rho0*( int_n2w2 / (us%s_to_T*freq)**2 )
502  if (en(i,j) >= 0.0) then
503  w0 = sqrt( en(i,j) / (ke_term + pe_term) )
504  else
505  call mom_error(warning, "wave_structure: En < 0.0; setting to W0 to 0.0")
506  print *, "En(i,j)=", en(i,j), " at ig=", ig, ", jg=", jg
507  w0 = 0.0
508  endif
509  ! Calculate actual vertical velocity profile and derivative
510  w_profile(:) = w0*w_strct(:)
511  ! dWdz_profile(:) = W0*u_strct(:)
512  ! Calculate average magnitude of actual horizontal velocity over a period
513  uavg_profile(:) = us%Z_to_L*abs(w0*u_strct(:)) * sqrt((freq**2 + f2) / (2.0*freq**2*kmag2))
514  else
515  w_profile(:) = 0.0
516  ! dWdz_profile(:) = 0.0
517  uavg_profile(:) = 0.0
518  endif
519 
520  ! Store values in control structure
521  cs%w_strct(i,j,1:nzm) = w_strct(1:nzm)
522  cs%u_strct(i,j,1:nzm) = u_strct(1:nzm)
523  cs%W_profile(i,j,1:nzm) = w_profile(1:nzm)
524  cs%Uavg_profile(i,j,1:nzm)= uavg_profile(1:nzm)
525  cs%z_depths(i,j,1:nzm) = us%Z_to_m*z_int(1:nzm)
526  cs%N2(i,j,1:nzm) = n2(1:nzm)
527  cs%num_intfaces(i,j) = nzm
528  else
529  ! If not enough layers, default to zero
530  nzm = kc+1
531  cs%w_strct(i,j,1:nzm) = 0.0
532  cs%u_strct(i,j,1:nzm) = 0.0
533  cs%W_profile(i,j,1:nzm) = 0.0
534  cs%Uavg_profile(i,j,1:nzm)= 0.0
535  cs%z_depths(i,j,1:nzm) = 0.0 ! could use actual values
536  cs%N2(i,j,1:nzm) = 0.0 ! could use with actual values
537  cs%num_intfaces(i,j) = nzm
538  endif ! kc >= 3 and kc > ModeNum + 1?
539  endif ! drxh_sum >= 0?
540  !else ! if at test point - delete later
541  ! return ! if at test point - delete later
542  !endif ! if at test point - delete later
543  endif ! mask2dT > 0.5?
544  else
545  ! if cn=0.0, default to zero
546  nzm = nz+1! could use actual values
547  cs%w_strct(i,j,1:nzm) = 0.0
548  cs%u_strct(i,j,1:nzm) = 0.0
549  cs%W_profile(i,j,1:nzm) = 0.0
550  cs%Uavg_profile(i,j,1:nzm)= 0.0
551  cs%z_depths(i,j,1:nzm) = 0.0 ! could use actual values
552  cs%N2(i,j,1:nzm) = 0.0 ! could use with actual values
553  cs%num_intfaces(i,j) = nzm
554  endif ; enddo ! if cn>0.0? ; i-loop
555  enddo ! j-loop
556 
557 end subroutine wave_structure
558 
559 !> Solves a tri-diagonal system Ax=y using either the standard
560 !! Thomas algorithm (TDMA_T) or its more stable variant that invokes the
561 !! "Hallberg substitution" (TDMA_H).
562 subroutine tridiag_solver(a, b, c, h, y, method, x)
563  real, dimension(:), intent(in) :: a !< lower diagonal with first entry equal to zero.
564  real, dimension(:), intent(in) :: b !< middle diagonal.
565  real, dimension(:), intent(in) :: c !< upper diagonal with last entry equal to zero.
566  real, dimension(:), intent(in) :: h !< vector of values that have already been added to b; used
567  !! for systems of the form (e.g. average layer thickness in vertical diffusion case):
568  !! [ -alpha(k-1/2) ] * e(k-1) +
569  !! [ alpha(k-1/2) + alpha(k+1/2) + h(k) ] * e(k) +
570  !! [ -alpha(k+1/2) ] * e(k+1) = y(k)
571  !! where a(k)=[-alpha(k-1/2)], b(k)=[alpha(k-1/2)+alpha(k+1/2) + h(k)],
572  !! and c(k)=[-alpha(k+1/2)]. Only used with TDMA_H method.
573  real, dimension(:), intent(in) :: y !< vector of known values on right hand side.
574  character(len=*), intent(in) :: method !< A string describing the algorithm to use
575  real, dimension(:), intent(out) :: x !< vector of unknown values to solve for.
576  ! Local variables
577  integer :: nrow ! number of rows in A matrix
578 ! real, allocatable, dimension(:,:) :: A_check ! for solution checking
579 ! real, allocatable, dimension(:) :: y_check ! for solution checking
580  real, allocatable, dimension(:) :: c_prime, y_prime, q, alpha
581  ! intermediate values for solvers
582  real :: Q_prime, beta ! intermediate values for solver
583  integer :: k ! row (e.g. interface) index
584  integer :: i,j
585 
586  nrow = size(y)
587  allocate(c_prime(nrow))
588  allocate(y_prime(nrow))
589  allocate(q(nrow))
590  allocate(alpha(nrow))
591 ! allocate(A_check(nrow,nrow))
592 ! allocate(y_check(nrow))
593 
594  if (method == 'TDMA_T') then
595  ! Standard Thomas algoritim (4th variant).
596  ! Note: Requires A to be non-singular for accuracy/stability
597  c_prime(:) = 0.0 ; y_prime(:) = 0.0
598  c_prime(1) = c(1)/b(1) ; y_prime(1) = y(1)/b(1)
599 
600  ! Forward sweep
601  do k=2,nrow-1
602  c_prime(k) = c(k)/(b(k)-a(k)*c_prime(k-1))
603  enddo
604  !print *, 'c_prime=', c_prime(1:nrow)
605  do k=2,nrow
606  y_prime(k) = (y(k)-a(k)*y_prime(k-1))/(b(k)-a(k)*c_prime(k-1))
607  enddo
608  !print *, 'y_prime=', y_prime(1:nrow)
609  x(nrow) = y_prime(nrow)
610 
611  ! Backward sweep
612  do k=nrow-1,1,-1
613  x(k) = y_prime(k)-c_prime(k)*x(k+1)
614  enddo
615  !print *, 'x=',x(1:nrow)
616 
617  ! Check results - delete later
618  !do j=1,nrow ; do i=1,nrow
619  ! if (i==j)then ; A_check(i,j) = b(i)
620  ! elseif (i==j+1)then ; A_check(i,j) = a(i)
621  ! elseif (i==j-1)then ; A_check(i,j) = c(i)
622  ! endif
623  !enddo ; enddo
624  !print *, 'A(2,1),A(2,2),A(1,2)=', A_check(2,1), A_check(2,2), A_check(1,2)
625  !y_check = matmul(A_check,x)
626  !if (all(y_check /= y))then
627  ! print *, "tridiag_solver: Uh oh, something's not right!"
628  ! print *, "y=", y
629  ! print *, "y_check=", y_check
630  !endif
631 
632  elseif (method == 'TDMA_H') then
633  ! Thomas algoritim (4th variant) w/ Hallberg substitution.
634  ! For a layered system where k is at interfaces, alpha{k+1/2} refers to
635  ! some property (e.g. inverse thickness for mode-structure problem) of the
636  ! layer below and alpha{k-1/2} refers to the layer above.
637  ! Here, alpha(k)=alpha{k+1/2} and alpha(k-1)=alpha{k-1/2}.
638  ! Strictly speaking, this formulation requires A to be a non-singular,
639  ! symmetric, diagonally dominant matrix, with h>0.
640  ! Need to add a check for these conditions.
641  do k=1,nrow-1
642  if (abs(a(k+1)-c(k)) > 1.e-10*(abs(a(k+1))+abs(c(k)))) then
643  call mom_error(warning, "tridiag_solver: matrix not symmetric; need symmetry when invoking TDMA_H")
644  endif
645  enddo
646  alpha = -c
647  ! Alpha of the bottom-most layer is not necessarily zero. Therefore,
648  ! back out the value from the provided b(nrow and h(nrow) values
649  alpha(nrow) = b(nrow)-h(nrow)-alpha(nrow-1)
650  ! Prime other variables
651  beta = 1/b(1)
652  y_prime(:) = 0.0 ; q(:) = 0.0
653  y_prime(1) = beta*y(1) ; q(1) = beta*alpha(1)
654  q_prime = 1-q(1)
655 
656  ! Forward sweep
657  do k=2,nrow-1
658  beta = 1/(h(k)+alpha(k-1)*q_prime+alpha(k))
659  if (isnan(beta))then ; print *, "Tridiag_solver: beta is NAN" ; endif
660  q(k) = beta*alpha(k)
661  y_prime(k) = beta*(y(k)+alpha(k-1)*y_prime(k-1))
662  q_prime = beta*(h(k)+alpha(k-1)*q_prime)
663  enddo
664  if ((h(nrow)+alpha(nrow-1)*q_prime+alpha(nrow)) == 0.0)then
665  call mom_error(fatal, "Tridiag_solver: this system is not stable.") ! ; overriding beta(nrow)
666  ! This has hard-coded dimensions: beta = 1/(1e-15) ! place holder for unstable systems - delete later
667  else
668  beta = 1/(h(nrow)+alpha(nrow-1)*q_prime+alpha(nrow))
669  endif
670  y_prime(nrow) = beta*(y(nrow)+alpha(nrow-1)*y_prime(nrow-1))
671  x(nrow) = y_prime(nrow)
672  ! Backward sweep
673  do k=nrow-1,1,-1
674  x(k) = y_prime(k)+q(k)*x(k+1)
675  enddo
676  !print *, 'yprime=',y_prime(1:nrow)
677  !print *, 'x=',x(1:nrow)
678  endif
679 
680  deallocate(c_prime,y_prime,q,alpha)
681 ! deallocate(A_check,y_check)
682 
683 end subroutine tridiag_solver
684 
685 !> Allocate memory associated with the wave structure module and read parameters.
686 subroutine wave_structure_init(Time, G, param_file, diag, CS)
687  type(time_type), intent(in) :: time !< The current model time.
688  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
689  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
690  !! parameters.
691  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
692  !! diagnostic output.
693  type(wave_structure_cs), pointer :: cs !< A pointer that is set to point to the
694  !! control structure for this module.
695 ! This include declares and sets the variable "version".
696 #include "version_variable.h"
697  character(len=40) :: mdl = "MOM_wave_structure" ! This module's name.
698  integer :: isd, ied, jsd, jed, nz
699 
700  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
701 
702  if (associated(cs)) then
703  call mom_error(warning, "wave_structure_init called with an "// &
704  "associated control structure.")
705  return
706  else ; allocate(cs) ; endif
707 
708  call get_param(param_file, mdl, "INTERNAL_TIDE_SOURCE_X", cs%int_tide_source_x, &
709  "X Location of generation site for internal tide", default=1.)
710  call get_param(param_file, mdl, "INTERNAL_TIDE_SOURCE_Y", cs%int_tide_source_y, &
711  "Y Location of generation site for internal tide", default=1.)
712 
713  cs%diag => diag
714 
715  ! Allocate memory for variable in control structure; note,
716  ! not all rows will be filled if layers get merged!
717  allocate(cs%w_strct(isd:ied,jsd:jed,nz+1))
718  allocate(cs%u_strct(isd:ied,jsd:jed,nz+1))
719  allocate(cs%W_profile(isd:ied,jsd:jed,nz+1))
720  allocate(cs%Uavg_profile(isd:ied,jsd:jed,nz+1))
721  allocate(cs%z_depths(isd:ied,jsd:jed,nz+1))
722  allocate(cs%N2(isd:ied,jsd:jed,nz+1))
723  allocate(cs%num_intfaces(isd:ied,jsd:jed))
724 
725  ! Write all relevant parameters to the model log.
726  call log_version(param_file, mdl, version, "")
727 
728 end subroutine wave_structure_init
729 
730 end module mom_wave_structure
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:80
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_wave_structure
Vertical structure functions for first baroclinic mode wave speed.
Definition: MOM_wave_structure.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_wave_structure::wave_structure_cs
The control structure for the MOM_wave_structure module.
Definition: MOM_wave_structure.F90:36
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_eos::calculate_density_derivs
Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition: MOM_EOS.F90:81
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:240