18 implicit none ;
private 20 #include <MOM_memory.h> 22 public benchmark_initialize_topography
23 public benchmark_initialize_thickness
24 public benchmark_init_temperature_salinity
34 subroutine benchmark_initialize_topography(D, G, param_file, max_depth, US)
36 real,
dimension(G%isd:G%ied,G%jsd:G%jed), &
39 real,
intent(in) :: max_depth
50 # include "version_variable.h" 51 character(len=40) :: mdl =
"benchmark_initialize_topography" 52 integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
53 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
54 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
56 call mom_mesg(
" benchmark_initialization.F90, benchmark_initialize_topography: setting topography", 5)
58 m_to_z = 1.0 ;
if (
present(us)) m_to_z = us%m_to_Z
61 call get_param(param_file, mdl,
"MINIMUM_DEPTH", min_depth, &
62 "The minimum depth of the ocean.", units=
"m", default=0.0, scale=m_to_z)
68 do j=js,je ;
do i=is,ie
69 x = (g%geoLonT(i,j)-g%west_lon) / g%len_lon
70 y = (g%geoLatT(i,j)-g%south_lat) / g%len_lat
72 d(i,j) = -d0 * ( y*(1.0 + 0.6*cos(4.0*pi*x)) &
74 + 0.05*cos(10.0*pi*x) - 0.7 )
75 if (d(i,j) > max_depth) d(i,j) = max_depth
76 if (d(i,j) < min_depth) d(i,j) = 0.
79 end subroutine benchmark_initialize_topography
85 subroutine benchmark_initialize_thickness(h, G, GV, US, param_file, eqn_of_state, &
86 P_Ref, just_read_params)
90 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
94 type(
eos_type),
pointer :: eqn_of_state
95 real,
intent(in) :: P_Ref
97 logical,
optional,
intent(in) :: just_read_params
100 real :: e0(szk_(gv)+1)
102 real :: e_pert(szk_(gv)+1)
104 real :: eta1D(szk_(gv)+1)
109 real :: thermocline_scale
110 real,
dimension(SZK_(GV)) :: &
115 real :: pres(szk_(gv))
125 # include "version_variable.h" 126 character(len=40) :: mdl =
"benchmark_initialize_thickness" 127 integer :: i, j, k, k1, is, ie, js, je, nz, itt
129 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
131 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
132 if (.not.just_read)
call log_version(param_file, mdl, version,
"")
133 call get_param(param_file, mdl,
"BENCHMARK_ML_DEPTH_IC", ml_depth, &
134 "Initial mixed layer depth in the benchmark test case.", &
135 units=
'm', default=50.0, scale=us%m_to_Z, do_not_log=just_read)
136 call get_param(param_file, mdl,
"BENCHMARK_THERMOCLINE_SCALE", thermocline_scale, &
137 "Initial thermocline depth scale in the benchmark test case.", &
138 default=500.0, units=
"m", scale=us%m_to_Z, do_not_log=just_read)
140 if (just_read)
return 142 call mom_mesg(
" benchmark_initialization.F90, benchmark_initialize_thickness: setting thickness", 5)
144 k1 = gv%nk_rho_varies + 1
151 pres(k) = p_ref ; s0(k) = 35.0
159 t0(k) = t0(k1) + (gv%Rlay(k) - rho_guess(k1)) / drho_dt(k1)
167 t0(k) = t0(k) + (gv%Rlay(k) - rho_guess(k)) / drho_dt(k)
172 i_ts = 1.0 / thermocline_scale
173 i_md = 1.0 / g%max_depth
174 do j=js,je ;
do i=is,ie
175 sst = 0.5*(t0(k1)+t0(nz)) - 0.9*0.5*(t0(k1)-t0(nz)) * &
176 cos(pi*(g%geoLatT(i,j)-g%south_lat)/(g%len_lat))
178 do k=1,nz ; e_pert(k) = 0.0 ;
enddo 184 eta1d(nz+1) = -g%bathyT(i,j)
187 t_int = 0.5*(t0(k) + t0(k-1))
188 t_frac = (t_int - t0(nz)) / (sst - t0(nz))
192 err = a_exp * exp(z*i_ts) + (1.0 - a_exp) * (z*i_md + 1.0) - t_frac
193 derr_dz = a_exp * i_ts * exp(z*i_ts) + (1.0 - a_exp) * i_md
194 z = z - err / derr_dz
199 eta1d(k) = e0(k) + e_pert(k)
201 if (eta1d(k) > -ml_depth) eta1d(k) = -ml_depth
203 if (eta1d(k) < eta1d(k+1) + gv%Angstrom_Z) &
204 eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
206 h(i,j,k) = max(gv%Z_to_H * (eta1d(k) - eta1d(k+1)), gv%Angstrom_H)
208 h(i,j,1) = max(gv%Z_to_H * (0.0 - eta1d(2)), gv%Angstrom_H)
212 end subroutine benchmark_initialize_thickness
215 subroutine benchmark_init_temperature_salinity(T, S, G, GV, US, param_file, &
216 eqn_of_state, P_Ref, just_read_params)
219 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: T
221 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: S
227 type(
eos_type),
pointer :: eqn_of_state
228 real,
intent(in) :: P_Ref
230 logical,
optional,
intent(in) :: just_read_params
233 real :: T0(szk_(g)), S0(szk_(g))
234 real :: pres(szk_(g))
235 real :: drho_dT(szk_(g))
236 real :: drho_dS(szk_(g))
237 real :: rho_guess(szk_(g))
242 character(len=40) :: mdl =
"benchmark_init_temperature_salinity" 243 integer :: i, j, k, k1, is, ie, js, je, nz, itt
245 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
247 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
249 if (just_read)
return 251 k1 = gv%nk_rho_varies + 1
254 pres(k) = p_ref ; s0(k) = 35.0
263 t0(k) = t0(k1) + (gv%Rlay(k) - rho_guess(k1)) / drho_dt(k1)
271 t0(k) = t0(k) + (gv%Rlay(k) - rho_guess(k)) / drho_dt(k)
275 do k=1,nz ;
do i=is,ie ;
do j=js,je
278 enddo ;
enddo ;
enddo 280 do i=is,ie ;
do j=js,je
281 sst = 0.5*(t0(k1)+t0(nz)) - 0.9*0.5*(t0(k1)-t0(nz)) * &
282 cos(pi*(g%geoLatT(i,j)-g%south_lat)/(g%len_lat))
288 end subroutine benchmark_init_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.
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.
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.
Initialization for the "bench mark" configuration.
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.