24 use mom_time_manager,
only : time_type, init_external_field, get_external_field_size, time_interp_external_init
29 implicit none ;
private
31 #include <MOM_memory.h>
35 module procedure set_up_ale_sponge_field_fixed
36 module procedure set_up_ale_sponge_field_varying
41 module procedure set_up_ale_sponge_vel_field_fixed
42 module procedure set_up_ale_sponge_vel_field_varying
50 module procedure initialize_ale_sponge_fixed
51 module procedure initialize_ale_sponge_varying
56 public get_ale_sponge_thicknesses, get_ale_sponge_nz_data
58 public rotate_ale_sponge, update_ale_sponge_field
70 real,
dimension(:,:,:),
pointer :: mask_in => null()
71 real,
dimension(:,:,:),
pointer :: p => null()
72 real,
dimension(:,:,:),
pointer :: h => null()
80 real,
dimension(:,:),
pointer :: mask_in => null()
81 real,
dimension(:,:),
pointer :: p => null()
82 real,
dimension(:,:),
pointer :: h => null()
107 integer,
pointer :: col_i(:) => null()
108 integer,
pointer :: col_j(:) => null()
109 integer,
pointer :: col_i_u(:) => null()
110 integer,
pointer :: col_j_u(:) => null()
111 integer,
pointer :: col_i_v(:) => null()
112 integer,
pointer :: col_j_v(:) => null()
114 real,
pointer :: iresttime_col(:) => null()
115 real,
pointer :: iresttime_col_u(:) => null()
116 real,
pointer :: iresttime_col_v(:) => null()
118 type(
p3d) :: var(max_fields_)
119 type(
p2d) :: ref_val(max_fields_)
132 logical :: remap_answers_2018
135 logical :: hor_regrid_answers_2018
139 logical :: time_varying_sponges
140 logical :: spongedataongrid
148 subroutine initialize_ale_sponge_fixed(Iresttime, G, param_file, CS, data_h, nz_data)
151 integer,
intent(in) :: nz_data
152 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: Iresttime
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.")
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.)
342 end subroutine initialize_ale_sponge_fixed
346 function get_ale_sponge_nz_data(CS)
349 integer :: get_ale_sponge_nz_data
351 if (
associated(cs))
then
352 get_ale_sponge_nz_data = cs%nz_data
354 get_ale_sponge_nz_data = 0
356 end function get_ale_sponge_nz_data
359 subroutine get_ale_sponge_thicknesses(G, data_h, sponge_mask, CS)
361 real,
allocatable,
dimension(:,:,:), &
362 intent(inout) :: data_h
363 logical,
dimension(SZI_(G),SZJ_(G)), &
364 intent(out) :: sponge_mask
368 integer :: c, i, j, k
370 if (
allocated(data_h))
call mom_error(fatal, &
371 "get_ALE_sponge_thicknesses called with an allocated data_h.")
373 if (.not.
associated(cs))
then
375 allocate(data_h(g%isd:g%ied,g%jsd:g%jed,1)) ; data_h(:,:,:) = -1.0
376 sponge_mask(:,:) = .false.
380 allocate(data_h(g%isd:g%ied,g%jsd:g%jed,cs%nz_data)) ; data_h(:,:,:) = -1.0
381 sponge_mask(:,:) = .false.
384 i = cs%col_i(c) ; j = cs%col_j(c)
385 sponge_mask(i,j) = .true.
387 data_h(i,j,k) = cs%Ref_h%p(k,c)
391 end subroutine get_ale_sponge_thicknesses
396 subroutine initialize_ale_sponge_varying(Iresttime, G, param_file, CS)
399 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: Iresttime
406 #include "version_variable.h"
407 character(len=40) :: mdl =
"MOM_sponge"
408 logical :: use_sponge
409 real,
allocatable,
dimension(:,:) :: Iresttime_u
410 real,
allocatable,
dimension(:,:) :: Iresttime_v
411 logical :: bndExtrapolation = .true.
412 logical :: default_2018_answers
413 logical :: spongeDataOngrid = .false.
414 integer :: i, j, k, col, total_sponge_cols, total_sponge_cols_u, total_sponge_cols_v
415 character(len=10) :: remapScheme
417 if (
associated(cs))
then
418 call mom_error(warning,
"initialize_ALE_sponge_varying called with an associated "// &
419 "control structure.")
424 call get_param(param_file, mdl,
"SPONGE", use_sponge, &
425 "If true, sponges may be applied anywhere in the domain. "//&
426 "The exact location and properties of those sponges are "//&
427 "specified from MOM_initialization.F90.", default=.false.)
428 if (.not.use_sponge)
return
430 call get_param(param_file, mdl,
"SPONGE_UV", cs%sponge_uv, &
431 "Apply sponges in u and v, in addition to tracers.", &
433 call get_param(param_file, mdl,
"REMAPPING_SCHEME", remapscheme, &
434 "This sets the reconstruction scheme used "//&
435 " for vertical remapping for all variables.", &
436 default=
"PLM", do_not_log=.true.)
437 call get_param(param_file, mdl,
"BOUNDARY_EXTRAPOLATION", bndextrapolation, &
438 "When defined, a proper high-order reconstruction "//&
439 "scheme is used within boundary cells rather "//&
440 "than PCM. E.g., if PPM is used for remapping, a "//&
441 "PPM reconstruction will also be used within boundary cells.", &
442 default=.false., do_not_log=.true.)
443 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
444 "This sets the default value for the various _2018_ANSWERS parameters.", &
446 call get_param(param_file, mdl,
"REMAPPING_2018_ANSWERS", cs%remap_answers_2018, &
447 "If true, use the order of arithmetic and expressions that recover the "//&
448 "answers from the end of 2018. Otherwise, use updated and more robust "//&
449 "forms of the same expressions.", default=default_2018_answers)
450 call get_param(param_file, mdl,
"SPONGE_DATA_ONGRID", cs%spongeDataOngrid, &
451 "When defined, the incoming sponge data are "//&
452 "assumed to be on the model grid " , &
454 cs%time_varying_sponges = .true.
456 cs%isc = g%isc ; cs%iec = g%iec ; cs%jsc = g%jsc ; cs%jec = g%jec
457 cs%isd = g%isd ; cs%ied = g%ied ; cs%jsd = g%jsd ; cs%jed = g%jed
458 cs%iscB = g%iscB ; cs%iecB = g%iecB; cs%jscB = g%jscB ; cs%jecB = g%jecB
461 cs%num_col = 0 ; cs%fldno = 0
462 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
463 if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) &
464 cs%num_col = cs%num_col + 1
466 if (cs%num_col > 0)
then
467 allocate(cs%Iresttime_col(cs%num_col)) ; cs%Iresttime_col = 0.0
468 allocate(cs%col_i(cs%num_col)) ; cs%col_i = 0
469 allocate(cs%col_j(cs%num_col)) ; cs%col_j = 0
472 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
473 if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0))
then
474 cs%col_i(col) = i ; cs%col_j(col) = j
475 cs%Iresttime_col(col) = iresttime(i,j)
480 total_sponge_cols = cs%num_col
481 call sum_across_pes(total_sponge_cols)
484 call initialize_remapping(cs%remap_cs, remapscheme, boundary_extrapolation=bndextrapolation, &
485 answers_2018=cs%remap_answers_2018)
486 call log_param(param_file, mdl,
"!Total sponge columns at h points", total_sponge_cols, &
487 "The total number of columns where sponges are applied at h points.", like_default=.true.)
488 if (cs%sponge_uv)
then
489 allocate(iresttime_u(g%isdB:g%iedB,g%jsd:g%jed)) ; iresttime_u(:,:) = 0.0
490 allocate(iresttime_v(g%isd:g%ied,g%jsdB:g%jedB)) ; iresttime_v(:,:) = 0.0
493 do j=cs%jsc,cs%jec;
do i=cs%iscB,cs%iecB
494 iresttime_u(i,j) = 0.5 * (iresttime(i,j) + iresttime(i+1,j))
495 if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0)) &
496 cs%num_col_u = cs%num_col_u + 1
498 if (cs%num_col_u > 0)
then
499 allocate(cs%Iresttime_col_u(cs%num_col_u)) ; cs%Iresttime_col_u = 0.0
500 allocate(cs%col_i_u(cs%num_col_u)) ; cs%col_i_u = 0
501 allocate(cs%col_j_u(cs%num_col_u)) ; cs%col_j_u = 0
504 do j=cs%jsc,cs%jec ;
do i=cs%iscB,cs%iecB
505 if ((iresttime_u(i,j)>0.0) .and. (g%mask2dCu(i,j)>0))
then
506 cs%col_i_u(col) = i ; cs%col_j_u(col) = j
507 cs%Iresttime_col_u(col) = iresttime_u(i,j)
513 total_sponge_cols_u = cs%num_col_u
514 call sum_across_pes(total_sponge_cols_u)
515 call log_param(param_file, mdl,
"!Total sponge columns at u points", total_sponge_cols_u, &
516 "The total number of columns where sponges are applied at u points.", like_default=.true.)
519 do j=cs%jscB,cs%jecB;
do i=cs%isc,cs%iec
520 iresttime_v(i,j) = 0.5 * (iresttime(i,j) + iresttime(i,j+1))
521 if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0)) &
522 cs%num_col_v = cs%num_col_v + 1
524 if (cs%num_col_v > 0)
then
525 allocate(cs%Iresttime_col_v(cs%num_col_v)) ; cs%Iresttime_col_v = 0.0
526 allocate(cs%col_i_v(cs%num_col_v)) ; cs%col_i_v = 0
527 allocate(cs%col_j_v(cs%num_col_v)) ; cs%col_j_v = 0
530 do j=cs%jscB,cs%jecB ;
do i=cs%isc,cs%iec
531 if ((iresttime_v(i,j)>0.0) .and. (g%mask2dCv(i,j)>0))
then
532 cs%col_i_v(col) = i ; cs%col_j_v(col) = j
533 cs%Iresttime_col_v(col) = iresttime_v(i,j)
538 total_sponge_cols_v = cs%num_col_v
539 call sum_across_pes(total_sponge_cols_v)
540 call log_param(param_file, mdl,
"!Total sponge columns at v points", total_sponge_cols_v, &
541 "The total number of columns where sponges are applied at v points.", like_default=.true.)
544 end subroutine initialize_ale_sponge_varying
548 subroutine init_ale_sponge_diags(Time, G, diag, CS)
549 type(time_type),
target,
intent(in) :: time
551 type(
diag_ctrl),
target,
intent(inout) :: diag
555 if (.not.
associated(cs))
return
559 end subroutine init_ale_sponge_diags
563 subroutine set_up_ale_sponge_field_fixed(sp_val, G, f_ptr, CS)
566 real,
dimension(SZI_(G),SZJ_(G),CS%nz_data), &
568 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
569 target,
intent(in) :: f_ptr
572 character(len=256) :: mesg
574 if (.not.
associated(cs))
return
576 cs%fldno = cs%fldno + 1
578 if (cs%fldno > max_fields_)
then
579 write(mesg,
'("Increase MAX_FIELDS_ to at least ",I3," in MOM_memory.h or decrease &
580 &the number of fields to be damped in the call to &
581 &initialize_ALE_sponge." )') cs%fldno
582 call mom_error(fatal,
"set_up_ALE_sponge_field: "//mesg)
586 allocate(cs%Ref_val(cs%fldno)%p(cs%nz_data,cs%num_col))
587 cs%Ref_val(cs%fldno)%p(:,:) = 0.0
590 cs%Ref_val(cs%fldno)%p(k,col) = sp_val(cs%col_i(col),cs%col_j(col),k)
594 cs%var(cs%fldno)%p => f_ptr
596 end subroutine set_up_ale_sponge_field_fixed
600 subroutine set_up_ale_sponge_field_varying(filename, fieldname, Time, G, GV, US, f_ptr, CS)
601 character(len=*),
intent(in) :: filename
603 character(len=*),
intent(in) :: fieldname
605 type(time_type),
intent(in) :: Time
609 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
610 target,
intent(in) :: f_ptr
614 real,
allocatable,
dimension(:,:,:) :: sp_val
615 real,
allocatable,
dimension(:,:,:) :: mask_z
616 real,
allocatable,
dimension(:),
target :: z_in, z_edges_in
617 real :: missing_value
619 integer :: isd,ied,jsd,jed
621 integer,
dimension(4) :: fld_sz
623 character(len=256) :: mesg
625 real,
dimension(:),
allocatable :: tmpT1d
626 real :: zTopOfCell, zBottomOfCell
629 if (.not.
associated(cs))
return
631 call time_interp_external_init()
632 isd = g%isd; ied = g%ied; jsd = g%jsd; jed = g%jed
633 cs%fldno = cs%fldno + 1
634 if (cs%fldno > max_fields_)
then
635 write(mesg,
'("Increase MAX_FIELDS_ to at least ",I3," in MOM_memory.h or decrease &
636 &the number of fields to be damped in the call to &
637 &initialize_ALE_sponge." )') cs%fldno
638 call mom_error(fatal,
"set_up_ALE_sponge_field: "//mesg)
642 if (cs%spongeDataOngrid)
then
643 cs%Ref_val(cs%fldno)%id = init_external_field(filename, fieldname,domain=g%Domain%mpp_domain)
645 cs%Ref_val(cs%fldno)%id = init_external_field(filename, fieldname)
648 fld_sz = get_external_field_size(cs%Ref_val(cs%fldno)%id)
650 cs%Ref_val(cs%fldno)%nz_data = nz_data
651 cs%Ref_val(cs%fldno)%num_tlevs = fld_sz(4)
654 allocate(cs%Ref_val(cs%fldno)%p(nz_data,cs%num_col))
655 cs%Ref_val(cs%fldno)%p(:,:) = 0.0
656 allocate( cs%Ref_val(cs%fldno)%h(nz_data,cs%num_col) )
657 cs%Ref_val(cs%fldno)%h(:,:) = 0.0
658 cs%var(cs%fldno)%p => f_ptr
661 end subroutine set_up_ale_sponge_field_varying
665 subroutine set_up_ale_sponge_vel_field_fixed(u_val, v_val, G, u_ptr, v_ptr, CS)
668 real,
dimension(SZIB_(G),SZJ_(G),CS%nz_data), &
670 real,
dimension(SZI_(G),SZJB_(G),CS%nz_data), &
672 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
target,
intent(in) :: u_ptr
673 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
target,
intent(in) :: v_ptr
676 character(len=256) :: mesg
678 if (.not.
associated(cs))
return
681 allocate(cs%Ref_val_u%p(cs%nz_data,cs%num_col_u))
682 cs%Ref_val_u%p(:,:) = 0.0
683 do col=1,cs%num_col_u
685 cs%Ref_val_u%p(k,col) = u_val(cs%col_i_u(col),cs%col_j_u(col),k)
689 allocate(cs%Ref_val_v%p(cs%nz_data,cs%num_col_v))
690 cs%Ref_val_v%p(:,:) = 0.0
691 do col=1,cs%num_col_v
693 cs%Ref_val_v%p(k,col) = v_val(cs%col_i_v(col),cs%col_j_v(col),k)
698 end subroutine set_up_ale_sponge_vel_field_fixed
702 subroutine set_up_ale_sponge_vel_field_varying(filename_u, fieldname_u, filename_v, fieldname_v, &
703 Time, G, US, CS, u_ptr, v_ptr)
704 character(len=*),
intent(in) :: filename_u
705 character(len=*),
intent(in) :: fieldname_u
706 character(len=*),
intent(in) :: filename_v
707 character(len=*),
intent(in) :: fieldname_v
708 type(time_type),
intent(in) :: Time
712 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
target,
intent(in) :: u_ptr
713 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
target,
intent(in) :: v_ptr
715 real,
allocatable,
dimension(:,:,:) :: u_val
716 real,
allocatable,
dimension(:,:,:) :: mask_u
717 real,
allocatable,
dimension(:,:,:) :: v_val
718 real,
allocatable,
dimension(:,:,:) :: mask_v
720 real,
allocatable,
dimension(:),
target :: z_in, z_edges_in
721 real :: missing_value
724 integer :: isd, ied, jsd, jed
725 integer :: isdB, iedB, jsdB, jedB
726 integer,
dimension(4) :: fld_sz
727 character(len=256) :: mesg
729 if (.not.
associated(cs))
return
731 isd = g%isd; ied = g%ied; jsd = g%jsd; jed = g%jed
732 isdb = g%isdB; iedb = g%iedB; jsdb = g%jsdB; jedb = g%jedB
736 cs%Ref_val_u%id = init_external_field(filename_u, fieldname_u)
738 fld_sz = get_external_field_size(cs%Ref_val_u%id)
739 cs%Ref_val_u%nz_data = fld_sz(3)
740 cs%Ref_val_u%num_tlevs = fld_sz(4)
741 cs%Ref_val_v%id = init_external_field(filename_v, fieldname_v)
743 fld_sz = get_external_field_size(cs%Ref_val_v%id)
744 cs%Ref_val_v%nz_data = fld_sz(3)
745 cs%Ref_val_v%num_tlevs = fld_sz(4)
746 allocate( u_val(isdb:iedb,jsd:jed, fld_sz(3)) )
747 allocate( mask_u(isdb:iedb,jsd:jed, fld_sz(3)) )
748 allocate( v_val(isd:ied,jsdb:jedb, fld_sz(3)) )
749 allocate( mask_v(isd:ied,jsdb:jedb, fld_sz(3)) )
755 missing_value, .true., .false., .false., m_to_z=us%m_to_Z, &
756 answers_2018=cs%hor_regrid_answers_2018)
763 missing_value, .true., .false., .false., m_to_z=us%m_to_Z, &
764 answers_2018=cs%hor_regrid_answers_2018)
766 allocate(cs%Ref_val_u%p(fld_sz(3),cs%num_col_u))
767 cs%Ref_val_u%p(:,:) = 0.0
768 do col=1,cs%num_col_u
770 cs%Ref_val_u%p(k,col) = u_val(cs%col_i_u(col),cs%col_j_u(col),k)
774 allocate(cs%Ref_val_v%p(fld_sz(3),cs%num_col_v))
775 cs%Ref_val_v%p(:,:) = 0.0
776 do col=1,cs%num_col_v
778 cs%Ref_val_v%p(k,col) = v_val(cs%col_i_v(col),cs%col_j_v(col),k)
783 end subroutine set_up_ale_sponge_vel_field_varying
787 subroutine apply_ale_sponge(h, dt, G, GV, US, CS, Time)
791 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
793 real,
intent(in) :: dt
796 type(time_type),
optional,
intent(in) :: time
801 real,
allocatable,
dimension(:) :: tmp_val2
802 real,
dimension(SZK_(G)) :: tmp_val1
803 real :: hu(szib_(g), szj_(g), szk_(g))
804 real :: hv(szi_(g), szjb_(g), szk_(g))
805 real,
allocatable,
dimension(:,:,:) :: sp_val
806 real,
allocatable,
dimension(:,:,:) :: mask_z
807 real,
dimension(:),
allocatable :: hsrc
809 real,
dimension(:),
allocatable :: tmpt1d
810 integer :: c, m, nkmb, i, j, k, is, ie, js, je, nz, nz_data
811 integer :: col, total_sponge_cols
812 real,
allocatable,
dimension(:),
target :: z_in, z_edges_in
813 real :: missing_value
814 real :: h_neglect, h_neglect_edge
815 real :: ztopofcell, zbottomofcell
818 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
819 if (.not.
associated(cs))
return
821 if (.not.cs%remap_answers_2018)
then
822 h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
823 elseif (gv%Boussinesq)
then
824 h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
826 h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
829 if (cs%time_varying_sponges)
then
830 if (.not.
present(time)) &
831 call mom_error(fatal,
"apply_ALE_sponge: No time information provided")
833 nz_data = cs%Ref_val(m)%nz_data
834 allocate(sp_val(g%isd:g%ied,g%jsd:g%jed,1:nz_data))
835 allocate(mask_z(g%isd:g%ied,g%jsd:g%jed,1:nz_data))
839 z_edges_in, missing_value, .true., .false., .false., &
840 spongeongrid=cs%SpongeDataOngrid, m_to_z=us%m_to_Z, &
841 answers_2018=cs%hor_regrid_answers_2018)
842 allocate( hsrc(nz_data) )
843 allocate( tmpt1d(nz_data) )
845 i = cs%col_i(c) ; j = cs%col_j(c)
846 cs%Ref_val(m)%p(1:nz_data,c) = sp_val(i,j,1:nz_data)
848 ztopofcell = 0. ; zbottomofcell = 0. ; npoints = 0; hsrc(:) = 0.0; tmpt1d(:) = -99.9
850 if (mask_z(cs%col_i(c),cs%col_j(c),k) == 1.0)
then
851 zbottomofcell = -min( z_edges_in(k+1), g%bathyT(cs%col_i(c),cs%col_j(c)) )
852 tmpt1d(k) = sp_val(cs%col_i(c),cs%col_j(c),k)
854 zbottomofcell = -g%bathyT(cs%col_i(c),cs%col_j(c))
855 tmpt1d(k) = tmpt1d(k-1)
859 hsrc(k) = ztopofcell - zbottomofcell
860 if (hsrc(k)>0.) npoints = npoints + 1
861 ztopofcell = zbottomofcell
864 hsrc(nz_data) = hsrc(nz_data) + ( ztopofcell + g%bathyT(cs%col_i(c),cs%col_j(c)) )
865 cs%Ref_val(m)%h(1:nz_data,c) = gv%Z_to_H*hsrc(1:nz_data)
866 cs%Ref_val(m)%p(1:nz_data,c) = tmpt1d(1:nz_data)
869 if (cs%Ref_val(m)%h(k,c) <= 0.001*gv%m_to_H) &
872 cs%Ref_val(m)%p(k,c) = cs%Ref_val(m)%p(k-1,c)
875 deallocate(sp_val, mask_z, hsrc, tmpt1d)
881 allocate(tmp_val2(nz_data))
886 i = cs%col_i(c) ; j = cs%col_j(c)
887 damp = dt * cs%Iresttime_col(c)
888 i1pdamp = 1.0 / (1.0 + damp)
889 tmp_val2(1:nz_data) = cs%Ref_val(m)%p(1:nz_data,c)
890 if (cs%time_varying_sponges)
then
891 call remapping_core_h(cs%remap_cs, nz_data, cs%Ref_val(m)%h(1:nz_data,c), tmp_val2, &
892 cs%nz, h(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
894 call remapping_core_h(cs%remap_cs,nz_data, cs%Ref_h%p(1:nz_data,c), tmp_val2, &
895 cs%nz, h(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
898 cs%var(m)%p(i,j,1:cs%nz) = i1pdamp * (cs%var(m)%p(i,j,1:cs%nz) + tmp_val1 * damp)
909 if (cs%sponge_uv)
then
911 do j=cs%jsc,cs%jec;
do i=cs%iscB,cs%iecB;
do k=1,nz
912 hu(i,j,k) = 0.5 * (h(i,j,k) + h(i+1,j,k))
913 enddo ;
enddo ;
enddo
914 if (cs%time_varying_sponges)
then
915 if (.not.
present(time)) &
916 call mom_error(fatal,
"apply_ALE_sponge: No time information provided")
918 nz_data = cs%Ref_val_u%nz_data
919 allocate(sp_val(g%isdB:g%iedB,g%jsd:g%jed,1:nz_data))
920 allocate(mask_z(g%isdB:g%iedB,g%jsd:g%jed,1:nz_data))
923 z_edges_in, missing_value, .true., .false., .false., &
924 m_to_z=us%m_to_Z, answers_2018=cs%hor_regrid_answers_2018)
930 i = cs%col_i(c) ; j = cs%col_j(c)
931 cs%Ref_val_u%p(1:nz_data,c) = sp_val(i,j,1:nz_data)
934 deallocate (sp_val, mask_z)
936 nz_data = cs%Ref_val_v%nz_data
937 allocate(sp_val(g%isd:g%ied,g%jsdB:g%jedB,1:nz_data))
938 allocate(mask_z(g%isd:g%ied,g%jsdB:g%jedB,1:nz_data))
941 z_edges_in, missing_value, .true., .false., .false., &
942 m_to_z=us%m_to_Z, answers_2018=cs%hor_regrid_answers_2018)
950 i = cs%col_i(c) ; j = cs%col_j(c)
951 cs%Ref_val_v%p(1:nz_data,c) = sp_val(i,j,1:nz_data)
954 deallocate (sp_val, mask_z)
961 i = cs%col_i_u(c) ; j = cs%col_j_u(c)
962 damp = dt * cs%Iresttime_col_u(c)
963 i1pdamp = 1.0 / (1.0 + damp)
964 if (cs%time_varying_sponges) nz_data = cs%Ref_val(m)%nz_data
965 tmp_val2(1:nz_data) = cs%Ref_val_u%p(1:nz_data,c)
966 if (cs%time_varying_sponges)
then
967 call remapping_core_h(cs%remap_cs, nz_data, cs%Ref_val_u%h(:,c), tmp_val2, &
968 cs%nz, hu(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
970 call remapping_core_h(cs%remap_cs, nz_data, cs%Ref_hu%p(:,c), tmp_val2, &
971 cs%nz, hu(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
974 cs%var_u%p(i,j,:) = i1pdamp * (cs%var_u%p(i,j,:) + tmp_val1 * damp)
978 do j=cs%jscB,cs%jecB;
do i=cs%isc,cs%iec;
do k=1,nz
979 hv(i,j,k) = 0.5 * (h(i,j,k) + h(i,j+1,k))
980 enddo ;
enddo ;
enddo
983 i = cs%col_i_v(c) ; j = cs%col_j_v(c)
984 damp = dt * cs%Iresttime_col_v(c)
985 i1pdamp = 1.0 / (1.0 + damp)
986 tmp_val2(1:nz_data) = cs%Ref_val_v%p(1:nz_data,c)
987 if (cs%time_varying_sponges)
then
988 call remapping_core_h(cs%remap_cs, cs%nz_data, cs%Ref_val_v%h(:,c), tmp_val2, &
989 cs%nz, hv(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
991 call remapping_core_h(cs%remap_cs, cs%nz_data, cs%Ref_hv%p(:,c), tmp_val2, &
992 cs%nz, hv(i,j,:), tmp_val1, h_neglect, h_neglect_edge)
995 cs%var_v%p(i,j,:) = i1pdamp * (cs%var_v%p(i,j,:) + tmp_val1 * damp)
1000 deallocate(tmp_val2)
1002 end subroutine apply_ale_sponge
1005 subroutine rotate_ale_sponge(sponge_in, G_in, sponge, G, turns, param_file)
1012 integer,
intent(in) :: turns
1022 real,
dimension(:,:),
allocatable :: iresttime_in, iresttime
1023 real,
dimension(:,:,:),
allocatable :: data_h_in, data_h
1024 real,
dimension(:,:,:),
allocatable :: sp_val_in, sp_val
1025 real,
dimension(:,:,:),
pointer :: sp_ptr => null()
1026 integer :: c, c_i, c_j
1027 integer :: k, nz_data
1029 logical :: fixed_sponge
1031 fixed_sponge = .not. sponge_in%time_varying_sponges
1034 allocate(iresttime_in(g_in%isd:g_in%ied, g_in%jsd:g_in%jed))
1035 allocate(iresttime(g%isd:g%ied, g%jsd:g%jed))
1036 iresttime_in(:,:) = 0.0
1038 if (fixed_sponge)
then
1039 nz_data = sponge_in%nz_data
1040 allocate(data_h_in(g_in%isd:g_in%ied, g_in%jsd:g_in%jed, nz_data))
1041 allocate(data_h(g%isd:g%ied, g%jsd:g%jed, nz_data))
1042 data_h_in(:,:,:) = 0.
1046 do c=1,sponge_in%num_col
1047 c_i = sponge_in%col_i(c)
1048 c_j = sponge_in%col_j(c)
1049 iresttime_in(c_i, c_j) = sponge_in%Iresttime_col(c)
1050 if (fixed_sponge)
then
1052 data_h(c_i, c_j, k) = sponge_in%Ref_h%p(k,c)
1058 if (fixed_sponge)
then
1060 call initialize_ale_sponge_fixed(iresttime, g, param_file, sponge, &
1063 call initialize_ale_sponge_varying(iresttime, g, param_file, sponge)
1066 deallocate(iresttime_in)
1067 deallocate(iresttime)
1068 if (fixed_sponge)
then
1069 deallocate(data_h_in)
1075 sponge%fldno = sponge_in%fldno
1077 if (fixed_sponge)
then
1078 allocate(sp_val_in(g_in%isd:g_in%ied, g_in%jsd:g_in%jed, nz_data))
1079 allocate(sp_val(g%isd:g%ied, g%jsd:g%jed, nz_data))
1082 do n=1,sponge_in%fldno
1084 sp_ptr => sponge_in%var(n)%p
1085 if (fixed_sponge)
then
1086 sp_val_in(:,:,:) = 0.0
1087 do c = 1, sponge_in%num_col
1088 c_i = sponge_in%col_i(c)
1089 c_j = sponge_in%col_j(c)
1091 sp_val_in(c_i, c_j, k) = sponge_in%Ref_val(n)%p(k,c)
1100 deallocate(sp_val_in)
1105 sponge%Ref_val(n)%id = sponge_in%Ref_val(n)%id
1106 sponge%Ref_val(n)%num_tlevs = sponge_in%Ref_val(n)%num_tlevs
1108 nz_data = sponge_in%Ref_val(n)%nz_data
1109 sponge%Ref_val(n)%nz_data = nz_data
1111 allocate(sponge%Ref_val(n)%p(nz_data, sponge_in%num_col))
1112 allocate(sponge%Ref_val(n)%h(nz_data, sponge_in%num_col))
1113 sponge%Ref_val(n)%p(:,:) = 0.0
1114 sponge%Ref_val(n)%h(:,:) = 0.0
1126 sponge%var(n)%p => sp_ptr
1131 if (
associated(sponge_in%var_u%p) .or.
associated(sponge_in%var_v%p)) &
1132 call mom_error(fatal,
"Rotation of ALE sponge velocities is not yet " &
1136 sponge%diag => sponge_in%diag
1139 end subroutine rotate_ale_sponge
1146 subroutine update_ale_sponge_field(sponge, p_old, G, GV, p_new)
1149 real,
dimension(:,:,:), &
1150 target,
intent(in) :: p_old
1153 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1154 target,
intent(in) :: p_new
1159 if (
associated(sponge%var(n)%p, p_old)) sponge%var(n)%p => p_new
1162 end subroutine update_ale_sponge_field
1168 subroutine ale_sponge_end(CS)
1174 if (.not.
associated(cs))
return
1176 if (
associated(cs%col_i))
deallocate(cs%col_i)
1177 if (
associated(cs%col_i_u))
deallocate(cs%col_i_u)
1178 if (
associated(cs%col_i_v))
deallocate(cs%col_i_v)
1179 if (
associated(cs%col_j))
deallocate(cs%col_j)
1180 if (
associated(cs%col_j_u))
deallocate(cs%col_j_u)
1181 if (
associated(cs%col_j_v))
deallocate(cs%col_j_v)
1183 if (
associated(cs%Iresttime_col))
deallocate(cs%Iresttime_col)
1184 if (
associated(cs%Iresttime_col_u))
deallocate(cs%Iresttime_col_u)
1185 if (
associated(cs%Iresttime_col_v))
deallocate(cs%Iresttime_col_v)
1188 if (
associated(cs%Ref_val(m)%p))
deallocate(cs%Ref_val(m)%p)
1193 end subroutine ale_sponge_end