MOM6
shelfwave_initialization Module Reference

Detailed Description

The module configures the model for the idealized shelfwave test case.

Data Types

type  shelfwave_obc_cs
 Control structure for shelfwave open boundaries. More...
 

Functions/Subroutines

logical function, public register_shelfwave_obc (param_file, CS, OBC_Reg)
 Add shelfwave to OBC registry. More...
 
subroutine, public shelfwave_obc_end (CS)
 Clean up the shelfwave OBC from registry. More...
 
subroutine, public shelfwave_initialize_topography (D, G, param_file, max_depth)
 Initialization of topography. More...
 
subroutine, public shelfwave_set_obc_data (OBC, CS, G, h, Time)
 This subroutine sets the properties of flow at open boundary conditions. More...
 

Variables

character(len=40) mdl = "shelfwave_initialization"
 

Function/Subroutine Documentation

◆ register_shelfwave_obc()

logical function, public shelfwave_initialization::register_shelfwave_obc ( type(param_file_type), intent(in)  param_file,
type(shelfwave_obc_cs), pointer  CS,
type(obc_registry_type), pointer  OBC_Reg 
)

Add shelfwave to OBC registry.

Parameters
[in]param_fileparameter file.
csshelfwave control structure.
obc_regOBC registry.

Definition at line 61 of file shelfwave_initialization.F90.

References mdl, and mom_error_handler::mom_error().

61  type(param_file_type), intent(in) :: param_file !< parameter file.
62  type(shelfwave_obc_cs), pointer :: cs !< shelfwave control structure.
63  type(obc_registry_type), pointer :: obc_reg !< OBC registry.
64  logical :: register_shelfwave_obc
65  ! Local variables
66  real :: kk, ll, pi, len_lat
67 
68  character(len=32) :: casename = "shelfwave" !< This case's name.
69 
70  pi = 4.0*atan(1.0)
71 
72  if (associated(cs)) then
73  call mom_error(warning, "register_shelfwave_OBC called with an "// &
74  "associated control structure.")
75  return
76  endif
77  allocate(cs)
78 
79  ! Register the tracer for horizontal advection & diffusion.
80  call register_obc(casename, param_file, obc_reg)
81  call get_param(param_file, mdl,"F_0",cs%f0, &
82  do_not_log=.true.)
83  call get_param(param_file, mdl,"LENLAT",len_lat, &
84  do_not_log=.true.)
85  call get_param(param_file, mdl,"SHELFWAVE_X_WAVELENGTH",cs%Lx, &
86  "Length scale of shelfwave in x-direction.",&
87  units="Same as x,y", default=100.)
88  call get_param(param_file, mdl,"SHELFWAVE_Y_LENGTH_SCALE",cs%Ly, &
89  "Length scale of exponential dropoff of topography\n"//&
90  "in the y-direction.", &
91  units="Same as x,y", default=50.)
92  call get_param(param_file, mdl,"SHELFWAVE_Y_MODE",cs%jj, &
93  "Cross-shore wave mode.", &
94  units="nondim", default=1.)
95  cs%alpha = 1. / cs%Ly
96  cs%ll = 2. * pi / cs%Lx
97  cs%kk = cs%jj * pi / len_lat
98  cs%omega = 2 * cs%alpha * cs%f0 * cs%ll / &
99  (cs%kk*cs%kk + cs%alpha*cs%alpha + cs%ll*cs%ll)
100  register_shelfwave_obc = .true.
101 
Here is the call graph for this function:

◆ shelfwave_initialize_topography()

