Skip to content

Commit

Permalink
+*Use thickness in ideal_age_tracer_column_physics
Browse files Browse the repository at this point in the history
  Pass mixed layer thickness, rather than mixed layer depth, to
ideal_age_tracer_column_physics, and determine the number of layers within
the mixed layer in count_BL_layers in thickness units rather than depth units.
To accommodate these changes there is now a call to convert_MLD_to_ML_thickness
inside of call_tracer_column_fns.  The thermo_var_ptrs type argument (tv) to
ideal_age_tracer_column_physics and count_BL_layers are no longer used and
have been removed.  All answers are bitwise identical in Boussinesq mode, but
in non-Boussinesq mode there are changes in the passive ideal age tracers at
the level of roundoff.  There are also changes in the units of arguments and
the number of arguments to a public interface.
  • Loading branch information
Hallberg-NOAA committed Apr 25, 2024
1 parent 6536d42 commit e1d9244
Show file tree
Hide file tree
Showing 2 changed files with 29 additions and 27 deletions.
27 changes: 17 additions & 10 deletions src/tracer/MOM_tracer_flow_control.F90
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,7 @@ module MOM_tracer_flow_control
use MOM_get_input, only : Get_MOM_input
use MOM_grid, only : ocean_grid_type
use MOM_hor_index, only : hor_index_type
use MOM_interface_heights, only : convert_MLD_to_ML_thickness
use MOM_CVMix_KPP, only : KPP_CS
use MOM_open_boundary, only : ocean_OBC_type
use MOM_restart, only : MOM_restart_CS
Expand Down Expand Up @@ -427,7 +428,7 @@ subroutine call_tracer_set_forcing(sfc_state, fluxes, day_start, day_interval, G
end subroutine call_tracer_set_forcing

!> This subroutine calls all registered tracer column physics subroutines.
subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV, US, tv, optics, CS, &
subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, mld, dt, G, GV, US, tv, optics, CS, &
debug, KPP_CSp, nonLocalTrans, evap_CFL_limit, minimum_forcing_depth)
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
Expand All @@ -444,7 +445,7 @@ subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV,
type(forcing), intent(in) :: fluxes !< A structure containing pointers to
!! any possible forcing fields.
!! Unused fields have NULL ptrs.
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: Hml !< Mixed layer depth [Z ~> m]
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: mld !< Mixed layer depth [Z ~> m]
real, intent(in) :: dt !< The amount of time covered by this
!! call [T ~> s]
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
Expand All @@ -464,6 +465,9 @@ subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV,
real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over
!! which fluxes can be applied [H ~> m or kg m-2]

! Local variables
real :: Hbl(SZI_(G),SZJ_(G)) !< Boundary layer thickness [H ~> m or kg m-2]

if (.not. associated(CS)) call MOM_error(FATAL, "call_tracer_column_fns: "// &
"Module must be initialized via call_tracer_register before it is used.")

