MOM6
shelfwave_initialization.F90
1 !> Configures the model for the idealized shelfwave test case.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_domains, only : sum_across_pes
8 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
10 use mom_grid, only : ocean_grid_type
11 use mom_open_boundary, only : ocean_obc_type, obc_none, obc_direction_w
12 use mom_open_boundary, only : obc_segment_type, register_obc
14 use mom_time_manager, only : time_type, time_type_to_real
16 
17 implicit none ; private
18 
19 #include <MOM_memory.h>
20 
21 character(len=40) :: mdl = "shelfwave_initialization" !< This module's name.
22 
23 ! The following routines are visible to the outside world
24 public shelfwave_initialize_topography
25 public shelfwave_set_obc_data
26 public register_shelfwave_obc, shelfwave_obc_end
27 
28 !> Control structure for shelfwave open boundaries.
29 type, public :: shelfwave_obc_cs ; private
30  real :: lx = 100.0 !< Long-shore length scale of bathymetry.
31  real :: ly = 50.0 !< Cross-shore length scale.
32  real :: f0 = 1.e-4 !< Coriolis parameter.
33  real :: jj = 1 !< Cross-shore wave mode.
34  real :: kk !< Parameter.
35  real :: ll !< Longshore wavenumber.
36  real :: alpha !< 1/Ly.
37  real :: omega !< Frequency.
38 end type shelfwave_obc_cs
39 
40 contains
41 
42 !> Add shelfwave to OBC registry.
43 function register_shelfwave_obc(param_file, CS, OBC_Reg)
44  type(param_file_type), intent(in) :: param_file !< parameter file.
45  type(shelfwave_obc_cs), pointer :: cs !< shelfwave control structure.
46  type(obc_registry_type), pointer :: obc_reg !< OBC registry.
47  logical :: register_shelfwave_obc
48  ! Local variables
49  real :: kk, ll, pi, len_lat
50 
51  character(len=32) :: casename = "shelfwave" !< This case's name.
52 
53  pi = 4.0*atan(1.0)
54 
55  if (associated(cs)) then
56  call mom_error(warning, "register_shelfwave_OBC called with an "// &
57  "associated control structure.")
58  return
59  endif
60  allocate(cs)
61 
62  ! Register the tracer for horizontal advection & diffusion.
63  call register_obc(casename, param_file, obc_reg)
64  call get_param(param_file, mdl,"F_0",cs%f0, &
65  do_not_log=.true.)
66  call get_param(param_file, mdl,"LENLAT",len_lat, &
67  do_not_log=.true.)
68  call get_param(param_file, mdl,"SHELFWAVE_X_WAVELENGTH",cs%Lx, &
69  "Length scale of shelfwave in x-direction.",&
70  units="Same as x,y", default=100.)
71  call get_param(param_file, mdl,"SHELFWAVE_Y_LENGTH_SCALE",cs%Ly, &
72  "Length scale of exponential dropoff of topography "//&
73  "in the y-direction.", &
74  units="Same as x,y", default=50.)
75  call get_param(param_file, mdl,"SHELFWAVE_Y_MODE",cs%jj, &
76  "Cross-shore wave mode.", &
77  units="nondim", default=1.)
78  cs%alpha = 1. / cs%Ly
79  cs%ll = 2. * pi / cs%Lx
80  cs%kk = cs%jj * pi / len_lat
81  cs%omega = 2 * cs%alpha * cs%f0 * cs%ll / &
82  (cs%kk*cs%kk + cs%alpha*cs%alpha + cs%ll*cs%ll)
83  register_shelfwave_obc = .true.
84 
85 end function register_shelfwave_obc
86 
87 !> Clean up the shelfwave OBC from registry.
88 subroutine shelfwave_obc_end(CS)
89  type(shelfwave_obc_cs), pointer :: cs !< shelfwave control structure.
90 
91  if (associated(cs)) then
92  deallocate(cs)
93  endif
94 end subroutine shelfwave_obc_end
95 
96 !> Initialization of topography.
97 subroutine shelfwave_initialize_topography( D, G, param_file, max_depth, US )
98  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
99  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
100  intent(out) :: d !< Ocean bottom depth in m or Z if US is present
101  type(param_file_type), intent(in) :: param_file !< Parameter file structure
102  real, intent(in) :: max_depth !< Maximum model depth in the units of D
103  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
104 
105  ! Local variables
106  real :: m_to_z ! A dimensional rescaling factor.
107  integer :: i, j
108  real :: y, rly, ly, h0
109 
110  m_to_z = 1.0 ; if (present(us)) m_to_z = us%m_to_Z
111 
112  call get_param(param_file, mdl,"SHELFWAVE_Y_LENGTH_SCALE",ly, &
113  default=50., do_not_log=.true.)
114  call get_param(param_file, mdl,"MINIMUM_DEPTH", h0, &
115  default=10., units="m", scale=m_to_z, do_not_log=.true.)
116 
117  rly = 0. ; if (ly>0.) rly = 1. / ly
118 
119  do j=g%jsc,g%jec ; do i=g%isc,g%iec
120  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
121  y = ( g%geoLatT(i,j) - g%south_lat )
122  d(i,j) = h0 * exp(2 * rly * y)
123  enddo ; enddo
124 
125 end subroutine shelfwave_initialize_topography
126 
127 !> This subroutine sets the properties of flow at open boundary conditions.
128 subroutine shelfwave_set_obc_data(OBC, CS, G, h, Time)
129  type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies
130  !! whether, where, and what open boundary
131  !! conditions are used.
132  type(shelfwave_obc_cs), pointer :: cs !< tidal bay control structure.
133  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
134  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< layer thickness.
135  type(time_type), intent(in) :: time !< model time.
136 
137  ! The following variables are used to set up the transport in the shelfwave example.
138  real :: my_amp, time_sec
139  real :: cos_wt, cos_ky, sin_wt, sin_ky, omega, alpha
140  real :: x, y, jj, kk, ll
141  character(len=40) :: mdl = "shelfwave_set_OBC_data" ! This subroutine's name.
142  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, n
143  integer :: isdb, iedb, jsdb, jedb
144  type(obc_segment_type), pointer :: segment => null()
145 
146  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
147  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
148  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
149 
150  if (.not.associated(obc)) return
151 
152  time_sec = time_type_to_real(time)
153  omega = cs%omega
154  alpha = cs%alpha
155  my_amp = 1.0
156  jj = cs%jj
157  kk = cs%kk
158  ll = cs%ll
159  do n = 1, obc%number_of_segments
160  segment => obc%segment(n)
161  if (.not. segment%on_pe) cycle
162  if (segment%direction /= obc_direction_w) cycle
163 
164  isdb = segment%HI%IsdB ; iedb = segment%HI%IedB
165  jsd = segment%HI%jsd ; jed = segment%HI%jed
166  do j=jsd,jed ; do i=isdb,iedb
167  x = g%geoLonCu(i,j) - g%west_lon
168  y = g%geoLatCu(i,j) - g%south_lat
169  sin_wt = sin(ll*x - omega*time_sec)
170  cos_wt = cos(ll*x - omega*time_sec)
171  sin_ky = sin(kk * y)
172  cos_ky = cos(kk * y)
173  segment%normal_vel_bt(i,j) = my_amp * exp(- alpha * y) * cos_wt * &
174  (alpha * sin_ky + kk * cos_ky)
175 ! segment%tangential_vel_bt(I,j) = my_amp * ll * exp(- alpha * y) * sin_wt * sin_ky
176 ! segment%vorticity_bt(I,j) = my_amp * exp(- alpha * y) * cos_wt * sin_ky&
177 ! (ll*ll + kk*kk + alpha*alpha)
178  enddo ; enddo
179  enddo
180 
181 end subroutine shelfwave_set_obc_data
182 
183 end module shelfwave_initialization
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.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
shelfwave_initialization
Configures the model for the idealized shelfwave test case.
Definition: shelfwave_initialization.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_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:272
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
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
shelfwave_initialization::shelfwave_obc_cs
Control structure for shelfwave open boundaries.
Definition: shelfwave_initialization.F90:29
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:186
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:107
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:22
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25