19 implicit none ;
private 21 #include <MOM_memory.h> 24 public sloshing_initialize_topography
25 public sloshing_initialize_thickness
26 public sloshing_initialize_temperature_salinity
31 subroutine sloshing_initialize_topography( D, G, param_file, max_depth )
33 real,
dimension(G%isd:G%ied,G%jsd:G%jed), &
36 real,
intent(in) :: max_depth
41 do i=g%isc,g%iec ;
do j=g%jsc,g%jec
45 end subroutine sloshing_initialize_topography
56 subroutine sloshing_initialize_thickness ( h, G, GV, US, param_file, just_read_params)
60 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
64 logical,
optional,
intent(in) :: just_read_params
67 real :: displ(szk_(g)+1)
68 real :: z_unif(szk_(g)+1)
69 real :: z_inter(szk_(g)+1)
72 real :: x1, y1, x2, y2
77 # include "version_variable.h" 78 character(len=40) :: mdl =
"sloshing_initialization" 80 integer :: i, j, k, is, ie, js, je, nx, nz
82 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
84 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
85 if (.not.just_read)
call log_version(param_file, mdl, version,
"")
86 call get_param(param_file, mdl,
"SLOSHING_IC_AMPLITUDE", a0, &
87 "Initial amplitude of sloshing internal interface height "//&
88 "displacements it the sloshing test case.", &
89 units=
'm', default=75.0, scale=us%m_to_Z, do_not_log=just_read)
90 call get_param(param_file, mdl,
"SLOSHING_IC_BUG", use_ic_bug, &
91 "If true, use code with a bug to set the sloshing initial conditions.", &
92 default=.false., do_not_log=just_read)
97 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
101 z_unif(k+1) = -
real(k)/
real(nz)
111 x1 = 0.30; y1 = 0.48; x2 = 0.70; y2 = 0.52
117 elseif ( (x > x1 ) .and. ( x < x2 ))
then 118 t = y1 + (y2-y1) * (x-x1) / (x2-x1)
120 t = y2 + (1.0-y2) * (x-x2) / (1.0-x2)
125 z_inter(k) = -t * g%max_depth
133 weight_z = - 4.0 * ( z_unif(k) + 0.5 )**2 + 1.0
135 x = g%geoLonT(i,j) / g%len_lon
137 displ(k) = a0 * cos(acos(-1.0)*x) + weight_z * us%m_to_Z
139 displ(k) = a0 * cos(acos(-1.0)*x) * weight_z
146 if ( k == nz+1 )
then 150 z_inter(k) = z_inter(k) + displ(k)
155 z_inter(nz+1) = -g%bathyT(i,j)
158 if ( z_inter(k) < (z_inter(k+1) + gv%Angstrom_Z) )
then 159 z_inter(k) = z_inter(k+1) + gv%Angstrom_Z
165 h(i,j,k) = gv%Z_to_H * (z_inter(k) - z_inter(k+1))
170 end subroutine sloshing_initialize_thickness
179 subroutine sloshing_initialize_temperature_salinity ( T, S, h, G, GV, param_file, &
180 eqn_of_state, just_read_params)
183 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(out) :: t
184 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(out) :: s
185 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(in) :: h
189 type(
eos_type),
pointer :: eqn_of_state
190 logical,
optional,
intent(in) :: just_read_params
193 integer :: i, j, k, is, ie, js, je, nz
194 real :: delta_s, delta_t
195 real :: s_ref, t_ref;
197 real :: s_range, t_range;
203 character(len=40) :: mdl =
"initialize_temp_salt_linear" 206 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
208 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
210 call get_param(param_file, mdl,
"S_REF", s_ref,
'Reference value for salinity', &
211 default=35.0, units=
'1e-3', do_not_log=just_read)
212 call get_param(param_file, mdl,
"T_REF", t_ref,
'Reference value for temperature', &
213 units=
'degC', fail_if_missing=.not.just_read, do_not_log=just_read)
216 call get_param(param_file, mdl,
"S_RANGE",s_range,
'Initial salinity range.', &
217 units=
'1e-3', default=2.0, do_not_log=just_read)
218 call get_param(param_file, mdl,
"T_RANGE",t_range,
'Initial temperature range', &
219 units=
'degC', default=0.0, do_not_log=just_read)
221 if (just_read)
return 231 deltah = g%max_depth / nz
232 do j=js,je ;
do i=is,ie
235 xi1 = xi0 + deltah / g%max_depth
236 s(i,j,k) = 34.0 + 0.5 * s_range * (xi0 + xi1)
242 delta_t = t_range / ( g%ke - 1.0 )
246 t(:,:,k) = t(:,:,k-1) + delta_t
249 t(:,:,g%ke/2 - (kdelta-1):g%ke/2 + kdelta) = 1.0
251 end subroutine sloshing_initialize_temperature_salinity
Ocean grid type. See mom_grid for details.
Calculates density of sea water from T, S and P.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Describes the horizontal ocean grid with only dynamic memory arrays.
The MOM6 facility to parse input files for runtime parameters.
Describes various unit conversion factors.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Type to carry basic tracer information.
Routines for error handling and I/O management.
This control structure holds memory and parameters for the MOM_sponge module.
Implements sponge regions in isopycnal mode.
Initialization for the "sloshing" internal waves configuration.
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Provides subroutines for quantities specific to the equation of state.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
A control structure for the equation of state.