Expand All @@ -488,12 +492,13 @@ subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV,
G, GV, US, CS%RGC_tracer_CSp, &
evap_CFL_limit=evap_CFL_limit, &
minimum_forcing_depth=minimum_forcing_depth)
if (CS%use_ideal_age) &
if (CS%use_ideal_age) then
call convert_MLD_to_ML_thickness(mld, h_new, Hbl, tv, G, GV)
call ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, &
G, GV, US, tv, CS%ideal_age_tracer_CSp, &
G, GV, US, CS%ideal_age_tracer_CSp, &
evap_CFL_limit=evap_CFL_limit, &
minimum_forcing_depth=minimum_forcing_depth, &
Hbl=Hml)
minimum_forcing_depth=minimum_forcing_depth, Hbl=Hbl)
endif
if (CS%use_regional_dyes) &
call dye_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, &
G, GV, US, tv, CS%dye_tracer_CSp, &
Expand Down Expand Up @@ -526,7 +531,7 @@ subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV,
if (US%QRZ_T_to_W_m2 /= 1.0) call MOM_error(FATAL, "MOM_generic_tracer_column_physics "//&
"has not been written to permit dimensionsal rescaling. Set all 4 of the "//&
"[QRZT]_RESCALE_POWER parameters to 0.")
call MOM_generic_tracer_column_physics(h_old, h_new, ea, eb, fluxes, Hml, dt, &
call MOM_generic_tracer_column_physics(h_old, h_new, ea, eb, fluxes, mld, dt, &
G, GV, US, CS%MOM_generic_tracer_CSp, tv, optics, &
evap_CFL_limit=evap_CFL_limit, &
minimum_forcing_depth=minimum_forcing_depth)
Expand Down Expand Up @@ -567,9 +572,11 @@ subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV,
if (CS%use_RGC_tracer) &
call RGC_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, &
G, GV, US, CS%RGC_tracer_CSp)
if (CS%use_ideal_age) &
if (CS%use_ideal_age) then
call convert_MLD_to_ML_thickness(mld, h_new, Hbl, tv, G, GV)
call ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, &
G, GV, US, tv, CS%ideal_age_tracer_CSp, Hbl=Hml)
G, GV, US, CS%ideal_age_tracer_CSp, Hbl=Hbl)
endif
if (CS%use_regional_dyes) &
call dye_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, &
G, GV, US, tv, CS%dye_tracer_CSp)
Expand All @@ -591,7 +598,7 @@ subroutine call_tracer_column_fns(h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV,
if (US%QRZ_T_to_W_m2 /= 1.0) call MOM_error(FATAL, "MOM_generic_tracer_column_physics "//&
"has not been written to permit dimensionsal rescaling. Set all 4 of the "//&
"[QRZT]_RESCALE_POWER parameters to 0.")
call MOM_generic_tracer_column_physics(h_old, h_new, ea, eb, fluxes, Hml, dt, &
call MOM_generic_tracer_column_physics(h_old, h_new, ea, eb, fluxes, mld, dt, &
G, GV, US, CS%MOM_generic_tracer_CSp, tv, optics)
endif
if (CS%use_pseudo_salt_tracer) &
Expand Down
29 changes: 12 additions & 17 deletions src/tracer/ideal_age_example.F90
Original file line number Diff line number Diff line change
Expand Up @@ -12,7 +12,6 @@ module ideal_age_example
use MOM_grid, only : ocean_grid_type
use MOM_hor_index, only : hor_index_type
use MOM_io, only : file_exists, MOM_read_data, slasher, vardesc, var_desc, query_vardesc
use MOM_interface_heights, only : thickness_to_dz
use MOM_open_boundary, only : ocean_OBC_type
use MOM_restart, only : query_initialized, set_initialized, MOM_restart_CS
use MOM_spatial_means, only : global_mass_int_EFP
Expand All @@ -22,7 +21,7 @@ module ideal_age_example
use MOM_tracer_diabatic, only : tracer_vertdiff, applyTracerBoundaryFluxesInOut
use MOM_tracer_Z_init, only : tracer_Z_init
use MOM_unit_scaling, only : unit_scale_type
use MOM_variables, only : surface, thermo_var_ptrs
use MOM_variables, only : surface
use MOM_verticalGrid, only : verticalGrid_type

implicit none ; private
Expand Down Expand Up @@ -297,7 +296,7 @@ subroutine initialize_ideal_age_tracer(restart, day, G, GV, US, h, diag, OBC, CS
end subroutine initialize_ideal_age_tracer

!> Applies diapycnal diffusion, aging and regeneration at the surface to the ideal age tracers
subroutine ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, US, tv, CS, &
subroutine ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, &
evap_CFL_limit, minimum_forcing_depth, Hbl)
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
Expand All @@ -317,14 +316,13 @@ subroutine ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G,
!! and tracer forcing fields. Unused fields have NULL ptrs.
real, intent(in) :: dt !< The amount of time covered by this call [T ~> s]
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
type(ideal_age_tracer_CS), pointer :: CS !< The control structure returned by a previous
!! call to register_ideal_age_tracer.
real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can
!! be fluxed out of the top layer in a timestep [nondim]
real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
!! fluxes can be applied [H ~> m or kg m-2]
real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: Hbl !< Boundary layer depth [Z ~> m]
real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: Hbl !< Boundary layer thickness [H ~> m or kg m-2]

! This subroutine applies diapycnal diffusion and any other column
! tracer physics or chemistry to the tracers from this file.
Expand All @@ -349,7 +347,7 @@ subroutine ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G,
endif

if (CS%use_real_BL_depth .and. present(Hbl)) then
call count_BL_layers(G, GV, h_old, Hbl, tv, BL_layers)
call count_BL_layers(G, GV, h_old, Hbl, BL_layers)
endif

if (.not.associated(CS)) return
Expand Down Expand Up @@ -578,30 +576,27 @@ subroutine ideal_age_example_end(CS)
endif
end subroutine ideal_age_example_end

subroutine count_BL_layers(G, GV, h, Hbl, tv, BL_layers)
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
subroutine count_BL_layers(G, GV, h, Hbl, BL_layers)
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: Hbl !< Boundary layer depth [Z ~> m]
type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: Hbl !< Boundary layer thickness [H ~> m or kg m-2]
real, dimension(SZI_(G),SZJ_(G)), intent(out) :: BL_layers !< Number of model layers in the boundary layer [nondim]

real :: dz(SZI_(G),SZK_(GV)) ! Height change across layers [Z ~> m]
real :: current_depth ! Distance from the free surface [Z ~> m]
real :: current_depth ! Distance from the free surface [H ~> m or kg m-2]
integer :: i, j, k, is, ie, js, je, nz, m, nk
character(len=255) :: msg
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec ; nz = GV%ke

BL_layers(:,:) = 0.
do j=js,je
call thickness_to_dz(h, tv, dz, j, G, GV)
do i=is,ie
current_depth = 0.
do k=1,nz
current_depth = current_depth + dz(i,k)
current_depth = current_depth + h(i,j,k)
if (Hbl(i,j) <= current_depth) then
BL_layers(i,j) = BL_layers(i,j) + (1.0 - (current_depth - Hbl(i,j)) / dz(i,k))
BL_layers(i,j) = BL_layers(i,j) + (1.0 - (current_depth - Hbl(i,j)) / h(i,j,k))
exit
else
BL_layers(i,j) = BL_layers(i,j) + 1.0
Expand Down

0 comments on commit e1d9244

Please sign in to comment.