Skip to content

Commit f111427

Browse files
Theresa MorrisonTheresa Morrison
Theresa Morrison
authored and
Theresa Morrison
committed
Add Initialization for Generic Tracer Sponges
Add a routine that will initialize the pointers for ALE sponges to be used with any of the generic tracers. The data for the nuding is set by the user with a number of new parameters.
1 parent d475de1 commit f111427

File tree

1 file changed

+147
-0
lines changed

1 file changed

+147
-0
lines changed

src/tracer/MOM_generic_tracer.F90

+147
Original file line numberDiff line numberDiff line change
@@ -30,6 +30,7 @@ module MOM_generic_tracer
3030
use g_tracer_utils, only: g_tracer_get_obc_segment_props
3131

3232
use MOM_ALE_sponge, only : set_up_ALE_sponge_field, ALE_sponge_CS
33+
use MOM_ALE_sponge, only : ALE_sponge_CS, initialize_ALE_sponge
3334
use MOM_coms, only : EFP_type, max_across_PEs, min_across_PEs, PE_here
3435
use MOM_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
3536
use MOM_diag_mediator, only : diag_ctrl, get_diag_time_end
@@ -308,6 +309,7 @@ subroutine initialize_MOM_generic_tracer(restart, day, G, GV, US, h, tv, param_f
308309

309310
character(len=128), parameter :: sub_name = 'initialize_MOM_generic_tracer'
310311
logical :: OK,obc_has
312+
logical :: do_use_GT_sponge
311313
integer :: i, j, k, isc, iec, jsc, jec, nk
312314
type(g_tracer_type), pointer :: g_tracer,g_tracer_next
313315
character(len=fm_string_len) :: g_tracer_name
@@ -438,6 +440,10 @@ subroutine initialize_MOM_generic_tracer(restart, day, G, GV, US, h, tv, param_f
438440
call g_tracer_set_common(G%isc,G%iec,G%jsc,G%jec,G%isd,G%ied,G%jsd,G%jed,&
439441
GV%ke,1,CS%diag%axesTL%handles,grid_tmask,grid_kmt,day)
440442

443+
call get_param(param_file, "initialize_sponges_file", "DO_SPONGE_GENERIC_TRACER", do_use_gt_sponge, &
444+
"If true, then some generic tracers may be nudged.", default=.false.)
445+
if (do_use_GT_sponge) call g_tracer_initialize_sponges(G, GV, US, CS, param_file, sponge_CSp, ALE_sponge_CSp, day)
446+
441447
! Register generic tracer modules diagnostics
442448

443449
#ifdef _USE_MOM6_DIAG
@@ -450,6 +456,147 @@ subroutine initialize_MOM_generic_tracer(restart, day, G, GV, US, h, tv, param_f
450456

451457
end subroutine initialize_MOM_generic_tracer
452458

459+
subroutine g_tracer_initialize_sponges(G, GV, US, CS, param_file, Layer_CSp, ALE_CSp, Time)
460+
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
461+
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
462+
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
463+
real, dimension(SZI_(G),SZJ_(G)) :: depth_tot !< The nominal total depth of the ocean [Z ~> m]
464+
type(MOM_generic_tracer_CS), pointer :: CS !< Pointer to the control structure for this module.
465+
type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters.
466+
type(sponge_CS), pointer :: Layer_CSp !< A pointer that is set to point to the control
467+
!! structure for this module (in layered mode).
468+
type(ALE_sponge_CS), pointer :: ALE_CSp !< A pointer that is set to point to the control
469+
!! structure for this module (in ALE mode).
470+
type(time_type), intent(in) :: Time !< Time at the start of the run segment. Time_in
471+
!! overrides any value set for Time.
472+
! Local variables
473+
real, allocatable, dimension(:,:,:) :: eta ! The target interface heights [Z ~> m].
474+
real, allocatable, dimension(:,:,:) :: dz ! The target interface thicknesses in height units [Z ~> m]
475+
real, allocatable, dimension(:,:,:) :: h ! The target interface thicknesses [H ~> m or kg m-2].
476+
477+
real, allocatable, dimension(:,:,:) :: tmp_GT ! A temporary array for reading sponge target temperatures
478+
! on the vertical grid of the input file [C ~> degC]
479+
480+
real :: Idamp(SZI_(G),SZJ_(G)) ! The sponge damping rate [T-1 ~> s-1]
481+
482+
integer :: i, j, k, is, ie, js, je, nz
483+
integer :: isd, ied, jsd, jed
484+
integer, dimension(4) :: siz
485+
integer :: nz_data ! The size of the sponge source grid
486+
character(len=40) :: tmp_var, Idamp_var, eta_var
487+
character(len=40) :: mdl = "initialize_sponges_file"
488+
character(len=200) :: damping_file, state_file, state_uv_file ! Strings for filenames
489+
character(len=200) :: filename, inputdir ! Strings for file/path and path.
490+
491+
character(len=200) :: g_tracer_name
492+
type(g_tracer_type), pointer :: g_tracer,g_tracer_next
493+
logical :: do_sponge_gt
494+
real, dimension(:,:,:), pointer :: g_tracer_ptr
495+
496+
logical :: use_ALE ! True if ALE is being used, False if in layered mode
497+
logical :: time_space_interp_sponge ! If true use sponge data that need to be interpolated in both
498+
! the horizontal dimension and in time prior to vertical remapping.
499+
500+
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec ; nz = GV%ke
501+
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
502+
503+
Idamp(:,:) = 0.0
504+
505+
call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
506+
inputdir = slasher(inputdir)
507+
call get_param(param_file, mdl, "SPONGE_GT_DAMPING_FILE", damping_file, &
508+
"The name of the file with the sponge damping rates.") !, &
509+
call get_param(param_file, mdl, "SPONGE_IDAMP_VAR", Idamp_var, &
510+
"The name of the inverse damping rate variable in "//&
511+
"SPONGE_DAMPING_FILE.", default="Idamp")
512+
call get_param(param_file, mdl, "USE_REGRIDDING", use_ALE, default=.false., do_not_log=.true.)
513+
call get_param(param_file, mdl, "INTERPOLATE_SPONGE_TIME_SPACE", time_space_interp_sponge, &
514+
"If True, perform on-the-fly regridding in lat-lon-time of sponge restoring data.", &
515+
default=.false.)
516+
517+
! Read in sponge damping rate for tracers
518+
filename = trim(inputdir)//trim(damping_file)
519+
call log_param(param_file, mdl, "INPUTDIR/SPONGE_DAMPING_FILE", filename)
520+
if (.not.file_exists(filename, G%Domain)) &
521+
call MOM_error(FATAL, " initialize_sponges: Unable to open "//trim(filename))
522+
523+
call MOM_read_data(filename, Idamp_var, Idamp(:,:), G%Domain, scale=US%T_to_s)
524+
525+
! Now register all of the fields which are damped in the sponge.
526+
! By default, momentum is advected vertically within the sponge, but
527+
! momentum is typically not damped within the sponge.
528+
529+
g_tracer=>CS%g_tracer_list
530+
do ! check for each generic tracer if it is nudged
531+
532+
call g_tracer_get_alias(g_tracer,g_tracer_name)
533+
534+
call get_param(param_file, mdl, "DO_SPONGE_GT_"//trim(g_tracer_name), do_sponge_gt, &
535+
"If true, then generic tracer "//trim(g_tracer_name)//" is nudged.", &
536+
default=.false.)
537+
538+
if (do_sponge_gt) then
539+
call get_param(param_file, mdl, "SPONGE_GT_"//trim(g_tracer_name)//"_FILE", state_file, &
540+
"The name of the file with the state to damp the generic ", &
541+
"tracer "//trim(g_tracer_name)//" toward.", default="sponge_"//trim(g_tracer_name)//".nc.")
542+
call get_param(param_file, mdl, "SPONGE_GT_"//trim(g_tracer_name)//"_VAR", tmp_var, &
543+
"The name of the variable to use in the sponge_GT file for generic ", &
544+
"tracer "//trim(g_tracer_name)//".", default=trim(g_tracer_name))
545+
call get_param(param_file, mdl, "SPONGE_GT_"//trim(g_tracer_name)//"_ETA_VAR", eta_var, &
546+
"The name of the interface height variable in "//&
547+
"SPONGE_GT_"//trim(g_tracer_name)//"_FILE.", default="ETA")
548+
549+
filename = trim(inputdir)//trim(state_file)
550+
call log_param(param_file, mdl, "INPUTDIR/SPONGE_GT"//trim(g_tracer_name)//"_FILE", filename)
551+
if (.not.file_exists(filename, G%Domain)) &
552+
call MOM_error(FATAL, " initialize_sponges: Unable to open "//trim(filename))
553+
554+
! get the pointer for this tracer
555+
call g_tracer_get_pointer(g_tracer,g_tracer_name,'field',g_tracer_ptr)
556+
557+
if (use_ALE) then ! ALE mode
558+
if (.not. time_space_interp_sponge) then
559+
!call field_size(filename,eta_var,siz,no_domain=.true.)
560+
if (siz(1) /= G%ieg-G%isg+1 .or. siz(2) /= G%jeg-G%jsg+1) &
561+
call MOM_error(FATAL,"initialize_sponge_file: Array size mismatch for sponge data.")
562+
nz_data = siz(3)-1
563+
allocate(eta(isd:ied,jsd:jed,nz_data+1))
564+
allocate(dz(isd:ied,jsd:jed,nz_data))
565+
call MOM_read_data(filename, eta_var, eta(:,:,:), G%Domain, scale=US%m_to_Z)
566+
do j=js,je ; do i=is,ie
567+
eta(i,j,nz_data+1) = -depth_tot(i,j)
568+
enddo ; enddo
569+
do k=nz_data,1,-1 ; do j=js,je ; do i=is,ie
570+
if (eta(i,j,K) < (eta(i,j,K+1) + GV%Angstrom_Z)) &
571+
eta(i,j,K) = eta(i,j,K+1) + GV%Angstrom_Z
572+
enddo ; enddo ; enddo
573+
do k=1,nz_data ; do j=js,je ; do i=is,ie
574+
dz(i,j,k) = eta(i,j,k)-eta(i,j,k+1)
575+
enddo; enddo ; enddo
576+
deallocate(eta)
577+
578+
allocate(tmp_GT(isd:ied,jsd:jed,nz_data))
579+
call MOM_read_data(filename, tmp_var, tmp_GT(:,:,:), G%Domain) !, scale=US%degC_to_C)
580+
581+
call set_up_ALE_sponge_field(tmp_GT, G, GV, g_tracer_ptr, &
582+
ALE_CSp, trim(g_tracer_name))
583+
deallocate(tmp_GT)
584+
deallocate(dz)
585+
586+
else
587+
call set_up_ALE_sponge_field(filename, tmp_var, Time, G, GV, US, g_tracer_ptr, &
588+
ALE_CSp, trim(g_tracer_name))
589+
endif
590+
endif
591+
endif
592+
593+
call g_tracer_get_next(g_tracer, g_tracer_next)
594+
if (.NOT. associated(g_tracer_next)) exit
595+
g_tracer=>g_tracer_next
596+
enddo
597+
598+
end subroutine g_tracer_initialize_sponges
599+
453600
!> Column physics for generic tracers.
454601
!! Get the coupler values for generic tracers that exchange with atmosphere
455602
!! Update generic tracer concentration fields from sources and sinks.

0 commit comments

Comments
 (0)