MOM6
Kelvin_initialization.F90
1 !> Configures the model for the Kelvin wave experiment.
2 !!
3 !! Kelvin = coastally-trapped Kelvin waves from the ROMS examples.
4 !! Initialize with level surfaces and drive the wave in at the west,
5 !! radiate out at the east.
7 
8 ! This file is part of MOM6. See LICENSE.md for the license.
9 
11 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
13 use mom_grid, only : ocean_grid_type
14 use mom_open_boundary, only : ocean_obc_type, obc_none
15 use mom_open_boundary, only : obc_segment_type, register_obc
16 use mom_open_boundary, only : obc_direction_n, obc_direction_e
17 use mom_open_boundary, only : obc_direction_s, obc_direction_w
21 use mom_time_manager, only : time_type, time_type_to_real
22 
23 implicit none ; private
24 
25 #include <MOM_memory.h>
26 
27 public kelvin_set_obc_data, kelvin_initialize_topography
28 public register_kelvin_obc, kelvin_obc_end
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 !> Control structure for Kelvin wave open boundaries.
36 type, public :: kelvin_obc_cs ; private
37  integer :: mode = 0 !< Vertical mode
38  real :: coast_angle = 0 !< Angle of coastline
39  real :: coast_offset1 = 0 !< Longshore distance to coastal angle
40  real :: coast_offset2 = 0 !< Longshore distance to coastal angle
41  real :: h0 = 0 !< Bottom depth
42  real :: f_0 !< Coriolis parameter
43  real :: rho_range !< Density range
44  real :: rho_0 !< Mean density
45 end type kelvin_obc_cs
46 
47 ! This include declares and sets the variable "version".
48 #include "version_variable.h"
49 
50 contains
51 
52 !> Add Kelvin wave to OBC registry.
53 function register_kelvin_obc(param_file, CS, OBC_Reg)
54  type(param_file_type), intent(in) :: param_file !< parameter file.
55  type(kelvin_obc_cs), pointer :: cs !< Kelvin wave control structure.
56  type(obc_registry_type), pointer :: obc_reg !< OBC registry.
57 
58  ! Local variables
59  logical :: register_kelvin_obc
60  character(len=40) :: mdl = "register_Kelvin_OBC" !< This subroutine's name.
61  character(len=32) :: casename = "Kelvin wave" !< This case's name.
62  character(len=200) :: config
63 
64  if (associated(cs)) then
65  call mom_error(warning, "register_Kelvin_OBC called with an "// &
66  "associated control structure.")
67  return
68  endif
69  allocate(cs)
70 
71  call log_version(param_file, mdl, version, "")
72  call get_param(param_file, mdl, "KELVIN_WAVE_MODE", cs%mode, &
73  "Vertical Kelvin wave mode imposed at upstream open boundary.", &
74  default=0)
75  call get_param(param_file, mdl, "F_0", cs%F_0, &
76  default=0.0, do_not_log=.true.)
77  call get_param(param_file, mdl, "TOPO_CONFIG", config, do_not_log=.true.)
78  if (trim(config) == "Kelvin") then
79  call get_param(param_file, mdl, "ROTATED_COAST_OFFSET_1", cs%coast_offset1, &
80  "The distance along the southern and northern boundaries "//&
81  "at which the coasts angle in.", &
82  units="km", default=100.0)
83  call get_param(param_file, mdl, "ROTATED_COAST_OFFSET_2", cs%coast_offset2, &
84  "The distance from the southern and northern boundaries "//&
85  "at which the coasts angle in.", &
86  units="km", default=10.0)
87  call get_param(param_file, mdl, "ROTATED_COAST_ANGLE", cs%coast_angle, &
88  "The angle of the southern bondary beyond X=ROTATED_COAST_OFFSET.", &
89  units="degrees", default=11.3)
90  cs%coast_angle = cs%coast_angle * (atan(1.0)/45.) ! Convert to radians
91  cs%coast_offset1 = cs%coast_offset1 * 1.e3 ! Convert to m
92  cs%coast_offset2 = cs%coast_offset2 * 1.e3 ! Convert to m
93  endif
94  if (cs%mode /= 0) then
95  call get_param(param_file, mdl, "DENSITY_RANGE", cs%rho_range, &
96  default=2.0, do_not_log=.true.)
97  call get_param(param_file, mdl, "RHO_0", cs%rho_0, &
98  default=1035.0, do_not_log=.true.)
99  call get_param(param_file, mdl, "MAXIMUM_DEPTH", cs%H0, &
100  default=1000.0, do_not_log=.true.)
101  endif
102 
103  ! Register the Kelvin open boundary.
104  call register_obc(casename, param_file, obc_reg)
105  register_kelvin_obc = .true.
106 
107 end function register_kelvin_obc
108 
109 !> Clean up the Kelvin wave OBC from registry.
110 subroutine kelvin_obc_end(CS)
111  type(kelvin_obc_cs), pointer :: cs !< Kelvin wave control structure.
112 
113  if (associated(cs)) then
114  deallocate(cs)
115  endif
116 end subroutine kelvin_obc_end
117 
118 ! -----------------------------------------------------------------------------
119 !> This subroutine sets up the Kelvin topography and land mask
120 subroutine kelvin_initialize_topography(D, G, param_file, max_depth, US)
121  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
122  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
123  intent(out) :: d !< Ocean bottom depth in m or Z if US is present
124  type(param_file_type), intent(in) :: param_file !< Parameter file structure
125  real, intent(in) :: max_depth !< Maximum model depth in the units of D
126  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
127 
128  ! Local variables
129  character(len=40) :: mdl = "Kelvin_initialize_topography" ! This subroutine's name.
130  real :: m_to_z ! A dimensional rescaling factor.
131  real :: min_depth ! The minimum and maximum depths [Z ~> m].
132  real :: pi ! 3.1415...
133  real :: coast_offset1, coast_offset2, coast_angle, right_angle
134  integer :: i, j
135 
136  call mom_mesg(" Kelvin_initialization.F90, Kelvin_initialize_topography: setting topography", 5)
137 
138  m_to_z = 1.0 ; if (present(us)) m_to_z = us%m_to_Z
139 
140  call get_param(param_file, mdl, "MINIMUM_DEPTH", min_depth, &
141  "The minimum depth of the ocean.", units="m", default=0.0, scale=m_to_z)
142  call get_param(param_file, mdl, "ROTATED_COAST_OFFSET_1", coast_offset1, &
143  default=100.0, do_not_log=.true.)
144  call get_param(param_file, mdl, "ROTATED_COAST_OFFSET_2", coast_offset2, &
145  default=10.0, do_not_log=.true.)
146  call get_param(param_file, mdl, "ROTATED_COAST_ANGLE", coast_angle, &
147  default=11.3, do_not_log=.true.)
148 
149  coast_angle = coast_angle * (atan(1.0)/45.) ! Convert to radians
150  right_angle = 2 * atan(1.0)
151 
152  do j=g%jsc,g%jec ; do i=g%isc,g%iec
153  d(i,j) = max_depth
154  ! Southern side
155  if ((g%geoLonT(i,j) - g%west_lon > coast_offset1) .AND. &
156  (atan2(g%geoLatT(i,j) - g%south_lat + coast_offset2, &
157  g%geoLonT(i,j) - g%west_lon - coast_offset1) < coast_angle)) &
158  d(i,j) = 0.5*min_depth
159  ! Northern side
160  if ((g%geoLonT(i,j) - g%west_lon < g%len_lon - coast_offset1) .AND. &
161  (atan2(g%len_lat + g%south_lat + coast_offset2 - g%geoLatT(i,j), &
162  g%len_lon + g%west_lon - coast_offset1 - g%geoLonT(i,j)) < coast_angle)) &
163  d(i,j) = 0.5*min_depth
164 
165  if (d(i,j) > max_depth) d(i,j) = max_depth
166  if (d(i,j) < min_depth) d(i,j) = 0.5*min_depth
167  enddo ; enddo
168 
169 end subroutine kelvin_initialize_topography
170 
171 !> This subroutine sets the properties of flow at open boundary conditions.
172 subroutine kelvin_set_obc_data(OBC, CS, G, GV, US, h, Time)
173  type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies
174  !! whether, where, and what open boundary
175  !! conditions are used.
176  type(kelvin_obc_cs), pointer :: cs !< Kelvin wave control structure.
177  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
178  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
179  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
180  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< layer thickness [H ~> m or kg m-2].
181  type(time_type), intent(in) :: time !< model time.
182 
183  ! The following variables are used to set up the transport in the Kelvin example.
184  real :: time_sec, cff
185  real :: n0 ! Brunt-Vaisala frequency [s-1]
186  real :: plx !< Longshore wave parameter
187  real :: pmz !< Vertical wave parameter
188  real :: lambda !< Offshore decay scale
189  real :: omega !< Wave frequency [s-1]
190  real :: pi
191  integer :: i, j, k, n, is, ie, js, je, isd, ied, jsd, jed, nz
192  integer :: isdb, iedb, jsdb, jedb
193  real :: fac, x, y, x1, y1
194  real :: val1, val2, sina, cosa
195  type(obc_segment_type), pointer :: segment => null()
196 
197  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
198  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
199  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
200 
201  if (.not.associated(obc)) call mom_error(fatal, 'Kelvin_initialization.F90: '// &
202  'Kelvin_set_OBC_data() was called but OBC type was not initialized!')
203 
204  time_sec = time_type_to_real(time)
205  pi = 4.0*atan(1.0)
206  fac = 1.0
207 
208  if (cs%mode == 0) then
209  omega = 2.0 * pi / (12.42 * 3600.0) ! M2 Tide period
210  val1 = us%m_to_Z * sin(omega * time_sec)
211  else
212  n0 = us%L_to_m*us%s_to_T * sqrt((cs%rho_range / cs%rho_0) * gv%g_Earth * (us%m_to_Z * cs%H0))
213  ! Two wavelengths in domain
214  plx = 4.0 * pi / g%len_lon
215  pmz = pi * cs%mode / cs%H0
216  lambda = pmz * cs%F_0 / n0
217  omega = cs%F_0 * plx / lambda
218 
219  ! lambda = PI * CS%mode * CS%F_0 / (CS%H0 * N0)
220  ! omega = (4.0 * CS%H0 * N0) / (CS%mode * G%len_lon)
221  endif
222 
223  sina = sin(cs%coast_angle)
224  cosa = cos(cs%coast_angle)
225  do n=1,obc%number_of_segments
226  segment => obc%segment(n)
227  if (.not. segment%on_pe) cycle
228  ! Apply values to the inflow end only.
229  if (segment%direction == obc_direction_e) cycle
230  if (segment%direction == obc_direction_n) cycle
231 
232  ! This should be somewhere else...
233  segment%Velocity_nudging_timescale_in = 1.0/(0.3*86400)
234 
235  if (segment%direction == obc_direction_w) then
236  isdb = segment%HI%IsdB ; iedb = segment%HI%IedB
237  jsd = segment%HI%jsd ; jed = segment%HI%jed
238  jsdb = segment%HI%JsdB ; jedb = segment%HI%JedB
239  do j=jsd,jed ; do i=isdb,iedb
240  x1 = 1000. * g%geoLonCu(i,j)
241  y1 = 1000. * g%geoLatCu(i,j)
242  x = (x1 - cs%coast_offset1) * cosa + y1 * sina
243  y = - (x1 - cs%coast_offset1) * sina + y1 * cosa
244  if (cs%mode == 0) then
245  ! Use inside bathymetry
246  cff = sqrt(gv%g_Earth * g%bathyT(i+1,j) )
247  val2 = fac * exp(- us%T_to_s*cs%F_0 * us%m_to_L*y / cff)
248  segment%eta(i,j) = val2 * cos(omega * time_sec)
249  segment%normal_vel_bt(i,j) = (val2 * (val1 * cff * cosa / &
250  (g%bathyT(i+1,j) )) )
251  if (segment%nudged) then
252  do k=1,nz
253  segment%nudged_normal_vel(i,j,k) = (val2 * (val1 * cff * cosa / &
254  (g%bathyT(i+1,j))) )
255  enddo
256  elseif (segment%specified) then
257  do k=1,nz
258  segment%normal_vel(i,j,k) = (val2 * (val1 * cff * cosa / &
259  (g%bathyT(i+1,j) )) )
260  segment%normal_trans(i,j,k) = segment%normal_vel(i,j,k) * h(i+1,j,k) * g%dyCu(i,j)
261  enddo
262  endif
263  else
264  ! Not rotated yet
265  segment%eta(i,j) = 0.0
266  segment%normal_vel_bt(i,j) = 0.0
267  if (segment%nudged) then
268  do k=1,nz
269  segment%nudged_normal_vel(i,j,k) = us%m_s_to_L_T * fac * lambda / cs%F_0 * &
270  exp(- lambda * y) * cos(pi * cs%mode * (k - 0.5) / nz) * &
271  cos(omega * time_sec)
272  enddo
273  elseif (segment%specified) then
274  do k=1,nz
275  segment%normal_vel(i,j,k) = us%m_s_to_L_T * fac * lambda / cs%F_0 * &
276  exp(- lambda * y) * cos(pi * cs%mode * (k - 0.5) / nz) * &
277  cos(omega * time_sec)
278  segment%normal_trans(i,j,k) = segment%normal_vel(i,j,k) * h(i+1,j,k) * g%dyCu(i,j)
279  enddo
280  endif
281  endif
282  enddo ; enddo
283  if (associated(segment%tangential_vel)) then
284  do j=jsdb+1,jedb-1 ; do i=isdb,iedb
285  x1 = 1000. * g%geoLonBu(i,j)
286  y1 = 1000. * g%geoLatBu(i,j)
287  x = (x1 - cs%coast_offset1) * cosa + y1 * sina
288  y = - (x1 - cs%coast_offset1) * sina + y1 * cosa
289  cff =sqrt(gv%g_Earth * g%bathyT(i+1,j) )
290  val2 = fac * exp(- us%T_to_s*cs%F_0 * us%m_to_L*y / cff)
291  if (cs%mode == 0) then ; do k=1,nz
292  segment%tangential_vel(i,j,k) = (val1 * val2 * cff * sina) / &
293  ( 0.5*(g%bathyT(i+1,j+1) + g%bathyT(i+1,j) ) )
294 
295  enddo ; endif
296  enddo ; enddo
297  endif
298  else ! Must be south
299  isd = segment%HI%isd ; ied = segment%HI%ied
300  jsdb = segment%HI%JsdB ; jedb = segment%HI%JedB
301  do j=jsdb,jedb ; do i=isd,ied
302  x1 = 1000. * g%geoLonCv(i,j)
303  y1 = 1000. * g%geoLatCv(i,j)
304  x = (x1 - cs%coast_offset1) * cosa + y1 * sina
305  y = - (x1 - cs%coast_offset1) * sina + y1 * cosa
306  if (cs%mode == 0) then
307  cff = sqrt(gv%g_Earth * g%bathyT(i,j+1) )
308  val2 = fac * exp(- 0.5 * (g%CoriolisBu(i,j) + g%CoriolisBu(i-1,j)) * us%m_to_L*y / cff)
309  segment%eta(i,j) = val2 * cos(omega * time_sec)
310  segment%normal_vel_bt(i,j) = us%L_T_to_m_s * (val1 * cff * sina / &
311  (g%bathyT(i,j+1) )) * val2
312  if (segment%nudged) then
313  do k=1,nz
314  segment%nudged_normal_vel(i,j,k) = us%L_T_to_m_s * (val1 * cff * sina / &
315  (g%bathyT(i,j+1) )) * val2
316  enddo
317  elseif (segment%specified) then
318  do k=1,nz
319  segment%normal_vel(i,j,k) = us%L_T_to_m_s * (val1 * cff * sina / &
320  (g%bathyT(i,j+1) )) * val2
321  segment%normal_trans(i,j,k) = segment%normal_vel(i,j,k) * h(i,j+1,k) * g%dxCv(i,j)
322  enddo
323  endif
324  else
325  ! Not rotated yet
326  segment%eta(i,j) = 0.0
327  segment%normal_vel_bt(i,j) = 0.0
328  if (segment%nudged) then
329  do k=1,nz
330  segment%nudged_normal_vel(i,j,k) = us%m_s_to_L_T*fac * lambda / cs%F_0 * &
331  exp(- lambda * y) * cos(pi * cs%mode * (k - 0.5) / nz) * cosa
332  enddo
333  elseif (segment%specified) then
334  do k=1,nz
335  segment%normal_vel(i,j,k) = us%m_s_to_L_T*fac * lambda / cs%F_0 * &
336  exp(- lambda * y) * cos(pi * cs%mode * (k - 0.5) / nz) * cosa
337  segment%normal_trans(i,j,k) = segment%normal_vel(i,j,k) * h(i,j+1,k) * g%dxCv(i,j)
338  enddo
339  endif
340  endif
341  enddo ; enddo
342  if (associated(segment%tangential_vel)) then
343  do j=jsdb,jedb ; do i=isdb+1,iedb-1
344  x1 = 1000. * g%geoLonBu(i,j)
345  y1 = 1000. * g%geoLatBu(i,j)
346  x = (x1 - cs%coast_offset1) * cosa + y1 * sina
347  y = - (x1 - cs%coast_offset1) * sina + y1 * cosa
348  cff = sqrt(gv%g_Earth * g%bathyT(i,j+1) )
349  val2 = fac * exp(- 0.5 * (g%CoriolisBu(i,j) + g%CoriolisBu(i-1,j)) * us%m_to_L*y / cff)
350  if (cs%mode == 0) then ; do k=1,nz
351  segment%tangential_vel(i,j,k) = ((val1 * val2 * cff * sina) / &
352  ( 0.5*((g%bathyT(i+1,j+1)) + g%bathyT(i,j+1))) )
353  enddo ; endif
354  enddo ; enddo
355  endif
356  endif
357  enddo
358 
359 end subroutine kelvin_set_obc_data
360 
361 end module kelvin_initialization
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
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_dyn_horgrid
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Definition: MOM_dyn_horgrid.F90:3
kelvin_initialization
Configures the model for the Kelvin wave experiment.
Definition: Kelvin_initialization.F90:6
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
kelvin_initialization::kelvin_obc_cs
Control structure for Kelvin wave open boundaries.
Definition: Kelvin_initialization.F90:36
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_open_boundary::obc_registry_type
Type to carry basic OBC information needed for updating values.
Definition: MOM_open_boundary.F90:337
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_open_boundary
Controls where open boundary conditions are applied.
Definition: MOM_open_boundary.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_open_boundary::ocean_obc_type
Open-boundary data.
Definition: MOM_open_boundary.F90:218
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_open_boundary::obc_segment_type
Open boundary segment data structure.
Definition: MOM_open_boundary.F90:113
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_dyn_horgrid::dyn_horgrid_type
Describes the horizontal ocean grid with only dynamic memory arrays.
Definition: MOM_dyn_horgrid.F90:23
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26