diff --git a/components/mpas-ocean/bld/build-namelist b/components/mpas-ocean/bld/build-namelist
index c7030bd5fc5a..2dfa6a9d1e92 100755
--- a/components/mpas-ocean/bld/build-namelist
+++ b/components/mpas-ocean/bld/build-namelist
@@ -560,6 +560,8 @@ add_default($nl, 'config_Redi_use_surface_taper');
add_default($nl, 'config_Redi_N2_based_taper_enable');
add_default($nl, 'config_Redi_N2_based_taper_min');
add_default($nl, 'config_Redi_N2_based_taper_limit_term1');
+add_default($nl, 'config_Redi_use_quasi_monotone_limiter');
+add_default($nl, 'config_Redi_quasi_monotone_safety_factor');
add_default($nl, 'config_Redi_min_layers_diag_terms');
add_default($nl, 'config_Redi_horizontal_taper');
add_default($nl, 'config_Redi_horizontal_ramp_min');
diff --git a/components/mpas-ocean/bld/build-namelist-section b/components/mpas-ocean/bld/build-namelist-section
index a47c8cc74656..98cfe66feea2 100644
--- a/components/mpas-ocean/bld/build-namelist-section
+++ b/components/mpas-ocean/bld/build-namelist-section
@@ -99,6 +99,8 @@ add_default($nl, 'config_Redi_use_surface_taper');
add_default($nl, 'config_Redi_N2_based_taper_enable');
add_default($nl, 'config_Redi_N2_based_taper_min');
add_default($nl, 'config_Redi_N2_based_taper_limit_term1');
+add_default($nl, 'config_Redi_use_quasi_monotone_limiter');
+add_default($nl, 'config_Redi_quasi_monotone_safety_factor');
add_default($nl, 'config_Redi_min_layers_diag_terms');
add_default($nl, 'config_Redi_horizontal_taper');
add_default($nl, 'config_Redi_horizontal_ramp_min');
diff --git a/components/mpas-ocean/bld/namelist_files/namelist_defaults_mpaso.xml b/components/mpas-ocean/bld/namelist_files/namelist_defaults_mpaso.xml
index f05e111ccd40..56a31a696922 100644
--- a/components/mpas-ocean/bld/namelist_files/namelist_defaults_mpaso.xml
+++ b/components/mpas-ocean/bld/namelist_files/namelist_defaults_mpaso.xml
@@ -149,6 +149,8 @@
.true.
0.1
.true.
+.true.
+0.9
6
15
'ramp'
diff --git a/components/mpas-ocean/bld/namelist_files/namelist_definition_mpaso.xml b/components/mpas-ocean/bld/namelist_files/namelist_definition_mpaso.xml
index 352303431ec7..f3da8bfe4c2d 100644
--- a/components/mpas-ocean/bld/namelist_files/namelist_definition_mpaso.xml
+++ b/components/mpas-ocean/bld/namelist_files/namelist_definition_mpaso.xml
@@ -447,6 +447,22 @@ Valid values: .true. or .false.
Default: Defined in namelist_defaults.xml
+
+If true, fluxes are reduced to prevent tracers from violating monotonicity. Cross-term fluxes are scaled toward zero to prevent tracers from under/overshooting the min/max values in adjacent cells and layers
+
+Valid values: .true. or .false.
+Default: Defined in namelist_defaults.xml
+
+
+
+A safety factor applied to flux scaling when monotonicity is violated. Smaller values scale fluxes toward zero more aggressively.
+
+Valid values: A value between 0 and 1
+Default: Defined in namelist_defaults.xml
+
+
Redi diagonal terms (2 and 3) are turned off from layer 1 through config_Redi_min_layers_diag_terms-1, and on from config_Redi_min_layers_diag_terms to nVertLevels. The Redi diagonal terms are not guaranteed to produce bounded tracer fields, and in practice produce growing temperatures in a few columns with fewer than 5 vertical cells. Redi is meant for isopycnal mixing in the deep ocean, so not applying Redi diagonal terms in very shallow regions is an acceptable solution.
diff --git a/components/mpas-ocean/src/Registry.xml b/components/mpas-ocean/src/Registry.xml
index 0c454fbca420..291a31c55a7c 100644
--- a/components/mpas-ocean/src/Registry.xml
+++ b/components/mpas-ocean/src/Registry.xml
@@ -326,6 +326,14 @@
description="If true, the N2 limiting is applied to the horizontal diffusion term"
possible_values=".true. or .false."
/>
+
+
maximumBnd(iTr, k, iCell))) then
- do i = 1, nEdgesOnCell(iCell)
- iEdge = edgesOnCell(i, iCell)
- redi_term2_edge(iTr, k, iEdge) = 0.0_RKIND
- end do
+ ! scale cross-fluxes down relative to size of bnds violation
+ over = max(tempTracer - maximumBnd(iTr, k, iCell), &
+ minimumBnd(iTr, k, iCell) - tempTracer)
+
+ diff = abs(tempTracer - tracers(iTr, k, iCell))
+
+ scal = min(1.0_RKIND, max(0.0_RKIND, &
+ limiter_safety_factor * ( &
+ 1.0_RKIND - over / (diff + 1.0E-16_RKIND))))
+
+ ! scale to zero if not using monotone limiter
+ if (.not. use_quasi_monotone) scal = 0.0_RKIND
+
+ ! cell flux scalings, to be interp. onto edge/levs next pass
+ redi_flux_scale(iTr, k, iCell) = scal
- redi_term3_topOfCell(iTr, k, iCell) = 0.0_RKIND
- redi_term3_topOfCell(iTr, k + 1, iCell) = 0.0_RKIND
if (isActiveTracer) then
rediLimiterCount(k, iCell) = rediLimiterCount(k, iCell) + 1
- endif
- endif
+ end if
+ end if
end do
end do
end do ! iCell
@@ -422,7 +499,8 @@ subroutine ocn_tracer_hmix_Redi_tend(meshPool, layerThickness, layerThickEdgeMea
!now go back and reapply all tendencies
!$omp parallel
- !$omp do schedule(runtime) private(invAreaCell, k, iTr, i, iEdge)
+ !$omp do schedule(runtime) &
+ !$omp private(invAreaCell, k, kUp, kDn, iTr, i, iEdge, jCell, scal, scalUp, scalDn)
do iCell = 1, nCells
invAreaCell = 1.0_RKIND/areaCell(iCell)
do k = minLevelCell(iCell), maxLevelCell(iCell)
@@ -432,21 +510,48 @@ subroutine ocn_tracer_hmix_Redi_tend(meshPool, layerThickness, layerThickEdgeMea
end do
if (maxLevelCell(iCell) .ge. config_Redi_min_layers_diag_terms) then
+ ! term 3: d/dz div( S grad psi )
do k = minLevelCell(iCell), maxLevelCell(iCell)
+ kUp = max(minLevelCell(iCell), k - 1)
+ kDn = min(maxLevelCell(iCell), k + 1)
do iTr = 1, ntracers
- tend(iTr, k, iCell) = tend(iTr, k, iCell) + rediKappaCell(iCell)*2.0_RKIND* &
- (RediKappaScaling(k, iCell)* &
- redi_term3_topOfCell(iTr, k, iCell) - &
- RediKappaScaling(k + 1, iCell)*redi_term3_topOfCell(iTr, k + 1, iCell))
+ ! harmonic-mean; not strictly monotone, but smooth
+ scalUp = 2.0_RKIND * &
+ redi_flux_scale(iTr, kUp, iCell) * &
+ redi_flux_scale(iTr, k, iCell) / ( &
+ redi_flux_scale(iTr, kUp, iCell) + &
+ redi_flux_scale(iTr, k, iCell) + 1.0E-16_RKIND)
+
+ scalDn = 2.0_RKIND * &
+ redi_flux_scale(iTr, kDn, iCell) * &
+ redi_flux_scale(iTr, k, iCell) / ( &
+ redi_flux_scale(iTr, kDn, iCell) + &
+ redi_flux_scale(iTr, k, iCell) + 1.0E-16_RKIND)
+
+ tend(iTr, k, iCell) = tend(iTr, k, iCell) + &
+ rediKappaCell(iCell) * ( &
+ RediKappaScaling(k, iCell)* &
+ scalUp*redi_term3_topOfCell(iTr, k, iCell) - &
+ RediKappaScaling(k + 1, iCell)* &
+ scalDn*redi_term3_topOfCell(iTr, k + 1, iCell))
end do
end do
+ ! term 2: div( h S d/dz psi )
do i = 1, nEdgesOnCell(iCell)
if (cellsOnCell(i, iCell) .eq. nCellsP1) cycle
iEdge = edgesOnCell(i, iCell)
+ jCell = cellsOnCell(i, iCell)
do k = minLevelEdgeBot(iEdge), maxLevelEdgeTop(iEdge)
do iTr = 1, ntracers
- tend(iTr, k, iCell) = tend(iTr, k, iCell) + edgeSignOnCell(i, iCell)*invAreaCell*redi_term2_edge(iTr, k, iEdge)
+ scal = 2.0_RKIND * &
+ redi_flux_scale(iTr, k, iCell) * &
+ redi_flux_scale(iTr, k, jCell) / ( &
+ redi_flux_scale(iTr, k, iCell) + &
+ redi_flux_scale(iTr, k, jCell) + 1.0E-16_RKIND)
+
+ tend(iTr, k, iCell) = tend(iTr, k, iCell) + &
+ scal*edgeSignOnCell(i, iCell)*invAreaCell*redi_term2_edge(iTr, k, iEdge)
end do
end do
end do
@@ -463,6 +568,9 @@ subroutine ocn_tracer_hmix_Redi_tend(meshPool, layerThickness, layerThickEdgeMea
deallocate (redi_term3)
deallocate (redi_term2_edge)
deallocate (redi_term3_topOfCell)
+ deallocate (minimumBnd)
+ deallocate (maximumBnd)
+ deallocate (redi_flux_scale)
call mpas_timer_stop("tracer redi")
@@ -545,6 +653,12 @@ subroutine ocn_tracer_hmix_Redi_init(domain, err)!{{{
call mpas_dmpar_finalize(domain%dminfo)
end if
+ ! quasi bounds-preserving limiter
+ use_quasi_monotone = config_Redi_use_quasi_monotone_limiter
+
+ ! safety factor for limiter: 0 <= fac <= 1
+ limiter_safety_factor = config_Redi_quasi_monotone_safety_factor
+
! Initialize horizontal taper
if (config_Redi_horizontal_taper == 'none' .or. &
config_Redi_horizontal_taper == 'RossbyRadius') then