MOM6
shelfwave_initialization Module Reference

Detailed Description

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, US)
 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"
 This module's name. More...
 

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 44 of file shelfwave_initialization.F90.

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 

References mdl, and mom_error_handler::mom_error().

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,
type(unit_scale_type), intent(in), optional  US 
)

Initialization of topography.

Parameters
[in]gThe dynamic horizontal grid type
[out]dOcean bottom depth in m or Z if US is present
[in]param_fileParameter file structure
[in]max_depthMaximum model depth in the units of D
[in]usA dimensional unit scaling type

Definition at line 98 of file shelfwave_initialization.F90.

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 

References mdl.

Referenced by mom_fixed_initialization::mom_initialize_topography().

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 89 of file shelfwave_initialization.F90.

89  type(shelfwave_OBC_CS), pointer :: CS !< shelfwave control structure.
90 
91  if (associated(cs)) then
92  deallocate(cs)
93  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 129 of file shelfwave_initialization.F90.

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) = g%US%m_s_to_L_T*my_amp * exp(- alpha * y) * cos_wt * &
174  (alpha * sin_ky + kk * cos_ky)
175 ! segment%tangential_vel_bt(I,j) = G%US%m_s_to_L_T*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 

Variable Documentation

◆ mdl

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

This module's name.

Definition at line 21 of file shelfwave_initialization.F90.

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

Referenced by register_shelfwave_obc(), and shelfwave_initialize_topography().