14 use regrid_consts,
only : coordinatemode, default_coordinate_mode
18 implicit none ;
private
20 character(len=40) :: mdl =
"adjustment_initialization"
22 #include <MOM_memory.h>
24 public adjustment_initialize_thickness
25 public adjustment_initialize_temperature_salinity
35 subroutine adjustment_initialize_thickness ( h, G, GV, US, param_file, just_read_params)
39 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
43 logical,
optional,
intent(in) :: just_read_params
48 real :: eta1d(szk_(g)+1)
53 real :: delta_s_strat, dsdz, delta_s, s_ref
54 real :: min_thickness, adjustment_width, adjustment_delta
55 real :: adjustment_deltas
56 real :: front_wave_amp, front_wave_length, front_wave_asym
57 real :: target_values(szk_(g)+1)
59 character(len=20) :: verticalcoordinate
61 #include "version_variable.h"
62 integer :: i, j, k, is, ie, js, je, nz
64 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
66 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
69 call mom_mesg(
"initialize_thickness_uniform: setting thickness")
72 if (.not.just_read)
call log_version(param_file, mdl, version,
"")
73 call get_param(param_file, mdl,
"S_REF", s_ref,
'Reference salinity', &
74 default=35.0, units=
'1e-3', do_not_log=just_read)
75 call get_param(param_file, mdl,
"MIN_THICKNESS",min_thickness,
'Minimum layer thickness', &
76 default=1.0e-3, units=
'm', scale=us%m_to_Z, do_not_log=just_read)
79 call get_param(param_file, mdl,
"REGRIDDING_COORDINATE_MODE",verticalcoordinate, &
80 default=default_coordinate_mode, do_not_log=just_read)
81 call get_param(param_file, mdl,
"ADJUSTMENT_WIDTH",adjustment_width, &
82 "Width of frontal zone", &
83 units=
"same as x,y", fail_if_missing=.not.just_read, do_not_log=just_read)
84 call get_param(param_file, mdl,
"DELTA_S_STRAT",delta_s_strat, &
85 "Top-to-bottom salinity difference of stratification", &
86 units=
"1e-3", fail_if_missing=.not.just_read, do_not_log=just_read)
87 call get_param(param_file, mdl,
"ADJUSTMENT_DELTAS",adjustment_deltas, &
88 "Salinity difference across front", &
89 units=
"1e-3", fail_if_missing=.not.just_read, do_not_log=just_read)
90 call get_param(param_file, mdl,
"FRONT_WAVE_AMP",front_wave_amp, &
91 "Amplitude of trans-frontal wave perturbation", &
92 units=
"same as x,y", default=0., do_not_log=just_read)
93 call get_param(param_file, mdl,
"FRONT_WAVE_LENGTH",front_wave_length, &
94 "Wave-length of trans-frontal wave perturbation", &
95 units=
"same as x,y", default=0., do_not_log=just_read)
96 call get_param(param_file, mdl,
"FRONT_WAVE_ASYM",front_wave_asym, &
97 "Amplitude of frontal asymmetric perturbation", &
98 units=
"same as x,y", default=0., do_not_log=just_read)
100 if (just_read)
return
110 dsdz = -delta_s_strat / g%max_depth
112 select case ( coordinatemode(verticalcoordinate) )
114 case ( regridding_layer, regridding_rho )
115 drho_ds = 1.0 * us%kg_m3_to_R
116 if (delta_s_strat /= 0.)
then
118 adjustment_delta = (adjustment_deltas / delta_s_strat) * g%max_depth
120 e0(k) = adjustment_delta - (g%max_depth + 2*adjustment_delta) * (real(k-1) / real(nz))
123 adjustment_delta = 2.*g%max_depth
125 e0(k) = -g%max_depth * (real(k-1) / real(nz))
129 target_values(1) = ( gv%Rlay(1) + 0.5*(gv%Rlay(1)-gv%Rlay(2)) )
130 target_values(nz+1) = ( gv%Rlay(nz) + 0.5*(gv%Rlay(nz)-gv%Rlay(nz-1)) )
132 target_values(1) = 0.0
133 target_values(nz+1) = 2.0 * gv%Rlay(1)
136 target_values(k) = target_values(k-1) + ( gv%Rlay(nz) - gv%Rlay(1) ) / (nz-1)
138 target_values(:) = target_values(:) - 1000.*us%kg_m3_to_R
139 do j=js,je ;
do i=is,ie
140 if (front_wave_length /= 0.)
then
141 y = ( 0.125 + g%geoLatT(i,j) / front_wave_length ) * ( 4. * acos(0.) )
142 yy = 2. * ( g%geoLatT(i,j) - 0.5 * g%len_lat ) / adjustment_width
143 yy = min(1.0, yy); yy = max(-1.0, yy)
144 yy = yy * 2. * acos( 0. )
145 y = front_wave_amp*sin(y) + front_wave_asym*sin(yy)
149 x = ( ( g%geoLonT(i,j) - 0.5 * g%len_lon ) + y ) / adjustment_width
150 x = min(1.0, x); x = max(-1.0, x)
152 delta_s = adjustment_deltas * 0.5 * (1. - sin( x ) )
154 if (drho_ds*dsdz /= 0.)
then
155 eta1d(k) = ( target_values(k) - drho_ds*( s_ref + delta_s ) ) / (drho_ds*dsdz)
157 eta1d(k) = e0(k) - (0.5*adjustment_delta) * sin( x )
159 eta1d(k) = max( eta1d(k), -g%max_depth )
160 eta1d(k) = min( eta1d(k), 0. )
162 eta1d(1) = 0.; eta1d(nz+1) = -g%max_depth
164 if (eta1d(k) > 0.)
then
165 eta1d(k) = max( eta1d(k+1) + min_thickness, 0. )
166 h(i,j,k) = gv%Z_to_H * max( eta1d(k) - eta1d(k+1), min_thickness )
167 elseif (eta1d(k) <= (eta1d(k+1) + min_thickness))
then
168 eta1d(k) = eta1d(k+1) + min_thickness
169 h(i,j,k) = gv%Z_to_H * min_thickness
171 h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
176 case ( regridding_zstar, regridding_sigma )
178 eta1d(k) = -g%max_depth * (real(k-1) / real(nz))
179 eta1d(k) = max(min(eta1d(k), 0.), -g%max_depth)
181 do j=js,je ;
do i=is,ie
183 h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
188 call mom_error(fatal,
"adjustment_initialize_thickness: "// &
189 "Unrecognized i.c. setup - set ADJUSTMENT_IC")
193 end subroutine adjustment_initialize_thickness
196 subroutine adjustment_initialize_temperature_salinity(T, S, h, G, GV, param_file, &
197 eqn_of_state, just_read_params)
200 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(out) :: t
201 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(out) :: s
202 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(in) :: h
205 type(
eos_type),
pointer :: eqn_of_state
206 logical,
optional,
intent(in) :: just_read_params
209 integer :: i, j, k, is, ie, js, je, nz
211 integer :: index_bay_z
214 real :: s_range, t_range
216 real :: xi0, xi1, dsdz, delta_s, delta_s_strat
217 real :: adjustment_width, adjustment_deltas
218 real :: front_wave_amp, front_wave_length, front_wave_asym
219 real :: eta1d(szk_(g)+1)
221 character(len=20) :: verticalcoordinate
223 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
225 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
228 call get_param(param_file, mdl,
"S_REF", s_ref,
'Reference salinity', &
229 default=35.0, units=
'1e-3', do_not_log=just_read)
230 call get_param(param_file, mdl,
"T_REF",t_ref,
'Reference temperature', units=
'C', &
231 fail_if_missing=.not.just_read, do_not_log=just_read)
232 call get_param(param_file, mdl,
"S_RANGE",s_range,
'Initial salinity range', units=
'1e-3', &
233 default=2.0, do_not_log=just_read)
234 call get_param(param_file, mdl,
"T_RANGE",t_range,
'Initial temperature range', units=
'C', &
235 default=0.0, do_not_log=just_read)
237 call get_param(param_file, mdl,
"REGRIDDING_COORDINATE_MODE",verticalcoordinate, &
238 default=default_coordinate_mode, do_not_log=just_read)
239 call get_param(param_file, mdl,
"ADJUSTMENT_WIDTH", adjustment_width, &
240 fail_if_missing=.not.just_read, do_not_log=.true.)
241 call get_param(param_file, mdl,
"ADJUSTMENT_DELTAS", adjustment_deltas, &
242 fail_if_missing=.not.just_read, do_not_log=.true.)
243 call get_param(param_file, mdl,
"DELTA_S_STRAT", delta_s_strat, &
244 fail_if_missing=.not.just_read, do_not_log=.true.)
245 call get_param(param_file, mdl,
"FRONT_WAVE_AMP", front_wave_amp, default=0., &
247 call get_param(param_file, mdl,
"FRONT_WAVE_LENGTH",front_wave_length, &
248 default=0., do_not_log=.true.)
249 call get_param(param_file, mdl,
"FRONT_WAVE_ASYM", front_wave_asym, default=0., &
252 if (just_read)
return
258 select case ( coordinatemode(verticalcoordinate) )
260 case ( regridding_zstar, regridding_sigma )
261 dsdz = -delta_s_strat / g%max_depth
262 do j=js,je ;
do i=is,ie
263 eta1d(nz+1) = -g%bathyT(i,j)
265 eta1d(k) = eta1d(k+1) + h(i,j,k)*gv%H_to_Z
267 if (front_wave_length /= 0.)
then
268 y = ( 0.125 + g%geoLatT(i,j) / front_wave_length ) * ( 4. * acos(0.) )
269 yy = 2. * ( g%geoLatT(i,j) - 0.5 * g%len_lat ) / front_wave_length
270 yy = min(1.0, yy); yy = max(-1.0, yy)
271 yy = yy * 2. * acos( 0. )
272 y = front_wave_amp*sin(y) + front_wave_asym*sin(yy)
276 x = ( ( g%geoLonT(i,j) - 0.5 * g%len_lon ) + y ) / adjustment_width
277 x = min(1.0, x); x = max(-1.0, x)
279 delta_s = adjustment_deltas * 0.5 * (1. - sin( x ) )
281 s(i,j,k) = s_ref + delta_s + 0.5 * ( eta1d(k)+eta1d(k+1) ) * dsdz
282 x = abs(s(i,j,k) - 0.5*real(nz-1)/real(nz)*s_range)/s_range*real(2*nz)
290 case ( regridding_layer, regridding_rho )
292 s(:,:,k) = s_ref + s_range * ( (real(k)-0.5) / real( nz ) )
299 call mom_error(fatal,
"adjustment_initialize_temperature_salinity: "// &
300 "Unrecognized i.c. setup - set ADJUSTMENT_IC")
304 end subroutine adjustment_initialize_temperature_salinity