subroutine, public shelfwave_initialization::shelfwave_initialize_topography ( real, dimension(g%isd:g%ied,g%jsd:g%jed), intent(out)  D,
type(dyn_horgrid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
real, intent(in)  max_depth 
)

Initialization of topography.

Parameters
[in]gThe dynamic horizontal grid type
[out]dOcean bottom depth in m
[in]param_fileParameter file structure
[in]max_depthMaximum depth of model in m

Definition at line 115 of file shelfwave_initialization.F90.

References mdl.

Referenced by mom_fixed_initialization::mom_initialize_topography().

115  ! Arguments
116  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
117  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
118  intent(out) :: d !< Ocean bottom depth in m
119  type(param_file_type), intent(in) :: param_file !< Parameter file structure
120  real, intent(in) :: max_depth !< Maximum depth of model in m
121 
122  ! Local variables
123  integer :: i, j
124  real :: y, rly, ly, h0
125 
126  call get_param(param_file, mdl,"SHELFWAVE_Y_LENGTH_SCALE",ly, &
127  default=50., do_not_log=.true.)
128  call get_param(param_file, mdl,"MINIMUM_DEPTH",h0, &
129  default=10., do_not_log=.true.)
130 
131  rly = 0. ; if (ly>0.) rly = 1. / ly
132  do i=g%isc,g%iec
133  do j=g%jsc,g%jec
134  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
135  y = ( g%geoLatT(i,j) - g%south_lat )
136  d(i,j) = h0 * exp(2 * rly * y)
137  enddo
138  enddo
139 
Here is the caller graph for this function:

◆ shelfwave_obc_end()

subroutine, public shelfwave_initialization::shelfwave_obc_end ( type(shelfwave_obc_cs), pointer  CS)

Clean up the shelfwave OBC from registry.

Parameters
csshelfwave control structure.

Definition at line 106 of file shelfwave_initialization.F90.

106  type(shelfwave_obc_cs), pointer :: cs !< shelfwave control structure.
107 
108  if (associated(cs)) then
109  deallocate(cs)
110  endif

◆ shelfwave_set_obc_data()

subroutine, public shelfwave_initialization::shelfwave_set_obc_data ( type(ocean_obc_type), pointer  OBC,
type(shelfwave_obc_cs), pointer  CS,
type(ocean_grid_type), intent(in)  G,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(time_type), intent(in)  Time 
)

This subroutine sets the properties of flow at open boundary conditions.

Parameters
obcThis open boundary condition type specifies whether, where, and what open boundary conditions are used.
cstidal bay control structure.
[in]gThe ocean's grid structure.
[in]hlayer thickness.
[in]timemodel time.

Definition at line 144 of file shelfwave_initialization.F90.

144  type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies
145  !! whether, where, and what open boundary
146  !! conditions are used.
147  type(shelfwave_obc_cs), pointer :: cs !< tidal bay control structure.
148  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
149  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< layer thickness.
150  type(time_type), intent(in) :: time !< model time.
151 
152  ! The following variables are used to set up the transport in the tidal_bay example.
153  real :: my_amp, time_sec
154  real :: cos_wt, cos_ky, sin_wt, sin_ky, omega, alpha
155  real :: x, y, jj, kk, ll
156  character(len=40) :: mdl = "shelfwave_set_OBC_data" ! This subroutine's name.
157  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, n
158  integer :: isdb, iedb, jsdb, jedb
159  type(obc_segment_type), pointer :: segment
160 
161  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
162  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
163  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
164 
165  if (.not.associated(obc)) return
166 
167  time_sec = time_type_to_real(time)
168  omega = cs%omega
169  alpha = cs%alpha
170  my_amp = 1.0
171  jj = cs%jj
172  kk = cs%kk
173  ll = cs%ll
174  do n = 1, obc%number_of_segments
175  segment => obc%segment(n)
176  if (.not. segment%on_pe) cycle
177  if (segment%direction /= obc_direction_w) cycle
178 
179  isdb = segment%HI%IsdB ; iedb = segment%HI%IedB
180  jsd = segment%HI%jsd ; jed = segment%HI%jed
181  do j=jsd,jed ; do i=isdb,iedb
182  x = g%geoLonCu(i,j) - g%west_lon
183  y = g%geoLatCu(i,j) - g%south_lat
184  sin_wt = sin(ll*x - omega*time_sec)
185  cos_wt = cos(ll*x - omega*time_sec)
186  sin_ky = sin(kk * y)
187  cos_ky = cos(kk * y)
188  segment%normal_vel_bt(i,j) = my_amp * exp(- alpha * y) * cos_wt * &
189  (alpha * sin_ky + kk * cos_ky)
190 ! segment%tangential_vel_bt(I,j) = my_amp * ll * exp(- alpha * y) * sin_wt * sin_ky
191 ! segment%vorticity_bt(I,j) = my_amp * exp(- alpha * y) * cos_wt * sin_ky&
192 ! (ll*ll + kk*kk + alpha*alpha)
193  enddo ; enddo
194  enddo
195 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19

Variable Documentation

◆ mdl

character(len=40) shelfwave_initialization::mdl = "shelfwave_initialization"
private

Definition at line 36 of file shelfwave_initialization.F90.

Referenced by register_shelfwave_obc(), and shelfwave_initialize_topography().

36 character(len=40) :: mdl = "shelfwave_initialization" ! This module's name.