This subroutine determines the number of points which are within sponges in this computational domain. Only points that have positive values of Iresttime and which mask2dT indicates are ocean points are included in the sponges. It also stores the target interface heights.
150 type(ocean_grid_type),
intent(in) :: g
151 integer,
intent(in) :: nz_data
152 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: iresttime
153 type(param_file_type),
intent(in) :: param_file
155 type(ale_sponge_cs),
pointer :: cs
157 real,
dimension(SZI_(G),SZJ_(G),nz_data),
intent(in) :: data_h
162 #include "version_variable.h" 163 character(len=40) :: mdl =
"MOM_sponge" 164 logical :: use_sponge
165 real,
allocatable,
dimension(:,:,:) :: data_hu
166 real,
allocatable,
dimension(:,:,:) :: data_hv
167 real,
allocatable,
dimension(:,:) :: iresttime_u
168 real,
allocatable,
dimension(:,:) :: iresttime_v
169 logical :: bndextrapolation = .true.
170 logical :: default_2018_answers
171 integer :: i, j, k, col, total_sponge_cols, total_sponge_cols_u, total_sponge_cols_v
172 character(len=10) :: remapscheme
173 if (
associated(cs))
then 174 call mom_error(warning,
"initialize_ALE_sponge_fixed called with an associated "// &
175 "control structure.")
180 call log_version(param_file, mdl, version,
"")
181 call get_param(param_file, mdl,
"SPONGE", use_sponge, &
182 "If true, sponges may be applied anywhere in the domain. "//&
183 "The exact location and properties of those sponges are "//&
184 "specified from MOM_initialization.F90.", default=.false.)
186 if (.not.use_sponge)
return 190 call get_param(param_file, mdl,
"SPONGE_UV", cs%sponge_uv, &
191 "Apply sponges in u and v, in addition to tracers.", &
194 call get_param(param_file, mdl,
"REMAPPING_SCHEME", remapscheme, &
195 "This sets the reconstruction scheme used "//&
196 " for vertical remapping for all variables.", &
197 default=
"PLM", do_not_log=.true.)
199 call get_param(param_file, mdl,
"BOUNDARY_EXTRAPOLATION", bndextrapolation, &
200 "When defined, a proper high-order reconstruction "//&
201 "scheme is used within boundary cells rather "//&
202 "than PCM. E.g., if PPM is used for remapping, a "//&
203 "PPM reconstruction will also be used within boundary cells.", &
204 default=.false., do_not_log=.true.)
205 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
206 "This sets the default value for the various _2018_ANSWERS parameters.", &
208 call get_param(param_file, mdl,
"REMAPPING_2018_ANSWERS", cs%remap_answers_2018, &
209 "If true, use the order of arithmetic and expressions that recover the "//&
210 "answers from the end of 2018. Otherwise, use updated and more robust "//&
211 "forms of the same expressions.", default=default_2018_answers)
212 call get_param(param_file, mdl,
"HOR_REGRID_2018_ANSWERS", cs%hor_regrid_answers_2018, &
213 "If true, use the order of arithmetic for horizonal regridding that recovers "//&
214 "the answers from the end of 2018. Otherwise, use rotationally symmetric "//&
215 "forms of the same expressions.", default=default_2018_answers)
217 cs%time_varying_sponges = .false.
219 cs%isc = g%isc ; cs%iec = g%iec ; cs%jsc = g%jsc ; cs%jec = g%jec
220 cs%isd = g%isd ; cs%ied = g%ied ; cs%jsd = g%jsd ; cs%jed = g%jed
221 cs%iscB = g%iscB ; cs%iecB = g%iecB; cs%jscB = g%jscB ; cs%jecB = g%jecB
224 cs%num_col = 0 ; cs%fldno = 0
225 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
226 if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) &
227 cs%num_col = cs%num_col + 1
230 if (cs%num_col > 0)
then 231 allocate(cs%Iresttime_col(cs%num_col)) ; cs%Iresttime_col = 0.0
232 allocate(cs%col_i(cs%num_col)) ; cs%col_i = 0
233 allocate(cs%col_j(cs%num_col)) ; cs%col_j = 0
236 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
237 if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0))
then 238 cs%col_i(col) = i ; cs%col_j(col) = j
239 cs%Iresttime_col(col) = iresttime(i,j)
245 allocate(cs%Ref_h%p(cs%nz_data,cs%num_col))
246 do col=1,cs%num_col ;
do k=1,cs%nz_data
247 cs%Ref_h%p(k,col) = data_h(cs%col_i(col),cs%col_j(col),k)
251 total_sponge_cols = cs%num_col
252 call sum_across_pes(total_sponge_cols)
255 call initialize_remapping(cs%remap_cs, remapscheme, boundary_extrapolation=bndextrapolation, &
256 answers_2018=cs%remap_answers_2018)
258 call log_param(param_file, mdl,
"!Total sponge columns at h points", total_sponge_cols, &
259 "The total number of columns where sponges are applied at h points.", like_default=.true.)
261 if (cs%sponge_uv)
then 263 allocate(data_hu(g%isdB:g%iedB,g%jsd:g%jed,nz_data)) ; data_hu(:,:,:) = 0.0
264 allocate(data_hv(g%isd:g%ied,g%jsdB:g%jedB,nz_data)) ; data_hv(:,:,:) = 0.0
265 allocate(iresttime_u(g%isdB:g%iedB,g%jsd:g%jed)) ; iresttime_u(:,:) = 0.0
266 allocate(iresttime_v(g%isd:g%ied,g%jsdB:g%jedB)) ; iresttime_v(:,:) = 0.0
270 do j=cs%jsc,cs%jec ;
do i=cs%iscB,cs%iecB
271 data_hu(i,j,:) = 0.5 * (data_h(i,j,:) + data_h(i+1,j,:))
272 iresttime_u(i,j) = 0.5 * (iresttime(i,j) + iresttime(i+1,j))
273 if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0)) &
274 cs%num_col_u = cs%num_col_u + 1
277 if (cs%num_col_u > 0)
then 279 allocate(cs%Iresttime_col_u(cs%num_col_u)) ; cs%Iresttime_col_u(:) = 0.0
280 allocate(cs%col_i_u(cs%num_col_u)) ; cs%col_i_u(:) = 0
281 allocate(cs%col_j_u(cs%num_col_u)) ; cs%col_j_u(:) = 0
285 do j=cs%jsc,cs%jec ;
do i=cs%iscB,cs%iecB
286 if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0))
then 287 cs%col_i_u(col) = i ; cs%col_j_u(col) = j
288 cs%Iresttime_col_u(col) = iresttime_u(i,j)
295 allocate(cs%Ref_hu%p(cs%nz_data,cs%num_col_u))
296 do col=1,cs%num_col_u ;
do k=1,cs%nz_data
297 cs%Ref_hu%p(k,col) = data_hu(cs%col_i_u(col),cs%col_j_u(col),k)
300 total_sponge_cols_u = cs%num_col_u
301 call sum_across_pes(total_sponge_cols_u)
302 call log_param(param_file, mdl,
"!Total sponge columns at u points", total_sponge_cols_u, &
303 "The total number of columns where sponges are applied at u points.", like_default=.true.)
307 do j=cs%jscB,cs%jecB;
do i=cs%isc,cs%iec
308 data_hv(i,j,:) = 0.5 * (data_h(i,j,:) + data_h(i,j+1,:))
309 iresttime_v(i,j) = 0.5 * (iresttime(i,j) + iresttime(i,j+1))
310 if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0)) &
311 cs%num_col_v = cs%num_col_v + 1
314 if (cs%num_col_v > 0)
then 316 allocate(cs%Iresttime_col_v(cs%num_col_v)) ; cs%Iresttime_col_v = 0.0
317 allocate(cs%col_i_v(cs%num_col_v)) ; cs%col_i_v = 0
318 allocate(cs%col_j_v(cs%num_col_v)) ; cs%col_j_v = 0
322 do j=cs%jscB,cs%jecB ;
do i=cs%isc,cs%iec
323 if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0))
then 324 cs%col_i_v(col) = i ; cs%col_j_v(col) = j
325 cs%Iresttime_col_v(col) = iresttime_v(i,j)
331 allocate(cs%Ref_hv%p(cs%nz_data,cs%num_col_v))
332 do col=1,cs%num_col_v ;
do k=1,cs%nz_data
333 cs%Ref_hv%p(k,col) = data_hv(cs%col_i_v(col),cs%col_j_v(col),k)
336 total_sponge_cols_v = cs%num_col_v
337 call sum_across_pes(total_sponge_cols_v)
338 call log_param(param_file, mdl,
"!Total sponge columns at v points", total_sponge_cols_v, &
339 "The total number of columns where sponges are applied at v points.", like_default=.true.)