20 implicit none ;
private 22 #include <MOM_memory.h> 25 character(len=40) ::
mdl =
"Rossby_front_2d_initialization" 27 #include "version_variable.h" 44 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
48 logical,
optional,
intent(in) :: just_read_params
51 integer :: i, j, k, is, ie, js, je, nz
52 real :: Tz, Dml, eta, stretch, h0
53 real :: min_thickness, T_range, dRho_dT
55 character(len=40) :: verticalCoordinate
57 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
59 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
62 call mom_mesg(
"Rossby_front_2d_initialization.F90, Rossby_front_initialize_thickness: setting thickness")
66 call get_param(param_file,
mdl,
"MIN_THICKNESS", min_thickness, &
67 'Minimum layer thickness',units=
'm',default=1.e-3, do_not_log=just_read)
68 call get_param(param_file,
mdl,
"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
69 default=default_coordinate_mode, do_not_log=just_read)
70 call get_param(param_file,
mdl,
"T_RANGE", t_range,
'Initial temperature range', &
71 units=
'C', default=0.0, do_not_log=just_read)
72 call get_param(param_file,
mdl,
"DRHO_DT", drho_dt, default=-0.2, do_not_log=.true.)
76 tz = t_range / g%max_depth
78 select case ( coordinatemode(verticalcoordinate) )
81 do j = g%jsc,g%jec ;
do i = g%isc,g%iec
82 dml =
hml( g, g%geoLatT(i,j) )
83 eta = -( -drho_dt / gv%Rho0 ) * tz * 0.5 * ( dml * dml )
84 stretch = ( ( g%max_depth + eta ) / g%max_depth )
85 h0 = ( g%max_depth /
real(nz) ) * stretch
92 do j = g%jsc,g%jec ;
do i = g%isc,g%iec
93 dml =
hml( g, g%geoLatT(i,j) )
94 eta = -( -drho_dt / gv%Rho0 ) * tz * 0.5 * ( dml * dml )
95 stretch = ( ( g%max_depth + eta ) / g%max_depth )
96 h0 = ( g%max_depth /
real(nz) ) * stretch
103 call mom_error(fatal,
"Rossby_front_initialize: "// &
104 "Unrecognized i.c. setup - set REGRIDDING_COORDINATE_MODE")
113 param_file, eqn_of_state, just_read_params)
115 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(out) :: T
116 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(out) :: S
117 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(in) :: h
119 type(
eos_type),
pointer :: eqn_of_state
120 logical,
optional,
intent(in) :: just_read_params
123 integer :: i, j, k, is, ie, js, je, nz
126 real :: y, zc, zi, dTdz
128 character(len=40) :: verticalCoordinate
131 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
133 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
135 call get_param(param_file,
mdl,
"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
136 default=default_coordinate_mode, do_not_log=just_read)
137 call get_param(param_file,
mdl,
"S_REF",s_ref,
'Reference salinity', units=
'1e-3', &
138 fail_if_missing=.not.just_read, do_not_log=just_read)
139 call get_param(param_file,
mdl,
"T_REF",t_ref,
'Reference temperature',units=
'C',&
140 fail_if_missing=.not.just_read, do_not_log=just_read)
141 call get_param(param_file,
mdl,
"T_RANGE",t_range,
'Initial temperature range',&
142 units=
'C', default=0.0, do_not_log=just_read)
144 if (just_read)
return 148 dtdz = t_range / g%max_depth
150 do j = g%jsc,g%jec ;
do i = g%isc,g%iec
154 zc = zi - 0.5*h(i,j,k)
155 zc = min( zc, -
hml(g, g%geoLatT(i,j)) )
156 t(i,j,k) = t_ref + dtdz * zc
167 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(out) :: u
168 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(out) :: v
169 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(in) :: h
173 logical,
optional,
intent(in) :: just_read_params
179 real :: dRho_dT, zi, zc, zm, f, Ty, Dml, hAtU
180 integer :: i, j, k, is, ie, js, je, nz
182 character(len=40) :: verticalCoordinate
184 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
186 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
188 call get_param(param_file,
mdl,
"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
189 default=default_coordinate_mode, do_not_log=just_read)
190 call get_param(param_file,
mdl,
"T_RANGE", t_range,
'Initial temperature range', &
191 units=
'C', default=0.0, do_not_log=just_read)
192 call get_param(param_file,
mdl,
"DRHO_DT", drho_dt, default=-0.2, do_not_log=.true.)
194 if (just_read)
return 199 do j = g%jsc,g%jec ;
do i = g%isc-1,g%iec+1
200 f = 0.5*( g%CoriolisBu(i,j) + g%CoriolisBu(i,j-1) )
201 dudt = 0.0 ;
if (abs(f) > 0.0) &
202 dudt = ( gv%g_Earth * drho_dt ) / ( f * gv%Rho0 )
203 dml =
hml( g, g%geoLatT(i,j) )
204 ty =
dtdy( g, t_range, g%geoLatT(i,j) )
207 hatu = 0.5*(h(i,j,k)+h(i+1,j,k))
210 zm = max( zc + dml, 0. )
211 u(i,j,k) = dudt * ty * zm
220 real function ypseudo( G, lat )
222 real,
intent(in) :: lat
227 ypseudo = ( ( lat - g%south_lat ) / g%len_lat ) - 0.5
233 real function hml( G, lat )
235 real,
intent(in) :: lat
237 real :: dHML, HMLmean
246 real function dtdy( G, dT, lat )
248 real,
intent(in) :: dT
249 real,
intent(in) :: lat
251 real :: PI, dHML, dHdy
257 dtdy = -( dt / g%max_depth ) * dhdy
integer, parameter regridding_layer
Layer mode.
subroutine, public read_axis_data(filename, axis_name, var)
character(len=40) mdl
This module's name.
integer function coordinatemode(string)
Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer...
integer, parameter regridding_sigma
Sigma coordinates.
real, parameter frontfractionalwidth
Width of front as fraction of domain.
Ocean grid type. See mom_grid for details.
Calculates density of sea water from T, S and P.
Provides the ocean grid type.
character(len= *), parameter default_coordinate_mode
Default coordinate mode.
real, parameter hmlmax
Deepest ML as fractional depth of ocean.
This module contains I/O framework code.
Initial conditions for the 2D Rossby front test.
subroutine, public rossby_front_initialize_velocity(u, v, h, G, GV, param_file, just_read_params)
Initialization of u and v in the Rossby front test.
subroutine, public calculate_density_derivs(T, S, pressure, drho_dT, drho_dS, start, npts, EOS)
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
logical function, public is_root_pe()
subroutine, public rossby_front_initialize_temperature_salinity(T, S, h, G, param_file, eqn_of_state, just_read_params)
Initialization of temperature and salinity in the Rossby front test.
subroutine, public rossby_front_initialize_thickness(h, G, GV, param_file, just_read_params)
Initialization of thicknesses in 2D Rossby front test.
subroutine, public mom_mesg(message, verb, all_print)
Provides subroutines for quantities specific to the equation of state.
Type for describing a variable, typically a tracer.
integer, parameter regridding_zstar
z* coordinates
real, parameter hmlmin
Shallowest ML as fractional depth of ocean.
real function hml(G, lat)
Analytic prescription of mixed layer depth in 2d Rossby front test.
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
Contains constants for interpreting input parameters that control regridding.
real function dtdy(G, dT, lat)
Analytic prescription of mixed layer temperature gradient in 2d Rossby front test.
subroutine, public mom_error(level, message, all_print)
integer, parameter regridding_rho
Target interface densities.
real function ypseudo(G, lat)
Pseudo coordinate across domain used by Hml() and dTdy() returns a coordinate from -PI/2 ...
A control structure for the equation of state.