MOM6
mom_entrain_diffusive Module Reference

Detailed Description

Diapycnal mixing and advection in isopycnal mode.

By Robert Hallberg, September 1997 - July 2000

This file contains the subroutines that implement diapycnal mixing and advection in isopycnal layers. The main subroutine, calculate_entrainment, returns the entrainment by each layer across the interfaces above and below it. These are calculated subject to the constraints that no layers can be driven to neg- ative thickness and that the each layer maintains its target density, using the scheme described in Hallberg (MWR 2000). There may or may not be a bulk mixed layer above the isopycnal layers. The solution is iterated until the change in the entrainment between successive iterations is less than some small tolerance.

The dual-stream entrainment scheme of MacDougall and Dewar (JPO 1997) is used for combined diapycnal advection and diffusion, modified as described in Hallberg (MWR 2000) to be solved implicitly in time. Any profile of diffusivities may be used. Diapycnal advection is fundamentally the residual of diapycnal diffusion, so the fully implicit upwind differencing scheme that is used is entirely appropriate. The downward buoyancy flux in each layer is determined from an implicit calculation based on the previously calculated flux of the layer above and an estim- ated flux in the layer below. This flux is subject to the foll- owing conditions: (1) the flux in the top and bottom layers are set by the boundary conditions, and (2) no layer may be driven below an Angstrom thickness. If there is a bulk mixed layer, the mixed and buffer layers are treated as Eulerian layers, whose thicknesses only change due to entrainment by the interior layers.

Data Types

type  entrain_diffusive_cs
 The control structure holding parametes for the MOM_entrain_diffusive module. More...
 

Functions/Subroutines

subroutine, public entrainment_diffusive (h, tv, fluxes, dt, G, GV, US, CS, ea, eb, kb_out, Kd_Lay, Kd_int)
 This subroutine calculates ea and eb, the rates at which a layer entrains from the layers above and below. The entrainment rates are proportional to the buoyancy flux in a layer and inversely proportional to the density differences between layers. The scheme that is used here is described in detail in Hallberg, Mon. Wea. Rev. 2000. More...
 
subroutine f_to_ent (F, h, kb, kmb, j, G, GV, CS, dsp1_ds, eakb, Ent_bl, ea, eb, do_i_in)
 This subroutine calculates the actual entrainments (ea and eb) and the amount of surface forcing that is applied to each layer if there is no bulk mixed layer. More...
 
subroutine set_ent_bl (h, dtKd_int, tv, kb, kmb, do_i, G, GV, US, CS, j, Ent_bl, Sref, h_bl)
 This subroutine sets the average entrainment across each of the interfaces between buffer layers within a timestep. It also causes thin and relatively light interior layers to be entrained by the deepest buffer layer. Also find the initial coordinate potential densities (Sref) of each layer. More...
 
subroutine determine_dskb (h_bl, Sref, Ent_bl, E_kb, is, ie, kmb, G, GV, limit, dSkb, ddSkb_dE, dSlay, ddSlay_dE, dS_anom_lim, do_i_in)
 This subroutine determines the reference density difference between the bottommost buffer layer and the first interior after the mixing between mixed and buffer layers and mixing with the layer below. Within the mixed and buffer layers, entrainment from the layer above is increased when it is necessary to keep the layers from developing a negative thickness; otherwise it equals Ent_bl. At each interface, the upward and downward fluxes average out to Ent_bl, unless entrainment by the layer below is larger than twice Ent_bl. The density difference across the first interior layer may also be returned. It could also be limited to avoid negative values or values that greatly exceed the density differences across an interface. Additionally, the partial derivatives of dSkb and dSlay with E_kb could also be returned. More...
 
subroutine f_kb_to_ea_kb (h_bl, Sref, Ent_bl, I_dSkbp1, F_kb, kmb, i, G, GV, CS, ea_kb, tol_in)
 Given an entrainment from below for layer kb, determine a consistent entrainment from above, such that dSkb * ea_kb = dSkbp1 * F_kb. The input value of ea_kb is both the maximum value that can be obtained and the first guess of the iterations. Ideally ea_kb should be an under-estimate. More...
 
subroutine determine_ea_kb (h_bl, dtKd_kb, Sref, I_dSkbp1, Ent_bl, ea_kbp1, min_eakb, max_eakb, kmb, is, ie, do_i, G, GV, CS, Ent, error, err_min_eakb0, err_max_eakb0, F_kb, dFdfm_kb)
 This subroutine determines the entrainment from above by the top interior layer (labeled kb elsewhere) given an entrainment by the layer below it, constrained to be within the provided bounds. More...
 
subroutine find_maxf_kb (h_bl, Sref, Ent_bl, I_dSkbp1, min_ent_in, max_ent_in, kmb, is, ie, G, GV, CS, maxF, ent_maxF, do_i_in, F_lim_maxent, F_thresh)
 Maximize F = ent*ds_kb*I_dSkbp1 in the range min_ent < ent < max_ent. More...
 
subroutine, public entrain_diffusive_init (Time, G, GV, US, param_file, diag, CS, just_read_params)
 This subroutine initializes the parameters and memory associated with the entrain_diffusive module. More...
 
subroutine, public entrain_diffusive_end (CS)
 This subroutine cleans up and deallocates any memory associated with the entrain_diffusive module. More...
 

Function/Subroutine Documentation

◆ determine_dskb()

subroutine mom_entrain_diffusive::determine_dskb ( real, dimension(szi_(g),szk_(g)), intent(in)  h_bl,
real, dimension(szi_(g),szk_(g)), intent(in)  Sref,
real, dimension(szi_(g),szk_(g)), intent(in)  Ent_bl,
real, dimension(szi_(g)), intent(in)  E_kb,
integer, intent(in)  is,
integer, intent(in)  ie,
integer, intent(in)  kmb,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
logical, intent(in)  limit,
real, dimension(szi_(g)), intent(inout)  dSkb,
real, dimension(szi_(g)), intent(inout), optional  ddSkb_dE,
real, dimension(szi_(g)), intent(inout), optional  dSlay,
real, dimension(szi_(g)), intent(inout), optional  ddSlay_dE,
real, dimension(szi_(g)), intent(inout), optional  dS_anom_lim,
logical, dimension(szi_(g)), intent(in), optional  do_i_in 
)
private

This subroutine determines the reference density difference between the bottommost buffer layer and the first interior after the mixing between mixed and buffer layers and mixing with the layer below. Within the mixed and buffer layers, entrainment from the layer above is increased when it is necessary to keep the layers from developing a negative thickness; otherwise it equals Ent_bl. At each interface, the upward and downward fluxes average out to Ent_bl, unless entrainment by the layer below is larger than twice Ent_bl. The density difference across the first interior layer may also be returned. It could also be limited to avoid negative values or values that greatly exceed the density differences across an interface. Additionally, the partial derivatives of dSkb and dSlay with E_kb could also be returned.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]h_blLayer thickness [H ~> m or kg m-2]
[in]srefReference potential density [R ~> kg m-3]
[in]ent_blThe average entrainment upward and downward across each interface around the buffer layers [H ~> m or kg m-2].
[in]e_kbThe entrainment by the top interior layer [H ~> m or kg m-2].
[in]isThe start of the i-index range to work on.
[in]ieThe end of the i-index range to work on.
[in]kmbThe number of mixed and buffer layers.
[in]limitIf true, limit dSkb and dSlay to avoid negative values.
[in,out]dskbThe limited potential density difference across the interface between the bottommost buffer layer and the topmost interior layer. [R ~> kg m-3] dSkb > 0.
[in,out]ddskb_deThe partial derivative of dSkb with E [R H-1 ~> kg m-4 or m-1].
[in,out]dslayThe limited potential density difference across the topmost interior layer. 0 < dSkb [R ~> kg m-3]
[in,out]ddslay_deThe partial derivative of dSlay with E [R H-1 ~> kg m-4 or m-1].
[in,out]ds_anom_limA limiting value to use for the density anomalies below the buffer layer [R ~> kg m-3].
[in]do_i_inIf present, determines which columns are worked on.

Definition at line 1201 of file MOM_entrain_diffusive.F90.

1201  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1202  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid
1203  !! structure.
1204  real, dimension(SZI_(G),SZK_(G)), intent(in) :: h_bl !< Layer thickness [H ~> m or kg m-2]
1205  real, dimension(SZI_(G),SZK_(G)), intent(in) :: Sref !< Reference potential density [R ~> kg m-3]
1206  real, dimension(SZI_(G),SZK_(G)), intent(in) :: Ent_bl !< The average entrainment upward and
1207  !! downward across each interface
1208  !! around the buffer layers [H ~> m or kg m-2].
1209  real, dimension(SZI_(G)), intent(in) :: E_kb !< The entrainment by the top interior
1210  !! layer [H ~> m or kg m-2].
1211  integer, intent(in) :: is !< The start of the i-index range to work on.
1212  integer, intent(in) :: ie !< The end of the i-index range to work on.
1213  integer, intent(in) :: kmb !< The number of mixed and buffer layers.
1214  logical, intent(in) :: limit !< If true, limit dSkb and dSlay to
1215  !! avoid negative values.
1216  real, dimension(SZI_(G)), intent(inout) :: dSkb !< The limited potential density
1217  !! difference across the interface
1218  !! between the bottommost buffer layer
1219  !! and the topmost interior layer. [R ~> kg m-3]
1220  !! dSkb > 0.
1221  real, dimension(SZI_(G)), optional, intent(inout) :: ddSkb_dE !< The partial derivative of dSkb
1222  !! with E [R H-1 ~> kg m-4 or m-1].
1223  real, dimension(SZI_(G)), optional, intent(inout) :: dSlay !< The limited potential density
1224  !! difference across the topmost
1225  !! interior layer. 0 < dSkb [R ~> kg m-3]
1226  real, dimension(SZI_(G)), optional, intent(inout) :: ddSlay_dE !< The partial derivative of dSlay
1227  !! with E [R H-1 ~> kg m-4 or m-1].
1228  real, dimension(SZI_(G)), optional, intent(inout) :: dS_anom_lim !< A limiting value to use for
1229  !! the density anomalies below the
1230  !! buffer layer [R ~> kg m-3].
1231  logical, dimension(SZI_(G)), optional, intent(in) :: do_i_in !< If present, determines which
1232  !! columns are worked on.
1233 
1234 ! Note that dSkb, ddSkb_dE, dSlay, ddSlay_dE, and dS_anom_lim are declared
1235 ! intent inout because they should not change where do_i_in is false.
1236 
1237 ! This subroutine determines the reference density difference between the
1238 ! bottommost buffer layer and the first interior after the mixing between mixed
1239 ! and buffer layers and mixing with the layer below. Within the mixed and buffer
1240 ! layers, entrainment from the layer above is increased when it is necessary to
1241 ! keep the layers from developing a negative thickness; otherwise it equals
1242 ! Ent_bl. At each interface, the upward and downward fluxes average out to
1243 ! Ent_bl, unless entrainment by the layer below is larger than twice Ent_bl.
1244 ! The density difference across the first interior layer may also be returned.
1245 ! It could also be limited to avoid negative values or values that greatly
1246 ! exceed the density differences across an interface.
1247 ! Additionally, the partial derivatives of dSkb and dSlay with E_kb could
1248 ! also be returned.
1249 
1250  ! Local variables
1251  real, dimension(SZI_(G),SZK_(G)) :: &
1252  b1, c1, & ! b1 and c1 are variables used by the tridiagonal solver.
1253  S, dS_dE, & ! The coordinate density [R ~> kg m-3] and its derivative with E.
1254  ea, dea_dE, & ! The entrainment from above and its derivative with E.
1255  eb, deb_dE ! The entrainment from below and its derivative with E.
1256  real :: deriv_dSkb(SZI_(G))
1257  real :: d1(SZI_(G)) ! d1 = 1.0-c1 is also used by the tridiagonal solver.
1258  real :: src ! A source term for dS_dR.
1259  real :: h1 ! The thickness in excess of the minimum that will remain
1260  ! after exchange with the layer below [H ~> m or kg m-2].
1261  logical, dimension(SZI_(G)) :: do_i
1262  real :: h_neglect ! A thickness that is so small it is usually lost
1263  ! in roundoff and can be neglected [H ~> m or kg m-2].
1264  real :: h_tr ! h_tr is h at tracer points with a tiny thickness
1265  ! added to ensure positive definiteness [H ~> m or kg m-2].
1266  real :: b_denom_1 ! The first term in the denominator of b1 [H ~> m or kg m-2].
1267  real :: rat
1268  real :: dS_kbp1, IdS_kbp1
1269  real :: deriv_dSLay
1270  real :: Inv_term ! [nondim]
1271  real :: f1, df1_drat ! Temporary variables [nondim].
1272  real :: z, dz_drat, f2, df2_dz, expz ! Temporary variables [nondim].
1273  real :: eps_dSLay, eps_dSkb ! Small nondimensional constants.
1274  integer :: i, k
1275 
1276  if (present(ddslay_de) .and. .not.present(dslay)) call mom_error(fatal, &
1277  "In deterimine_dSkb, ddSLay_dE may only be present if dSlay is.")
1278 
1279  h_neglect = gv%H_subroundoff
1280 
1281  do i=is,ie
1282  ea(i,kmb+1) = e_kb(i) ; dea_de(i,kmb+1) = 1.0
1283  s(i,kmb+1) = sref(i,kmb+1) ; ds_de(i,kmb+1) = 0.0
1284  b1(i,kmb+1) = 0.0
1285  d1(i) = 1.0
1286  do_i(i) = .true.
1287  enddo
1288  if (present(do_i_in)) then
1289  do i=is,ie ; do_i(i) = do_i_in(i) ; enddo
1290  endif
1291  do k=kmb,1,-1 ; do i=is,ie
1292  if (do_i(i)) then
1293  ! The do_i test here is only for efficiency.
1294  ! Determine the entrainment from below for each buffer layer.
1295  if (2.0*ent_bl(i,k+1) > ea(i,k+1)) then
1296  eb(i,k) = 2.0*ent_bl(i,k+1) - ea(i,k+1) ; deb_de(i,k) = -dea_de(i,k+1)
1297  else
1298  eb(i,k) = 0.0 ; deb_de(i,k) = 0.0
1299  endif
1300 
1301  ! Determine the entrainment from above for each buffer layer.
1302  h1 = (h_bl(i,k) - gv%Angstrom_H) + (eb(i,k) - ea(i,k+1))
1303  if (h1 >= 0.0) then
1304  ea(i,k) = ent_bl(i,k) ; dea_de(i,k) = 0.0
1305  elseif (ent_bl(i,k) + 0.5*h1 >= 0.0) then
1306  ea(i,k) = ent_bl(i,k) - 0.5*h1
1307  dea_de(i,k) = 0.5*(dea_de(i,k+1) - deb_de(i,k))
1308  else
1309  ea(i,k) = -h1
1310  dea_de(i,k) = dea_de(i,k+1) - deb_de(i,k)
1311  endif
1312  else
1313  ea(i,k) = 0.0 ; dea_de(i,k) = 0.0 ; eb(i,k) = 0.0 ; deb_de(i,k) = 0.0
1314  endif
1315 
1316  ! This is the first-pass of a tridiagonal solver for S.
1317  h_tr = h_bl(i,k) + h_neglect
1318  c1(i,k) = ea(i,k+1) * b1(i,k+1)
1319  b_denom_1 = (h_tr + d1(i)*eb(i,k))
1320  b1(i,k) = 1.0 / (b_denom_1 + ea(i,k))
1321  d1(i) = b_denom_1 * b1(i,k)
1322 
1323  s(i,k) = (h_tr*sref(i,k) + eb(i,k)*s(i,k+1)) * b1(i,k)
1324  enddo ; enddo
1325  do k=2,kmb ; do i=is,ie
1326  s(i,k) = s(i,k) + c1(i,k-1)*s(i,k-1)
1327  enddo ; enddo
1328 
1329  if (present(ddskb_de) .or. present(ddslay_de)) then
1330  ! These two tridiagonal solvers cannot be combined because the solutions for
1331  ! S are required as a source for dS_dE.
1332  do k=kmb,2,-1 ; do i=is,ie
1333  if (do_i(i) .and. (dea_de(i,k) - deb_de(i,k) > 0.0)) then
1334  src = (((s(i,k+1) - sref(i,k)) * (h_bl(i,k) + h_neglect) + &
1335  (s(i,k+1) - s(i,k-1)) * ea(i,k)) * deb_de(i,k) - &
1336  ((sref(i,k) - s(i,k-1)) * h_bl(i,k) + &
1337  (s(i,k+1) - s(i,k-1)) * eb(i,k)) * dea_de(i,k)) / &
1338  ((h_bl(i,k) + h_neglect + ea(i,k)) + eb(i,k))
1339  else ; src = 0.0 ; endif
1340  ds_de(i,k) = (src + eb(i,k)*ds_de(i,k+1)) * b1(i,k)
1341  enddo ; enddo
1342  do i=is,ie
1343  if (do_i(i) .and. (deb_de(i,1) < 0.0)) then
1344  src = (((s(i,2) - sref(i,1)) * (h_bl(i,1) + h_neglect)) * deb_de(i,1)) / &
1345  (h_bl(i,1) + h_neglect + eb(i,1))
1346  else ; src = 0.0 ; endif
1347  ds_de(i,1) = (src + eb(i,1)*ds_de(i,2)) * b1(i,1)
1348  enddo
1349  do k=2,kmb ; do i=is,ie
1350  ds_de(i,k) = ds_de(i,k) + c1(i,k-1)*ds_de(i,k-1)
1351  enddo ; enddo
1352  endif
1353 
1354  ! Now, apply any limiting and return the requested variables.
1355 
1356  eps_dskb = 1.0e-6 ! Should be a small, nondimensional, positive number.
1357  if (.not.limit) then
1358  do i=is,ie ; if (do_i(i)) then
1359  dskb(i) = sref(i,kmb+1) - s(i,kmb)
1360  endif ; enddo
1361  if (present(ddskb_de)) then ; do i=is,ie ; if (do_i(i)) then
1362  ddskb_de(i) = -1.0*ds_de(i,kmb)
1363  endif ; enddo ; endif
1364 
1365  if (present(dslay)) then ; do i=is,ie ; if (do_i(i)) then
1366  dslay(i) = 0.5 * (sref(i,kmb+2) - s(i,kmb))
1367  endif ; enddo ; endif
1368  if (present(ddslay_de)) then ; do i=is,ie ; if (do_i(i)) then
1369  ddslay_de(i) = -0.5*ds_de(i,kmb)
1370  endif ; enddo ; endif
1371  else
1372  do i=is,ie ; if (do_i(i)) then
1373  ! Need to ensure that 0 < dSkb <= S_kb - Sbl
1374  if (sref(i,kmb+1) - s(i,kmb) < eps_dskb*(sref(i,kmb+2) - sref(i,kmb+1))) then
1375  dskb(i) = eps_dskb * (sref(i,kmb+2) - sref(i,kmb+1)) ; deriv_dskb(i) = 0.0
1376  else
1377  dskb(i) = sref(i,kmb+1) - s(i,kmb) ; deriv_dskb(i) = -1.0
1378  endif
1379  if (present(ddskb_de)) ddskb_de(i) = deriv_dskb(i)*ds_de(i,kmb)
1380  endif ; enddo
1381 
1382  if (present(dslay)) then
1383  dz_drat = 1000.0 ! The limit of large dz_drat the same as choosing a
1384  ! Heaviside function.
1385  eps_dslay = 1.0e-10 ! Should be ~= GV%Angstrom_H / sqrt(Kd*dt)
1386  do i=is,ie ; if (do_i(i)) then
1387  ds_kbp1 = sref(i,kmb+2) - sref(i,kmb+1)
1388  ids_kbp1 = 1.0 / (sref(i,kmb+2) - sref(i,kmb+1))
1389  rat = (sref(i,kmb+1) - s(i,kmb)) * ids_kbp1
1390  ! Need to ensure that 0 < dSLay <= 2*dSkb
1391  if (rat < 0.5) then
1392  ! The coefficients here are chosen so that at rat = 0.5, the value (1.5)
1393  ! and first derivative (-0.5) match with the "typical" case (next).
1394  ! The functional form here is arbitrary.
1395  ! f1 provides a reasonable profile that matches the value and derivative
1396  ! of the "typical" case at rat = 0.5, and has a maximum of less than 2.
1397  inv_term = 1.0 / (1.0-rat)
1398  f1 = 2.0 - 0.125*(inv_term**2)
1399  df1_drat = - 0.25*(inv_term**3)
1400 
1401  ! f2 ensures that dSLay goes to 0 rapidly if rat is significantly
1402  ! negative.
1403  z = dz_drat * rat + 4.0 ! The 4 here gives f2(0) = 0.982.
1404  if (z >= 18.0) then ; f2 = 1.0 ; df2_dz = 0.0
1405  elseif (z <= -58.0) then ; f2 = eps_dslay ; df2_dz = 0.0
1406  else
1407  expz = exp(z) ; inv_term = 1.0 / (1.0 + expz)
1408  f2 = (eps_dslay + expz) * inv_term
1409  df2_dz = (1.0 - eps_dslay) * expz * inv_term**2
1410  endif
1411 
1412  dslay(i) = dskb(i) * f1 * f2
1413  deriv_dslay = deriv_dskb(i) * (f1 * f2) - (dskb(i)*ids_kbp1) * &
1414  (df1_drat*f2 + f1 * dz_drat * df2_dz)
1415  elseif (dskb(i) <= 3.0*ds_kbp1) then
1416  ! This is the "typical" case.
1417  dslay(i) = 0.5 * (dskb(i) + ds_kbp1)
1418  deriv_dslay = 0.5 * deriv_dskb(i) ! = -0.5
1419  else
1420  dslay(i) = 2.0*ds_kbp1
1421  deriv_dslay = 0.0
1422  endif
1423  if (present(ddslay_de)) ddslay_de(i) = deriv_dslay*ds_de(i,kmb)
1424  endif ; enddo
1425  endif ! present(dSlay)
1426  endif ! Not limited.
1427 
1428  if (present(ds_anom_lim)) then ; do i=is,ie ; if (do_i(i)) then
1429  ds_anom_lim(i) = max(0.0, eps_dskb * (sref(i,kmb+2) - sref(i,kmb+1)) - &
1430  (sref(i,kmb+1) - s(i,kmb)) )
1431  endif ; enddo ; endif
1432 

◆ determine_ea_kb()

subroutine mom_entrain_diffusive::determine_ea_kb ( real, dimension(szi_(g),szk_(g)), intent(in)  h_bl,
real, dimension(szi_(g)), intent(in)  dtKd_kb,
real, dimension(szi_(g),szk_(g)), intent(in)  Sref,
real, dimension(szi_(g)), intent(in)  I_dSkbp1,
real, dimension(szi_(g),szk_(g)), intent(in)  Ent_bl,
real, dimension(szi_(g)), intent(in)  ea_kbp1,
real, dimension(szi_(g)), intent(in)  min_eakb,
real, dimension(szi_(g)), intent(in)  max_eakb,
integer, intent(in)  kmb,
integer, intent(in)  is,
integer, intent(in)  ie,
logical, dimension(szi_(g)), intent(in)  do_i,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(entrain_diffusive_cs), pointer  CS,
real, dimension( g %isd: g %ied), intent(inout)  Ent,
real, dimension( g %isd: g %ied), intent(out), optional  error,
real, dimension( g %isd: g %ied), intent(in), optional  err_min_eakb0,
real, dimension( g %isd: g %ied), intent(in), optional  err_max_eakb0,
real, dimension( g %isd: g %ied), intent(out), optional  F_kb,
real, dimension( g %isd: g %ied), intent(out), optional  dFdfm_kb 
)
private

This subroutine determines the entrainment from above by the top interior layer (labeled kb elsewhere) given an entrainment by the layer below it, constrained to be within the provided bounds.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]h_blLayer thickness, with the top interior layer at k-index kmb+1 [H ~> m or kg m-2].
[in]srefThe coordinate reference potential density, with the value of the topmost interior layer at layer kmb+1 [R ~> kg m-3].
[in]ent_blThe average entrainment upward and downward across each interface around the buffer layers [H ~> m or kg m-2].
[in]i_dskbp1The inverse of the difference in reference potential density across the base of the uppermost interior layer [R-1 ~> m3 kg-1].
[in]dtkd_kbThe diapycnal diffusivity in the top interior layer times the time step [H2 ~> m2 or kg2 m-4].
[in]ea_kbp1The entrainment from above by layer kb+1 [H ~> m or kg m-2].
[in]min_eakbThe minimum permissible rate of entrainment [H ~> m or kg m-2].
[in]max_eakbThe maximum permissible rate of entrainment [H ~> m or kg m-2].
[in]kmbThe number of mixed and buffer layers.
[in]isThe start of the i-index range to work on.
[in]ieThe end of the i-index range to work on.
[in]do_iA logical variable indicating which i-points to work on.
csThis module's control structure.
[in,out]entThe entrainment rate of the uppermost interior layer [H ~> m or kg m-2]. The input value is the first guess.
[out]errorThe error (locally defined in this routine) associated with the returned solution.
[in]err_min_eakb0The errors (locally defined) associated with min_eakb when ea_kbp1 = 0, returned from a previous call to this fn.
[in]err_max_eakb0The errors (locally defined) associated with min_eakb when ea_kbp1 = 0, returned from a previous call to this fn.
[out]f_kbThe entrainment from below by the uppermost interior layer corresponding to the returned value of Ent [H ~> m or kg m-2].
[out]dfdfm_kbThe partial derivative of F_kb with ea_kbp1 [nondim].

Definition at line 1574 of file MOM_entrain_diffusive.F90.

1574  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1575  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1576  real, dimension(SZI_(G),SZK_(G)), intent(in) :: h_bl !< Layer thickness, with the top interior
1577  !! layer at k-index kmb+1 [H ~> m or kg m-2].
1578  real, dimension(SZI_(G),SZK_(G)), intent(in) :: Sref !< The coordinate reference potential
1579  !! density, with the value of the
1580  !! topmost interior layer at layer
1581  !! kmb+1 [R ~> kg m-3].
1582  real, dimension(SZI_(G),SZK_(G)), intent(in) :: Ent_bl !< The average entrainment upward and
1583  !! downward across each interface around
1584  !! the buffer layers [H ~> m or kg m-2].
1585  real, dimension(SZI_(G)), intent(in) :: I_dSkbp1 !< The inverse of the difference in
1586  !! reference potential density across
1587  !! the base of the uppermost interior
1588  !! layer [R-1 ~> m3 kg-1].
1589  real, dimension(SZI_(G)), intent(in) :: dtKd_kb !< The diapycnal diffusivity in the top
1590  !! interior layer times the time step
1591  !! [H2 ~> m2 or kg2 m-4].
1592  real, dimension(SZI_(G)), intent(in) :: ea_kbp1 !< The entrainment from above by layer
1593  !! kb+1 [H ~> m or kg m-2].
1594  real, dimension(SZI_(G)), intent(in) :: min_eakb !< The minimum permissible rate of
1595  !! entrainment [H ~> m or kg m-2].
1596  real, dimension(SZI_(G)), intent(in) :: max_eakb !< The maximum permissible rate of
1597  !! entrainment [H ~> m or kg m-2].
1598  integer, intent(in) :: kmb !< The number of mixed and buffer layers.
1599  integer, intent(in) :: is !< The start of the i-index range to work on.
1600  integer, intent(in) :: ie !< The end of the i-index range to work on.
1601  logical, dimension(SZI_(G)), intent(in) :: do_i !< A logical variable indicating which
1602  !! i-points to work on.
1603  type(entrain_diffusive_CS), pointer :: CS !< This module's control structure.
1604  real, dimension(SZI_(G)), intent(inout) :: Ent !< The entrainment rate of the uppermost
1605  !! interior layer [H ~> m or kg m-2].
1606  !! The input value is the first guess.
1607  real, dimension(SZI_(G)), optional, intent(out) :: error !< The error (locally defined in this
1608  !! routine) associated with the returned
1609  !! solution.
1610  real, dimension(SZI_(G)), optional, intent(in) :: err_min_eakb0 !< The errors (locally defined)
1611  !! associated with min_eakb when ea_kbp1 = 0,
1612  !! returned from a previous call to this fn.
1613  real, dimension(SZI_(G)), optional, intent(in) :: err_max_eakb0 !< The errors (locally defined)
1614  !! associated with min_eakb when ea_kbp1 = 0,
1615  !! returned from a previous call to this fn.
1616  real, dimension(SZI_(G)), optional, intent(out) :: F_kb !< The entrainment from below by the
1617  !! uppermost interior layer
1618  !! corresponding to the returned
1619  !! value of Ent [H ~> m or kg m-2].
1620  real, dimension(SZI_(G)), optional, intent(out) :: dFdfm_kb !< The partial derivative of F_kb with
1621  !! ea_kbp1 [nondim].
1622 
1623 ! This subroutine determines the entrainment from above by the top interior
1624 ! layer (labeled kb elsewhere) given an entrainment by the layer below it,
1625 ! constrained to be within the provided bounds.
1626 
1627  ! Local variables
1628  real, dimension(SZI_(G)) :: &
1629  dS_kb, & ! The coordinate-density difference between the
1630  ! layer kb and deepest buffer layer, limited to
1631  ! ensure that it is positive [R ~> kg m-3].
1632  ds_lay, & ! The coordinate-density difference across layer
1633  ! kb, limited to ensure that it is positive and not
1634  ! too much bigger than dS_kb or dS_kbp1 [R ~> kg m-3].
1635  ddskb_de, ddslay_de, & ! The derivatives of dS_kb and dS_Lay with E
1636  ! [R H-1 ~> kg m-4 or m-1].
1637  derror_de, & ! The derivative of err with E [H ~> m or kg m-2].
1638  err, & ! The "error" whose zero is being sought [H2 ~> m2 or kg2 m-4].
1639  e_min, e_max, & ! The minimum and maximum values of E [H ~> m or kg m-2].
1640  error_mine, error_maxe ! err when E = E_min or E = E_max [H2 ~> m2 or kg2 m-4].
1641  real :: err_est ! An estimate of what err will be [H2 ~> m2 or kg2 m-4].
1642  real :: eL ! 1 or 0, depending on whether increases in E lead
1643  ! to decreases in the entrainment from below by the
1644  ! deepest buffer layer [nondim].
1645  real :: fa ! Temporary variable used to calculate err [nondim].
1646  real :: fk ! Temporary variable used to calculate err [H2 ~> m2 or kg2 m-4].
1647  real :: fm, fr ! Temporary variables used to calculate err [H ~> m or kg m-2].
1648  real :: tolerance ! The tolerance within which E must be converged [H ~> m or kg m-2].
1649  real :: E_prev ! The previous value of E [H ~> m or kg m-2].
1650  logical, dimension(SZI_(G)) :: false_position ! If true, the false position
1651  ! method might be used for the next iteration.
1652  logical, dimension(SZI_(G)) :: redo_i ! If true, more work is needed on this column.
1653  logical :: do_any
1654  real :: large_err ! A large error measure [H2 ~> m2 or kg2 m-4].
1655  integer :: i, it
1656  integer, parameter :: MAXIT = 30
1657 
1658  if (.not.cs%bulkmixedlayer) then
1659  call mom_error(fatal, "determine_Ea_kb should not be called "//&
1660  "unless BULKMIXEDLAYER is defined.")
1661  endif
1662  tolerance = cs%Tolerance_Ent
1663  large_err = gv%m_to_H**2 * 1.0e30
1664 
1665  do i=is,ie ; redo_i(i) = do_i(i) ; enddo
1666 
1667  do i=is,ie ; if (do_i(i)) then
1668  ! The first guess of Ent was the value from the previous iteration.
1669 
1670  ! These were previously calculated and provide good limits and estimates
1671  ! of the errors there. By construction the errors increase with R*ea_kbp1.
1672  e_min(i) = min_eakb(i) ; e_max(i) = max_eakb(i)
1673  error_mine(i) = -large_err ; error_maxe(i) = large_err
1674  false_position(i) = .true. ! Used to alternate between false_position and
1675  ! bisection when Newton's method isn't working.
1676  if (present(err_min_eakb0)) error_mine(i) = err_min_eakb0(i) - e_min(i) * ea_kbp1(i)
1677  if (present(err_max_eakb0)) error_maxe(i) = err_max_eakb0(i) - e_max(i) * ea_kbp1(i)
1678 
1679  if ((error_maxe(i) <= 0.0) .or. (error_mine(i) >= 0.0)) then
1680  ! The root is not bracketed and one of the limiting values should be used.
1681  if (error_maxe(i) <= 0.0) then
1682  ! The errors decrease with E*ea_kbp1, so E_max is the best solution.
1683  ent(i) = e_max(i) ; err(i) = error_maxe(i)
1684  else ! error_minE >= 0 is equivalent to ea_kbp1 = 0.0.
1685  ent(i) = e_min(i) ; err(i) = error_mine(i)
1686  endif
1687  derror_de(i) = 0.0
1688  redo_i(i) = .false.
1689  endif
1690  endif ; enddo ! End of i-loop
1691 
1692  do it = 1,maxit
1693  do_any = .false. ; do i=is,ie ; if (redo_i(i)) do_any = .true. ; enddo
1694  if (.not.do_any) exit
1695  call determine_dskb(h_bl, sref, ent_bl, ent, is, ie, kmb, g, gv, .true., ds_kb, &
1696  ddskb_de, ds_lay, ddslay_de, do_i_in = redo_i)
1697  do i=is,ie ; if (redo_i(i)) then
1698  ! The correct root is bracketed between E_min and E_max.
1699  ! Note the following limits: Ent >= 0 ; fa > 1 ; fk > 0
1700  el = 0.0 ; if (2.0*ent_bl(i,kmb+1) >= ent(i)) el = 1.0
1701  fa = (1.0 + el) + ds_kb(i)*i_dskbp1(i)
1702  fk = dtkd_kb(i) * (ds_lay(i)/ds_kb(i))
1703  fm = (ea_kbp1(i) - h_bl(i,kmb+1)) + el*2.0*ent_bl(i,kmb+1)
1704  if (fm > -gv%Angstrom_H) fm = fm + gv%Angstrom_H ! This could be smooth if need be.
1705  err(i) = (fa * ent(i)**2 - fm * ent(i)) - fk
1706  derror_de(i) = ((2.0*fa + (ddskb_de(i)*i_dskbp1(i))*ent(i))*ent(i) - fm) - &
1707  dtkd_kb(i) * (ddslay_de(i)*ds_kb(i) - ddskb_de(i)*ds_lay(i))/(ds_kb(i)**2)
1708 
1709  if (err(i) == 0.0) then
1710  redo_i(i) = .false. ; cycle
1711  elseif (err(i) > 0.0) then
1712  e_max(i) = ent(i) ; error_maxe(i) = err(i)
1713  else
1714  e_min(i) = ent(i) ; error_mine(i) = err(i)
1715  endif
1716 
1717  e_prev = ent(i)
1718  if ((it == 1) .or. (derror_de(i) <= 0.0)) then
1719  ! Assuming that the coefficients of the quadratic equation are correct
1720  ! will usually give a very good first guess. Also, if derror_dE < 0.0,
1721  ! R is on the wrong side of the approximate parabola. In either case,
1722  ! try assuming that the error is approximately a parabola and solve.
1723  fr = sqrt(fm**2 + 4.0*fa*fk)
1724  if (fm >= 0.0) then
1725  ent(i) = (fm + fr) / (2.0 * fa)
1726  else
1727  ent(i) = (2.0 * fk) / (fr - fm)
1728  endif
1729  ! But make sure that the root stays bracketed, bisecting if needed.
1730  if ((ent(i) > e_max(i)) .or. (ent(i) < e_min(i))) &
1731  ent(i) = 0.5*(e_max(i) + e_min(i))
1732  elseif (((e_max(i)-ent(i))*derror_de(i) > -err(i)) .and. &
1733  ((ent(i)-e_min(i))*derror_de(i) > err(i)) ) then
1734  ! Use Newton's method for the next estimate, provided it will
1735  ! remain bracketed between Rmin and Rmax.
1736  ent(i) = ent(i) - err(i) / derror_de(i)
1737  elseif (false_position(i) .and. &
1738  (error_maxe(i) - error_mine(i) < 0.9*large_err)) then
1739  ! Use the false postion method if there are decent error estimates.
1740  ent(i) = e_min(i) + (e_max(i)-e_min(i)) * &
1741  (-error_mine(i)/(error_maxe(i) - error_mine(i)))
1742  false_position(i) = .false.
1743  else ! Bisect as a last resort or if the false position method was used last.
1744  ent(i) = 0.5*(e_max(i) + e_min(i))
1745  false_position(i) = .true.
1746  endif
1747 
1748  if (abs(e_prev - ent(i)) < tolerance) then
1749  err_est = err(i) + (ent(i) - e_prev) * derror_de(i)
1750  if ((it > 1) .or. (err_est*err(i) <= 0.0) .or. &
1751  (abs(err_est) < abs(tolerance*derror_de(i)))) redo_i(i) = .false.
1752  endif
1753 
1754  endif ; enddo ! End of i-loop
1755  enddo ! End of iterations to determine Ent(i).
1756 
1757  ! Update the value of dS_kb for consistency with Ent.
1758  if (present(f_kb) .or. present(dfdfm_kb)) &
1759  call determine_dskb(h_bl, sref, ent_bl, ent, is, ie, kmb, g, gv, .true., &
1760  ds_kb, do_i_in = do_i)
1761 
1762  if (present(f_kb)) then ; do i=is,ie ; if (do_i(i)) then
1763  f_kb(i) = ent(i) * (ds_kb(i) * i_dskbp1(i))
1764  endif ; enddo ; endif
1765  if (present(error)) then ; do i=is,ie ; if (do_i(i)) then
1766  error(i) = err(i)
1767  endif ; enddo ; endif
1768  if (present(dfdfm_kb)) then ; do i=is,ie ; if (do_i(i)) then
1769  ! derror_dE and ddSkb_dE are _not_ recalculated here, since dFdfm_kb is
1770  ! only used in Newton's method, and slightly increasing the accuracy of the
1771  ! estimate is unlikely to speed convergence.
1772  if (derror_de(i) > 0.0) then
1773  dfdfm_kb(i) = ((ds_kb(i) + ent(i) * ddskb_de(i)) * i_dskbp1(i)) * &
1774  (ent(i) / derror_de(i))
1775  else ! Use Adcroft's division by 0 convention.
1776  dfdfm_kb(i) = 0.0
1777  endif
1778  endif ; enddo ; endif
1779 

◆ entrain_diffusive_end()

subroutine, public mom_entrain_diffusive::entrain_diffusive_end ( type(entrain_diffusive_cs), pointer  CS)

This subroutine cleans up and deallocates any memory associated with the entrain_diffusive module.

Parameters
csA pointer to the control structure for this module that will be deallocated.

Definition at line 2148 of file MOM_entrain_diffusive.F90.

2148  type(entrain_diffusive_CS), pointer :: CS !< A pointer to the control structure for this
2149  !! module that will be deallocated.
2150  if (associated(cs)) deallocate(cs)
2151 

◆ entrain_diffusive_init()

subroutine, public mom_entrain_diffusive::entrain_diffusive_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(entrain_diffusive_cs), pointer  CS,
logical, intent(in), optional  just_read_params 
)

This subroutine initializes the parameters and memory associated with the entrain_diffusive module.

Parameters
[in]timeThe current model time.
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[in]param_fileA structure to parse for run-time parameters.
[in,out]diagA structure that is used to regulate diagnostic output.
csA pointer that is set to point to the control structure.
[in]just_read_paramsIf present and true, this call will only read parameters logging them or registering any diagnostics

Definition at line 2081 of file MOM_entrain_diffusive.F90.

2081  type(time_type), intent(in) :: Time !< The current model time.
2082  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
2083  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
2084  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
2085  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
2086  !! parameters.
2087  type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate diagnostic
2088  !! output.
2089  type(entrain_diffusive_CS), pointer :: CS !< A pointer that is set to point to the control
2090  !! structure.
2091  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
2092  !! only read parameters logging them or registering
2093  !! any diagnostics
2094 
2095  ! Local variables
2096  real :: dt ! The dynamics timestep, used here in the default for TOLERANCE_ENT, in MKS units [s]
2097  real :: Kd ! A diffusivity used in the default for TOLERANCE_ENT, in MKS units [m2 s-1]
2098  logical :: just_read ! If true, just read parameters but do nothing else.
2099  ! This include declares and sets the variable "version".
2100 # include "version_variable.h"
2101  character(len=40) :: mdl = "MOM_entrain_diffusive" ! This module's name.
2102 
2103  if (associated(cs)) then
2104  call mom_error(warning, "entrain_diffusive_init called with an associated "// &
2105  "control structure.")
2106  return
2107  endif
2108  allocate(cs)
2109 
2110  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
2111 
2112  cs%diag => diag
2113 
2114  cs%bulkmixedlayer = (gv%nkml > 0)
2115 
2116 ! Set default, read and log parameters
2117  if (.not.just_read) call log_version(param_file, mdl, version, "")
2118  call get_param(param_file, mdl, "MAX_ENT_IT", cs%max_ent_it, &
2119  "The maximum number of iterations that may be used to "//&
2120  "calculate the interior diapycnal entrainment.", default=5, do_not_log=just_read)
2121  ! In this module, KD is only used to set the default for TOLERANCE_ENT. [m2 s-1]
2122  call get_param(param_file, mdl, "KD", kd, default=0.0)
2123  call get_param(param_file, mdl, "DT", dt, &
2124  "The (baroclinic) dynamics time step.", units = "s", &
2125  fail_if_missing=.true., do_not_log=just_read)
2126  call get_param(param_file, mdl, "TOLERANCE_ENT", cs%Tolerance_Ent, &
2127  "The tolerance with which to solve for entrainment values.", &
2128  units="m", default=max(100.0*gv%Angstrom_m,1.0e-4*sqrt(dt*kd)), scale=gv%m_to_H, &
2129  do_not_log=just_read)
2130 
2131  cs%Rho_sig_off = 1000.0*us%kg_m3_to_R
2132 
2133  if (.not.just_read) then
2134  cs%id_Kd = register_diag_field('ocean_model', 'Kd_effective', diag%axesTL, time, &
2135  'Diapycnal diffusivity as applied', 'm2 s-1', conversion=us%Z2_T_to_m2_s)
2136  cs%id_diff_work = register_diag_field('ocean_model', 'diff_work', diag%axesTi, time, &
2137  'Work actually done by diapycnal diffusion across each interface', &
2138  'W m-2', conversion=us%RZ3_T3_to_W_m2)
2139  endif
2140 
2141  if (just_read) deallocate(cs)
2142 

◆ entrainment_diffusive()

subroutine, public mom_entrain_diffusive::entrainment_diffusive ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
type(forcing), intent(in)  fluxes,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(entrain_diffusive_cs), pointer  CS,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(out)  ea,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(out)  eb,
integer, dimension( g %isd: g %ied, g %jsd: g %jed), intent(inout), optional  kb_out,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in), optional  Kd_Lay,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke+1), intent(in), optional  Kd_int 
)

This subroutine calculates ea and eb, the rates at which a layer entrains from the layers above and below. The entrainment rates are proportional to the buoyancy flux in a layer and inversely proportional to the density differences between layers. The scheme that is used here is described in detail in Hallberg, Mon. Wea. Rev. 2000.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[in]hLayer thicknesses [H ~> m or kg m-2].
[in]tvA structure containing pointers to any available thermodynamic fields. Absent fields have NULL ptrs.
[in]fluxesA structure of surface fluxes that may be used.
[in]dtThe time increment [T ~> s].
csThe control structure returned by a previous call to entrain_diffusive_init.
[out]eaThe amount of fluid entrained from the layer
[out]ebThe amount of fluid entrained from the layer
[in,out]kb_outThe index of the lightest layer denser than
[in]kd_layThe diapycnal diffusivity of layers
[in]kd_intThe diapycnal diffusivity of interfaces

Definition at line 52 of file MOM_entrain_diffusive.F90.

52  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
53  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
54  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
55  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
56  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
57  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers to any available
58  !! thermodynamic fields. Absent fields have NULL
59  !! ptrs.
60  type(forcing), intent(in) :: fluxes !< A structure of surface fluxes that may
61  !! be used.
62  real, intent(in) :: dt !< The time increment [T ~> s].
63  type(entrain_diffusive_CS), pointer :: CS !< The control structure returned by a previous
64  !! call to entrain_diffusive_init.
65  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
66  intent(out) :: ea !< The amount of fluid entrained from the layer
67  !! above within this time step [H ~> m or kg m-2].
68  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
69  intent(out) :: eb !< The amount of fluid entrained from the layer
70  !! below within this time step [H ~> m or kg m-2].
71  integer, dimension(SZI_(G),SZJ_(G)), &
72  optional, intent(inout) :: kb_out !< The index of the lightest layer denser than
73  !! the buffer layer.
74  ! At least one of the two following arguments must be present.
75  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
76  optional, intent(in) :: Kd_Lay !< The diapycnal diffusivity of layers
77  !! [Z2 T-1 ~> m2 s-1].
78  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
79  optional, intent(in) :: Kd_int !< The diapycnal diffusivity of interfaces
80  !! [Z2 T-1 ~> m2 s-1].
81 
82 ! This subroutine calculates ea and eb, the rates at which a layer entrains
83 ! from the layers above and below. The entrainment rates are proportional to
84 ! the buoyancy flux in a layer and inversely proportional to the density
85 ! differences between layers. The scheme that is used here is described in
86 ! detail in Hallberg, Mon. Wea. Rev. 2000.
87 
88  real, dimension(SZI_(G),SZK_(G)) :: &
89  dtKd ! The layer diapycnal diffusivity times the time step [H2 ~> m2 or kg2 m-4].
90  real, dimension(SZI_(G),SZK_(G)+1) :: &
91  dtKd_int ! The diapycnal diffusivity at the interfaces times the time step [H2 ~> m2 or kg2 m-4]
92  real, dimension(SZI_(G),SZK_(G)) :: &
93  F, & ! The density flux through a layer within a time step divided by the
94  ! density difference across the interface below the layer [H ~> m or kg m-2].
95  maxf, & ! maxF is the maximum value of F that will not deplete all of the
96  ! layers above or below a layer within a timestep [H ~> m or kg m-2].
97  minf, & ! minF is the minimum flux that should be expected in the absence of
98  ! interactions between layers [H ~> m or kg m-2].
99  fprev, &! The previous estimate of F [H ~> m or kg m-2].
100  dfdfm, &! The partial derivative of F with respect to changes in F of the
101  ! neighboring layers. [nondim]
102  h_guess ! An estimate of the layer thicknesses after entrainment, but
103  ! before the entrainments are adjusted to drive the layer
104  ! densities toward their target values [H ~> m or kg m-2].
105  real, dimension(SZI_(G),SZK_(G)+1) :: &
106  Ent_bl ! The average entrainment upward and downward across
107  ! each interface around the buffer layers [H ~> m or kg m-2].
108  real, allocatable, dimension(:,:,:) :: &
109  Kd_eff, & ! The effective diffusivity that actually applies to each
110  ! layer after the effects of boundary conditions are
111  ! considered [Z2 T-1 ~> m2 s-1].
112  diff_work ! The work actually done by diffusion across each
113  ! interface [R Z3 T-3 ~> W m-2]. Sum vertically for the total work.
114 
115  real :: hm, fm, fr, fk ! Work variables with units of H, H, H, and H2.
116 
117  real :: b1(SZI_(G)) ! b1 and c1 are variables used by the
118  real :: c1(SZI_(G),SZK_(G)) ! tridiagonal solver.
119 
120  real, dimension(SZI_(G)) :: &
121  htot, & ! The total thickness above or below a layer [H ~> m or kg m-2].
122  Rcv, & ! Value of the coordinate variable (potential density)
123  ! based on the simulated T and S and P_Ref [R ~> kg m-3].
124  pres, & ! Reference pressure (P_Ref) [R L2 T-2 ~> Pa].
125  eakb, & ! The entrainment from above by the layer below the buffer
126  ! layer (i.e. layer kb) [H ~> m or kg m-2].
127  ea_kbp1, & ! The entrainment from above by layer kb+1 [H ~> m or kg m-2].
128  eb_kmb, & ! The entrainment from below by the deepest buffer layer [H ~> m or kg m-2].
129  ds_kb, & ! The reference potential density difference across the
130  ! interface between the buffer layers and layer kb [R ~> kg m-3].
131  ds_anom_lim, &! The amount by which dS_kb is reduced when limits are
132  ! applied [R ~> kg m-3].
133  i_dskbp1, & ! The inverse of the potential density difference across the
134  ! interface below layer kb [R-1 ~> m3 kg-1].
135  dtkd_kb, & ! The diapycnal diffusivity in layer kb times the time step
136  ! [H2 ~> m2 or kg2 m-4].
137  maxf_correct, & ! An amount by which to correct maxF due to excessive
138  ! surface heat loss [H ~> m or kg m-2].
139  zeros, & ! An array of all zeros. (Usually used with [H ~> m or kg m-2].)
140  max_eakb, & ! The maximum value of eakb that might be realized [H ~> m or kg m-2].
141  min_eakb, & ! The minimum value of eakb that might be realized [H ~> m or kg m-2].
142  err_max_eakb0, & ! The value of error returned by determine_Ea_kb
143  err_min_eakb0, & ! when eakb = min_eakb and max_eakb and ea_kbp1 = 0.
144  err_eakb0, & ! A value of error returned by determine_Ea_kb.
145  f_kb, & ! The value of F in layer kb, or equivalently the entrainment
146  ! from below by layer kb [H ~> m or kg m-2].
147  dfdfm_kb, & ! The partial derivative of F with fm [nondim]. See dFdfm.
148  maxf_kb, & ! The maximum value of F_kb that might be realized [H ~> m or kg m-2].
149  eakb_maxf, & ! The value of eakb that gives F_kb=maxF_kb [H ~> m or kg m-2].
150  f_kb_maxent ! The value of F_kb when eakb = max_eakb [H ~> m or kg m-2].
151  real, dimension(SZI_(G),SZK_(G)) :: &
152  Sref, & ! The reference potential density of the mixed and buffer layers,
153  ! and of the two lightest interior layers (kb and kb+1) copied
154  ! into layers kmb+1 and kmb+2 [R ~> kg m-3].
155  h_bl ! The thicknesses of the mixed and buffer layers, and of the two
156  ! lightest interior layers (kb and kb+1) copied into layers kmb+1
157  ! and kmb+2 [H ~> m or kg m-2].
158 
159  real, dimension(SZI_(G),SZK_(G)) :: &
160  ds_dsp1, & ! The coordinate variable (sigma-2) difference across an
161  ! interface divided by the difference across the interface
162  ! below it. [nondim]
163  dsp1_ds, & ! The inverse coordinate variable (sigma-2) difference
164  ! across an interface times the difference across the
165  ! interface above it. [nondim]
166  i2p2dsp1_ds, & ! 1 / (2 + 2 * ds_k+1 / ds_k). [nondim]
167  grats ! 2*(2 + ds_k+1 / ds_k + ds_k / ds_k+1) =
168  ! 4*ds_Lay*(1/ds_k + 1/ds_k+1). [nondim]
169 
170  real :: dRHo ! The change in locally referenced potential density between
171  ! the layers above and below an interface [R ~> kg m-3].
172  real :: g_2dt ! 0.5 * G_Earth / dt, times unit conversion factors
173  ! [m3 H-2 s-2 T-1 ~> m s-3 or m7 kg-2 s-3].
174  real, dimension(SZI_(G)) :: &
175  pressure, & ! The pressure at an interface [R L2 T-2 ~> Pa].
176  T_eos, S_eos, & ! The potential temperature and salinity at which to
177  ! evaluate dRho_dT and dRho_dS [degC] and [ppt].
178  drho_dt, drho_ds ! The partial derivatives of potential density with temperature and
179  ! salinity, [R degC-1 ~> kg m-3 degC-1] and [R ppt-1 ~> kg m-3 ppt-1].
180 
181  real :: tolerance ! The tolerance within which E must be converged [H ~> m or kg m-2].
182  real :: Angstrom ! The minimum layer thickness [H ~> m or kg m-2].
183  real :: h_neglect ! A thickness that is so small it is usually lost
184  ! in roundoff and can be neglected [H ~> m or kg m-2].
185  real :: F_cor ! A correction to the amount of F that is used to
186  ! entrain from the layer above [H ~> m or kg m-2].
187  real :: Kd_here ! The effective diapycnal diffusivity [H2 s-1 ~> m2 s-1 or kg2 m-4 s-1].
188  real :: h_avail ! The thickness that is available for entrainment [H ~> m or kg m-2].
189  real :: dS_kb_eff ! The value of dS_kb after limiting is taken into account.
190  real :: Rho_cor ! The depth-integrated potential density anomaly that
191  ! needs to be corrected for [H R ~> kg m-2 or kg2 m-5].
192  real :: ea_cor ! The corrective adjustment to eakb [H ~> m or kg m-2].
193  real :: h1 ! The layer thickness after entrainment through the
194  ! interface below is taken into account [H ~> m or kg m-2].
195  real :: Idt ! The inverse of the time step [T-1 ~> s-1].
196 
197  logical :: do_any
198  logical :: do_entrain_eakb ! True if buffer layer is entrained
199  logical :: do_i(SZI_(G)), did_i(SZI_(G)), reiterate
200  integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
201  integer :: it, i, j, k, is, ie, js, je, nz, K2, kmb
202  integer :: kb(SZI_(G)) ! The value of kb in row j.
203  integer :: kb_min ! The minimum value of kb in the current j-row.
204  integer :: kb_min_act ! The minimum active value of kb in the current j-row.
205  integer :: is1, ie1 ! The minimum and maximum active values of i in the current j-row.
206  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
207  angstrom = gv%Angstrom_H
208  h_neglect = gv%H_subroundoff
209 
210  if (.not. associated(cs)) call mom_error(fatal, &
211  "MOM_entrain_diffusive: Module must be initialized before it is used.")
212 
213  if (.not.(present(kd_lay) .or. present(kd_int))) call mom_error(fatal, &
214  "MOM_entrain_diffusive: Either Kd_Lay or Kd_int must be present in call.")
215 
216  if ((.not.cs%bulkmixedlayer .and. .not.associated(fluxes%buoy)) .and. &
217  (associated(fluxes%lprec) .or. associated(fluxes%evap) .or. &
218  associated(fluxes%sens) .or. associated(fluxes%sw))) then
219  if (is_root_pe()) call mom_error(note, "Calculate_Entrainment: &
220  &The code to handle evaporation and precipitation without &
221  &a bulk mixed layer has not been implemented.")
222  if (is_root_pe()) call mom_error(fatal, &
223  "Either define BULKMIXEDLAYER in MOM_input or use fluxes%buoy &
224  &and a linear equation of state to drive the model.")
225  endif
226 
227  tolerance = cs%Tolerance_Ent
228  kmb = gv%nk_rho_varies
229  k2 = max(kmb+1,2) ; kb_min = k2
230  if (.not. cs%bulkmixedlayer) then
231  kb(:) = 1
232  ! These lines fill in values that are arbitrary, but needed because
233  ! they are used to normalize the buoyancy flux in layer nz.
234  do i=is,ie ; ds_dsp1(i,nz) = 2.0 ; dsp1_ds(i,nz) = 0.5 ; enddo
235  else
236  kb(:) = 0
237  do i=is,ie ; ds_dsp1(i,nz) = 0.0 ; dsp1_ds(i,nz) = 0.0 ; enddo
238  endif
239 
240  if (cs%id_diff_work > 0) allocate(diff_work(g%isd:g%ied,g%jsd:g%jed,nz+1))
241  if (cs%id_Kd > 0) allocate(kd_eff(g%isd:g%ied,g%jsd:g%jed,nz))
242 
243  if (associated(tv%eqn_of_state)) then
244  pres(:) = tv%P_Ref
245  else
246  pres(:) = 0.0
247  endif
248  eosdom(:) = eos_domain(g%HI)
249 
250  !$OMP parallel do default(none) shared(is,ie,js,je,nz,Kd_Lay,G,GV,US,dt,CS,h,tv, &
251  !$OMP kmb,Angstrom,fluxes,K2,h_neglect,tolerance, &
252  !$OMP ea,eb,Kd_int,Kd_eff,EOSdom,diff_work,g_2dt, kb_out) &
253  !$OMP firstprivate(kb,ds_dsp1,dsp1_ds,pres,kb_min) &
254  !$OMP private(dtKd,dtKd_int,do_i,Ent_bl,dtKd_kb,h_bl, &
255  !$OMP I2p2dsp1_ds,grats,htot,max_eakb,I_dSkbp1, &
256  !$OMP zeros,maxF_kb,maxF,ea_kbp1,eakb,Sref, &
257  !$OMP maxF_correct,do_any,do_entrain_eakb, &
258  !$OMP err_min_eakb0,err_max_eakb0,eakb_maxF, &
259  !$OMP min_eakb,err_eakb0,F,minF,hm,fk,F_kb_maxent,&
260  !$OMP F_kb,is1,ie1,kb_min_act,dFdfm_kb,b1,dFdfm, &
261  !$OMP Fprev,fm,fr,c1,reiterate,eb_kmb,did_i, &
262  !$OMP h_avail,h_guess,dS_kb,Rcv,F_cor,dS_kb_eff, &
263  !$OMP Rho_cor,ea_cor,h1,Idt,Kd_here,pressure, &
264  !$OMP T_eos,S_eos,dRho_dT,dRho_dS,dRho,dS_anom_lim)
265  do j=js,je
266  do i=is,ie ; kb(i) = 1 ; enddo
267 
268  if (present(kd_lay)) then
269  do k=1,nz ; do i=is,ie
270  dtkd(i,k) = gv%Z_to_H**2 * (dt * kd_lay(i,j,k))
271  enddo ; enddo
272  if (present(kd_int)) then
273  do k=1,nz+1 ; do i=is,ie
274  dtkd_int(i,k) = gv%Z_to_H**2 * (dt * kd_int(i,j,k))
275  enddo ; enddo
276  else
277  do k=2,nz ; do i=is,ie
278  dtkd_int(i,k) = gv%Z_to_H**2 * (0.5 * dt * (kd_lay(i,j,k-1) + kd_lay(i,j,k)))
279  enddo ; enddo
280  endif
281  else ! Kd_int must be present, or there already would have been an error.
282  do k=1,nz ; do i=is,ie
283  dtkd(i,k) = gv%Z_to_H**2 * (0.5 * dt * (kd_int(i,j,k)+kd_int(i,j,k+1)))
284  enddo ; enddo
285  dO k=1,nz+1 ; do i=is,ie
286  dtkd_int(i,k) = gv%Z_to_H**2 * (dt * kd_int(i,j,k))
287  enddo ; enddo
288  endif
289 
290  do i=is,ie ; do_i(i) = (g%mask2dT(i,j) > 0.5) ; enddo
291  do i=is,ie ; ds_dsp1(i,nz) = 0.0 ; enddo
292  do i=is,ie ; dsp1_ds(i,nz) = 0.0 ; enddo
293 
294  do k=2,nz-1 ; do i=is,ie
295  ds_dsp1(i,k) = gv%g_prime(k) / gv%g_prime(k+1)
296  enddo ; enddo
297 
298  if (cs%bulkmixedlayer) then
299  ! This subroutine determines the averaged entrainment across each
300  ! interface and causes thin and relatively light interior layers to be
301  ! entrained by the deepest buffer layer. This also determines kb.
302  call set_ent_bl(h, dtkd_int, tv, kb, kmb, do_i, g, gv, us, cs, j, ent_bl, sref, h_bl)
303 
304  do i=is,ie
305  dtkd_kb(i) = 0.0 ; if (kb(i) < nz) dtkd_kb(i) = dtkd(i,kb(i))
306  enddo
307  else
308  do i=is,ie ; ent_bl(i,kmb+1) = 0.0 ; enddo
309  endif
310 
311  do k=2,nz-1 ; do i=is,ie
312  dsp1_ds(i,k) = 1.0 / ds_dsp1(i,k)
313  i2p2dsp1_ds(i,k) = 0.5/(1.0+dsp1_ds(i,k))
314  grats(i,k) = 2.0*(2.0+(dsp1_ds(i,k)+ds_dsp1(i,k)))
315  enddo ; enddo
316 
317 ! Determine the maximum flux, maxF, for each of the isopycnal layers.
318 ! Also determine when the fluxes start entraining
319 ! from various buffer or mixed layers, where appropriate.
320  if (cs%bulkmixedlayer) then
321  kb_min = nz
322  do i=is,ie
323  htot(i) = h(i,j,1) - angstrom
324  enddo
325  do k=2,kmb ; do i=is,ie
326  htot(i) = htot(i) + (h(i,j,k) - angstrom)
327  enddo ; enddo
328  do i=is,ie
329  max_eakb(i) = max(ent_bl(i,kmb+1) + 0.5*htot(i), htot(i))
330  i_dskbp1(i) = 1.0 / (sref(i,kmb+2) - sref(i,kmb+1))
331  zeros(i) = 0.0
332  enddo
333 
334  ! Find the maximum amount of entrainment from below that the top
335  ! interior layer could exhibit (maxF(i,kb)), approximating that
336  ! entrainment by (eakb*max(dS_kb/dSkbp1,0)). eakb is in the range
337  ! from 0 to max_eakb.
338  call find_maxf_kb(h_bl, sref, ent_bl, i_dskbp1, zeros, max_eakb, kmb, &
339  is, ie, g, gv, cs, maxf_kb, eakb_maxf, do_i, f_kb_maxent)
340  do i=is,ie ; if (kb(i) <= nz) then
341  maxf(i,kb(i)) = max(0.0, maxf_kb(i), f_kb_maxent(i))
342  if ((maxf_kb(i) > f_kb_maxent(i)) .and. (eakb_maxf(i) >= htot(i))) &
343  max_eakb(i) = eakb_maxf(i)
344  endif ; enddo
345 
346  do i=is,ie ; ea_kbp1(i) = 0.0 ; eakb(i) = max_eakb(i) ; enddo
347  call determine_ea_kb(h_bl, dtkd_kb, sref, i_dskbp1, ent_bl, ea_kbp1, &
348  max_eakb, max_eakb, kmb, is, ie, do_i, g, gv, cs, eakb, &
349  error=err_max_eakb0, f_kb=f_kb)
350 
351  ! The maximum value of F(kb) between htot and max_eakb determines
352  ! what maxF(kb+1) should be.
353  do i=is,ie ; min_eakb(i) = min(htot(i), max_eakb(i)) ; enddo
354  call find_maxf_kb(h_bl, sref, ent_bl, i_dskbp1, min_eakb, max_eakb, &
355  kmb, is, ie, g, gv, cs, f_kb_maxent, do_i_in = do_i)
356 
357  do i=is,ie
358  do_entrain_eakb = .false.
359  ! If error_max_eakb0 < 0, then buffer layers are always all entrained
360  if (do_i(i)) then ; if (err_max_eakb0(i) < 0.0) then
361  do_entrain_eakb = .true.
362  endif ; endif
363 
364  if (do_entrain_eakb) then
365  eakb(i) = max_eakb(i) ; min_eakb(i) = max_eakb(i)
366  else
367  eakb(i) = 0.0 ; min_eakb(i) = 0.0
368  endif
369  enddo
370 
371  ! Find the amount of entrainment of the buffer layers that would occur
372  ! if there were no entrainment by the deeper interior layers. Also find
373  ! how much entrainment of the deeper layers would occur.
374  call determine_ea_kb(h_bl, dtkd_kb, sref, i_dskbp1, ent_bl, ea_kbp1, &
375  zeros, max_eakb, kmb, is, ie, do_i, g, gv, cs, min_eakb, &
376  error=err_min_eakb0, f_kb=f_kb, err_max_eakb0=err_max_eakb0)
377  ! Error_min_eakb0 should be ~0 unless error_max_eakb0 < 0.
378  do i=is,ie ; if ((kb(i)<nz) .and. (kb_min>kb(i))) kb_min = kb(i) ; enddo
379  else
380  ! Without a bulk mixed layer, surface fluxes are applied in this
381  ! subroutine. (Otherwise, they are handled in mixedlayer.)
382  ! Initially the maximum flux in layer zero is given by the surface
383  ! buoyancy flux. It will later be limited if the surface flux is
384  ! too large. Here buoy is the surface buoyancy flux.
385  do i=is,ie
386  maxf(i,1) = 0.0
387  htot(i) = h(i,j,1) - angstrom
388  enddo
389  if (associated(fluxes%buoy)) then ; do i=is,ie
390  maxf(i,1) = gv%Z_to_H * (dt*fluxes%buoy(i,j)) / gv%g_prime(2)
391  enddo ; endif
392  endif
393 
394 ! The following code calculates the maximum flux, maxF, for the interior
395 ! layers.
396  do k=kb_min,nz-1 ; do i=is,ie
397  if ((k == kb(i)+1) .and. cs%bulkmixedlayer) then
398  maxf(i,k) = ds_dsp1(i,k)*(f_kb_maxent(i) + htot(i))
399  htot(i) = htot(i) + (h(i,j,k) - angstrom)
400  elseif (k > kb(i)) then
401  maxf(i,k) = ds_dsp1(i,k)*(maxf(i,k-1) + htot(i))
402  htot(i) = htot(i) + (h(i,j,k) - angstrom)
403  endif
404  enddo ; enddo
405  do i=is,ie
406  maxf(i,nz) = 0.0
407  if (.not.cs%bulkmixedlayer) then
408  maxf_correct(i) = max(0.0, -(maxf(i,nz-1) + htot(i)))
409  endif
410  htot(i) = h(i,j,nz) - angstrom
411  enddo
412  if (.not.cs%bulkmixedlayer) then
413  do_any = .false. ; do i=is,ie ; if (maxf_correct(i) > 0.0) do_any = .true. ; enddo
414  if (do_any) then
415  do k=nz-1,1,-1 ; do i=is,ie
416  maxf(i,k) = maxf(i,k) + maxf_correct(i)
417  maxf_correct(i) = maxf_correct(i) * dsp1_ds(i,k)
418  enddo ; enddo
419  endif
420  endif
421  do k=nz-1,kb_min,-1 ; do i=is,ie ; if (do_i(i)) then
422  if (k >= kb(i)) then
423  maxf(i,k) = min(maxf(i,k),dsp1_ds(i,k+1)*maxf(i,k+1) + htot(i))
424  htot(i) = htot(i) + (h(i,j,k) - angstrom)
425  endif
426  if (k == kb(i)) then
427  if ((maxf(i,k) < f_kb(i)) .or. (maxf(i,k) < maxf_kb(i)) &
428  .and. (eakb_maxf(i) <= max_eakb(i))) then
429  ! In this case, too much was being entrained by the topmost interior
430  ! layer, even with the minimum initial estimate. The buffer layer
431  ! will always entrain the maximum amount.
432  f_kb(i) = maxf(i,k)
433  if ((f_kb(i) <= maxf_kb(i)) .and. (eakb_maxf(i) <= max_eakb(i))) then
434  eakb(i) = eakb_maxf(i)
435  else
436  eakb(i) = max_eakb(i)
437  endif
438  call f_kb_to_ea_kb(h_bl, sref, ent_bl, i_dskbp1, f_kb, kmb, i, &
439  g, gv, cs, eakb, angstrom)
440  if ((eakb(i) < max_eakb(i)) .or. (eakb(i) < min_eakb(i))) then
441  call determine_ea_kb(h_bl, dtkd_kb, sref, i_dskbp1, ent_bl, zeros, &
442  eakb, eakb, kmb, i, i, do_i, g, gv, cs, eakb, &
443  error=err_eakb0)
444  if (eakb(i) < max_eakb(i)) then
445  max_eakb(i) = eakb(i) ; err_max_eakb0(i) = err_eakb0(i)
446  endif
447  if (eakb(i) < min_eakb(i)) then
448  min_eakb(i) = eakb(i) ; err_min_eakb0(i) = err_eakb0(i)
449  endif
450  endif
451  endif
452  endif
453  endif ; enddo ; enddo
454  if (.not.cs%bulkmixedlayer) then
455  do i=is,ie
456  maxf(i,1) = min(maxf(i,1),dsp1_ds(i,2)*maxf(i,2) + htot(i))
457  enddo
458  endif
459 
460 ! The following code provides an initial estimate of the flux in
461 ! each layer, F. The initial guess for the layer diffusive flux is
462 ! the smaller of a forward discretization or the maximum diffusive
463 ! flux starting from zero thickness in one time step without
464 ! considering adjacent layers.
465  do i=is,ie
466  f(i,1) = maxf(i,1)
467  f(i,nz) = maxf(i,nz) ; minf(i,nz) = 0.0
468  enddo
469  do k=nz-1,k2,-1
470  do i=is,ie
471  if ((k==kb(i)) .and. (do_i(i))) then
472  eakb(i) = min_eakb(i)
473  minf(i,k) = 0.0
474  elseif ((k>kb(i)) .and. (do_i(i))) then
475 ! Here the layer flux is estimated, assuming no entrainment from
476 ! the surrounding layers. The estimate is a forward (steady) flux,
477 ! limited by the maximum flux for a layer starting with zero
478 ! thickness. This is often a good guess and leads to few iterations.
479  hm = h(i,j,k) + h_neglect
480  ! Note: Tried sqrt((0.5*ds_dsp1(i,k))*dtKd(i,k)) for the second limit,
481  ! but it usually doesn't work as well.
482  f(i,k) = min(maxf(i,k), sqrt(ds_dsp1(i,k)*dtkd(i,k)), &
483  0.5*(ds_dsp1(i,k)+1.0) * (dtkd(i,k) / hm))
484 
485 ! Calculate the minimum flux that can be expected if there is no entrainment
486 ! from the neighboring layers. The 0.9 is used to give used to give a 10%
487 ! known error tolerance.
488  fk = dtkd(i,k) * grats(i,k)
489  minf(i,k) = min(maxf(i,k), &
490  0.9*(i2p2dsp1_ds(i,k) * fk / (hm + sqrt(hm*hm + fk))))
491  if (k==kb(i)) minf(i,k) = 0.0 ! BACKWARD COMPATIBILITY - DELETE LATER?
492  else
493  f(i,k) = 0.0
494  minf(i,k) = 0.0
495  endif
496  enddo ! end of i loop
497  enddo ! end of k loop
498 
499  ! This is where the fluxes are actually calculated.
500 
501  is1 = ie+1 ; ie1 = is-1
502  do i=is,ie ; if (do_i(i)) then ; is1 = i ; exit ; endif ; enddo
503  do i=ie,is,-1 ; if (do_i(i)) then ; ie1 = i ; exit ; endif ; enddo
504 
505  if (cs%bulkmixedlayer) then
506  kb_min_act = nz
507  do i=is,ie
508  if (do_i(i) .and. (kb(i) < kb_min_act)) kb_min_act = kb(i)
509  enddo
510  ! Solve for the entrainment rate from above in the topmost interior
511  ! layer, eakb, such that
512  ! eakb*dS_implicit = dt*Kd*dS_layer_implicit / h_implicit.
513  do i=is1,ie1
514  ea_kbp1(i) = 0.0
515  if (do_i(i) .and. (kb(i) < nz)) &
516  ea_kbp1(i) = dsp1_ds(i,kb(i)+1)*f(i,kb(i)+1)
517  enddo
518  call determine_ea_kb(h_bl, dtkd_kb, sref, i_dskbp1, ent_bl, ea_kbp1, min_eakb, &
519  max_eakb, kmb, is1, ie1, do_i, g, gv, cs, eakb, f_kb=f_kb, &
520  err_max_eakb0=err_max_eakb0, err_min_eakb0=err_min_eakb0, &
521  dfdfm_kb=dfdfm_kb)
522  else
523  kb_min_act = kb_min
524  endif
525 
526  do it=0,cs%max_ent_it-1
527  do i=is1,ie1 ; if (do_i(i)) then
528  if (.not.cs%bulkmixedlayer) f(i,1) = min(f(i,1),maxf(i,1))
529  b1(i) = 1.0
530  endif ; enddo ! end of i loop
531 
532  ! F_kb has already been found for this iteration, either above or at
533  ! the end of the code for the previous iteration.
534  do k=kb_min_act,nz-1 ; do i=is1,ie1 ; if (do_i(i) .and. (k>=kb(i))) then
535  ! Calculate the flux in layer k.
536  if (cs%bulkmixedlayer .and. (k==kb(i))) then
537  f(i,k) = f_kb(i)
538  dfdfm(i,k) = dfdfm_kb(i)
539  else ! k > kb(i)
540  fprev(i,k) = f(i,k)
541  fm = (f(i,k-1) - h(i,j,k)) + dsp1_ds(i,k+1)*f(i,k+1)
542  fk = grats(i,k)*dtkd(i,k)
543  fr = sqrt(fm*fm + fk)
544 
545  if (fm>=0) then
546  f(i,k) = min(maxf(i,k), i2p2dsp1_ds(i,k) * (fm+fr))
547  else
548  f(i,k) = min(maxf(i,k), i2p2dsp1_ds(i,k) * (fk / (-1.0*fm+fr)))
549  endif
550 
551  if ((f(i,k) >= maxf(i,k)) .or. (fr == 0.0)) then
552  dfdfm(i,k) = 0.0
553  else
554  dfdfm(i,k) = i2p2dsp1_ds(i,k) * ((fr + fm) / fr)
555  endif
556 
557  if (k > k2) then
558  ! This is part of a tridiagonal solver for the actual flux.
559  c1(i,k) = dfdfm(i,k-1)*(dsp1_ds(i,k)*b1(i))
560  b1(i) = 1.0 / (1.0 - c1(i,k)*dfdfm(i,k))
561  f(i,k) = min(b1(i)*(f(i,k)-fprev(i,k)) + fprev(i,k), maxf(i,k))
562  if (f(i,k) >= maxf(i,k)) dfdfm(i,k) = 0.0
563  endif
564  endif
565  endif ; enddo ; enddo
566 
567  do k=nz-2,kb_min_act,-1 ; do i=is1,ie1
568  if (do_i(i) .and. (k > kb(i))) &
569  f(i,k) = min((f(i,k)+c1(i,k+1)*(f(i,k+1)-fprev(i,k+1))),maxf(i,k))
570  enddo ; enddo
571 
572  if (cs%bulkmixedlayer) then
573  do i=is1,ie1
574  if (do_i(i) .and. (kb(i) < nz)) then
575  ! F will be increased to minF later.
576  ea_kbp1(i) = dsp1_ds(i,kb(i)+1)*max(f(i,kb(i)+1), minf(i,kb(i)+1))
577  else
578  ea_kbp1(i) = 0.0
579  endif
580  enddo
581  call determine_ea_kb(h_bl, dtkd_kb, sref, i_dskbp1, ent_bl, ea_kbp1, min_eakb, &
582  max_eakb, kmb, is1, ie1, do_i, g, gv, cs, eakb, f_kb=f_kb, &
583  err_max_eakb0=err_max_eakb0, err_min_eakb0=err_min_eakb0, &
584  dfdfm_kb=dfdfm_kb)
585  do i=is1,ie1
586  if (do_i(i) .and. (kb(i) < nz)) f(i,kb(i)) = f_kb(i)
587  enddo
588  endif
589 
590 ! Determine whether to do another iteration.
591  if (it < cs%max_ent_it-1) then
592 
593  reiterate = .false.
594  if (cs%bulkmixedlayer) then ; do i=is1,ie1 ; if (do_i(i)) then
595  eb_kmb(i) = max(2.0*ent_bl(i,kmb+1) - eakb(i), 0.0)
596  endif ; enddo ; endif
597  do i=is1,ie1
598  did_i(i) = do_i(i) ; do_i(i) = .false.
599  enddo
600  do k=kb_min_act,nz-1 ; do i=is1,ie1
601  if (did_i(i) .and. (k >= kb(i))) then
602  if (f(i,k) < minf(i,k)) then
603  f(i,k) = minf(i,k)
604  do_i(i) = .true. ; reiterate = .true.
605  elseif (k > kb(i)) then
606  if ((abs(f(i,k) - fprev(i,k)) > tolerance) .or. &
607  ((h(i,j,k) + ((1.0+dsp1_ds(i,k))*f(i,k) - &
608  (f(i,k-1) + dsp1_ds(i,k+1)*f(i,k+1)))) < 0.9*angstrom)) then
609  do_i(i) = .true. ; reiterate = .true.
610  endif
611  else ! (k == kb(i))
612 ! A more complicated test is required for the layer beneath the buffer layer,
613 ! since its flux may be partially used to entrain directly from the mixed layer.
614 ! Negative fluxes should not occur with the bulk mixed layer.
615  if (h(i,j,k) + ((f(i,k) + eakb(i)) - &
616  (eb_kmb(i) + dsp1_ds(i,k+1)*f(i,k+1))) < 0.9*angstrom) then
617  do_i(i) = .true. ; reiterate = .true.
618  endif
619  endif
620  endif
621  enddo ; enddo
622  if (.not.reiterate) exit
623  endif ! end of if (it < CS%max_ent_it-1)
624  enddo ! end of it loop
625 ! This is the end of the section that might be iterated.
626 
627 
628  if (it == (cs%max_ent_it)) then
629  ! Limit the flux so that the layer below is not depleted.
630  ! This should only be applied to the last iteration.
631  do i=is1,ie1 ; if (do_i(i)) then
632  f(i,nz-1) = max(f(i,nz-1), min(minf(i,nz-1), 0.0))
633  if (kb(i) >= nz-1) then ; ea_kbp1(i) = 0.0 ; endif
634  endif ; enddo
635  do k=nz-2,kb_min_act,-1 ; do i=is1,ie1 ; if (do_i(i)) then
636  if (k>kb(i)) then
637  f(i,k) = min(max(minf(i,k),f(i,k)), (dsp1_ds(i,k+1)*f(i,k+1) + &
638  max(((f(i,k+1)-dsp1_ds(i,k+2)*f(i,k+2)) + &
639  (h(i,j,k+1) - angstrom)), 0.5*(h(i,j,k+1)-angstrom))))
640  elseif (k==kb(i)) then
641  ea_kbp1(i) = dsp1_ds(i,k+1)*f(i,k+1)
642  h_avail = dsp1_ds(i,k+1)*f(i,k+1) + max(0.5*(h(i,j,k+1)-angstrom), &
643  ((f(i,k+1)-dsp1_ds(i,k+2)*f(i,k+2)) + (h(i,j,k+1) - angstrom)))
644  if ((f(i,k) > 0.0) .and. (f(i,k) > h_avail)) then
645  f_kb(i) = max(0.0, h_avail)
646  f(i,k) = f_kb(i)
647  if ((f_kb(i) < maxf_kb(i)) .and. (eakb_maxf(i) <= eakb(i))) &
648  eakb(i) = eakb_maxf(i)
649  call f_kb_to_ea_kb(h_bl, sref, ent_bl, i_dskbp1, f_kb, kmb, i, &
650  g, gv, cs, eakb)
651  endif
652  endif
653  endif ; enddo ; enddo
654 
655 
656  if (cs%bulkmixedlayer) then ; do i=is1,ie1
657  if (do_i(i) .and. (kb(i) < nz)) then
658  h_avail = eakb(i) + max(0.5*(h_bl(i,kmb+1)-angstrom), &
659  (f_kb(i)-ea_kbp1(i)) + (h_bl(i,kmb+1)-angstrom))
660  ! Ensure that 0 < eb_kmb < h_avail.
661  ent_bl(i,kmb+1) = min(ent_bl(i,kmb+1),0.5*(eakb(i) + h_avail))
662 
663  eb_kmb(i) = max(2.0*ent_bl(i,kmb+1) - eakb(i), 0.0)
664  endif
665  enddo ; endif
666 
667  ! Limit the flux so that the layer above is not depleted.
668  do k=kb_min_act+1,nz-1 ; do i=is1,ie1 ; if (do_i(i)) then
669  if ((.not.cs%bulkmixedlayer) .or. (k > kb(i)+1)) then
670  f(i,k) = min(f(i,k), ds_dsp1(i,k)*( ((f(i,k-1) + &
671  dsp1_ds(i,k-1)*f(i,k-1)) - f(i,k-2)) + (h(i,j,k-1) - angstrom)))
672  f(i,k) = max(f(i,k),min(minf(i,k),0.0))
673  elseif (k == kb(i)+1) then
674  f(i,k) = min(f(i,k), ds_dsp1(i,k)*( ((f(i,k-1) + eakb(i)) - &
675  eb_kmb(i)) + (h(i,j,k-1) - angstrom)))
676  f(i,k) = max(f(i,k),min(minf(i,k),0.0))
677  endif
678  endif ; enddo ; enddo
679  endif ! (it == (CS%max_ent_it))
680 
681  call f_to_ent(f, h, kb, kmb, j, g, gv, cs, dsp1_ds, eakb, ent_bl, ea, eb)
682 
683  ! Calculate the layer thicknesses after the entrainment to constrain the
684  ! corrective fluxes.
685  if (associated(tv%eqn_of_state)) then
686  do i=is,ie
687  h_guess(i,1) = (h(i,j,1) - angstrom) + (eb(i,j,1) - ea(i,j,2))
688  h_guess(i,nz) = (h(i,j,nz) - angstrom) + (ea(i,j,nz) - eb(i,j,nz-1))
689  if (h_guess(i,1) < 0.0) h_guess(i,1) = 0.0
690  if (h_guess(i,nz) < 0.0) h_guess(i,nz) = 0.0
691  enddo
692  do k=2,nz-1 ; do i=is,ie
693  h_guess(i,k) = (h(i,j,k) - angstrom) + ((ea(i,j,k) - eb(i,j,k-1)) + &
694  (eb(i,j,k) - ea(i,j,k+1)))
695  if (h_guess(i,k) < 0.0) h_guess(i,k) = 0.0
696  enddo ; enddo
697  if (cs%bulkmixedlayer) then
698  call determine_dskb(h_bl, sref, ent_bl, eakb, is, ie, kmb, g, gv, &
699  .true., ds_kb, ds_anom_lim=ds_anom_lim)
700  do k=nz-1,kb_min,-1
701  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pres, rcv, tv%eqn_of_state, eosdom)
702  do i=is,ie
703  if ((k>kb(i)) .and. (f(i,k) > 0.0)) then
704  ! Within a time step, a layer may entrain no more than its
705  ! thickness for correction. This limitation should apply
706  ! extremely rarely, but precludes undesirable behavior.
707  ! Note: Corrected a sign/logic error & factor of 2 error, and
708  ! the layers tracked the target density better, mostly due to
709  ! the factor of 2 error.
710  f_cor = h(i,j,k) * min(1.0 , max(-ds_dsp1(i,k), &
711  (gv%Rlay(k) - rcv(i)) / (gv%Rlay(k+1)-gv%Rlay(k))) )
712 
713  ! Ensure that (1) Entrainments are positive, (2) Corrections in
714  ! a layer cannot deplete the layer itself (very generously), and
715  ! (3) a layer can take no more than a quarter the mass of its
716  ! neighbor.
717  if (f_cor > 0.0) then
718  f_cor = min(f_cor, 0.9*f(i,k), ds_dsp1(i,k)*0.5*h_guess(i,k), &
719  0.25*h_guess(i,k+1))
720  else
721  f_cor = -min(-f_cor, 0.9*f(i,k), 0.5*h_guess(i,k), &
722  0.25*ds_dsp1(i,k)*h_guess(i,k-1) )
723  endif
724 
725  ea(i,j,k) = ea(i,j,k) - dsp1_ds(i,k)*f_cor
726  eb(i,j,k) = eb(i,j,k) + f_cor
727  elseif ((k==kb(i)) .and. (f(i,k) > 0.0)) then
728  ! Rho_cor is the density anomaly that needs to be corrected,
729  ! taking into account that the true potential density of the
730  ! deepest buffer layer is not exactly what is returned as dS_kb.
731  ds_kb_eff = 2.0*ds_kb(i) - ds_anom_lim(i) ! Could be negative!!!
732  rho_cor = h(i,j,k) * (gv%Rlay(k)-rcv(i)) + eakb(i)*ds_anom_lim(i)
733 
734  ! Ensure that -.9*eakb < ea_cor < .9*eakb
735  if (abs(rho_cor) < abs(0.9*eakb(i)*ds_kb_eff)) then
736  ea_cor = -rho_cor / ds_kb_eff
737  else
738  ea_cor = sign(0.9*eakb(i),-rho_cor*ds_kb_eff)
739  endif
740 
741  if (ea_cor > 0.0) then
742  ! Ensure that -F_cor < 0.5*h_guess
743  ea_cor = min(ea_cor, 0.5*(max_eakb(i) - eakb(i)), &
744  0.5*h_guess(i,k) / (ds_kb(i) * i_dskbp1(i)))
745  else
746  ! Ensure that -ea_cor < 0.5*h_guess & F_cor < 0.25*h_guess(k+1)
747  ea_cor = -min(-ea_cor, 0.5*h_guess(i,k), &
748  0.25*h_guess(i,k+1) / (ds_kb(i) * i_dskbp1(i)))
749  endif
750 
751  ea(i,j,k) = ea(i,j,k) + ea_cor
752  eb(i,j,k) = eb(i,j,k) - (ds_kb(i) * i_dskbp1(i)) * ea_cor
753  elseif (k < kb(i)) then
754  ! Repetative, unless ea(kb) has been corrected.
755  ea(i,j,k) = ea(i,j,k+1)
756  endif
757  enddo
758  enddo
759  do k=kb_min-1,k2,-1 ; do i=is,ie
760  ea(i,j,k) = ea(i,j,k+1)
761  enddo ; enddo
762 
763  ! Repetative, unless ea(kb) has been corrected.
764  k=kmb
765  do i=is,ie
766  ! Do not adjust eb through the base of the buffer layers, but it
767  ! may be necessary to change entrainment from above.
768  h1 = (h(i,j,k) - angstrom) + (eb(i,j,k) - ea(i,j,k+1))
769  ea(i,j,k) = max(ent_bl(i,k), ent_bl(i,k)-0.5*h1, -h1)
770  enddo
771  do k=kmb-1,2,-1 ; do i=is,ie
772  ! Determine the entrainment from below for each buffer layer.
773  eb(i,j,k) = max(2.0*ent_bl(i,k+1) - ea(i,j,k+1), 0.0)
774 
775  ! Determine the entrainment from above for each buffer layer.
776  h1 = (h(i,j,k) - angstrom) + (eb(i,j,k) - ea(i,j,k+1))
777  ea(i,j,k) = max(ent_bl(i,k), ent_bl(i,k)-0.5*h1, -h1)
778  enddo ; enddo
779  do i=is,ie
780  eb(i,j,1) = max(2.0*ent_bl(i,2) - ea(i,j,2), 0.0)
781  enddo
782 
783  else ! not bulkmixedlayer
784  do k=k2,nz-1;
785  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pres, rcv, tv%eqn_of_state, eosdom)
786  do i=is,ie ; if (f(i,k) > 0.0) then
787  ! Within a time step, a layer may entrain no more than
788  ! its thickness for correction. This limitation should
789  ! apply extremely rarely, but precludes undesirable
790  ! behavior.
791  f_cor = h(i,j,k) * min(dsp1_ds(i,k) , max(-1.0, &
792  (gv%Rlay(k) - rcv(i)) / (gv%Rlay(k+1)-gv%Rlay(k))) )
793 
794  ! Ensure that (1) Entrainments are positive, (2) Corrections in
795  ! a layer cannot deplete the layer itself (very generously), and
796  ! (3) a layer can take no more than a quarter the mass of its
797  ! neighbor.
798  if (f_cor >= 0.0) then
799  f_cor = min(f_cor, 0.9*f(i,k), 0.5*dsp1_ds(i,k)*h_guess(i,k), &
800  0.25*h_guess(i,k+1))
801  else
802  f_cor = -min(-f_cor, 0.9*f(i,k), 0.5*h_guess(i,k), &
803  0.25*ds_dsp1(i,k)*h_guess(i,k-1) )
804  endif
805 
806  ea(i,j,k) = ea(i,j,k) - dsp1_ds(i,k)*f_cor
807  eb(i,j,k) = eb(i,j,k) + f_cor
808  endif ; enddo
809  enddo
810  endif
811 
812  endif ! associated(tv%eqn_of_state))
813 
814  if (cs%id_Kd > 0) then
815  idt = gv%H_to_Z**2 / dt
816  do k=2,nz-1 ; do i=is,ie
817  if (k<kb(i)) then ; kd_here = 0.0 ; else
818  kd_here = f(i,k) * ( h(i,j,k) + ((ea(i,j,k) - eb(i,j,k-1)) + &
819  (eb(i,j,k) - ea(i,j,k+1))) ) / (i2p2dsp1_ds(i,k) * grats(i,k))
820  endif
821 
822  kd_eff(i,j,k) = max(dtkd(i,k), kd_here)*idt
823  enddo ; enddo
824  do i=is,ie
825  kd_eff(i,j,1) = dtkd(i,1)*idt
826  kd_eff(i,j,nz) = dtkd(i,nz)*idt
827  enddo
828  endif
829 
830  if (cs%id_diff_work > 0) then
831  g_2dt = 0.5 * gv%H_to_Z**2*us%L_to_Z**2 * (gv%g_Earth / dt)
832  do i=is,ie ; diff_work(i,j,1) = 0.0 ; diff_work(i,j,nz+1) = 0.0 ; enddo
833  if (associated(tv%eqn_of_state)) then
834  if (associated(fluxes%p_surf)) then
835  do i=is,ie ; pressure(i) = fluxes%p_surf(i,j) ; enddo
836  else
837  do i=is,ie ; pressure(i) = 0.0 ; enddo
838  endif
839  do k=2,nz
840  do i=is,ie ; pressure(i) = pressure(i) + (gv%g_Earth*gv%H_to_RZ)*h(i,j,k-1) ; enddo
841  do i=is,ie
842  if (k==kb(i)) then
843  t_eos(i) = 0.5*(tv%T(i,j,kmb) + tv%T(i,j,k))
844  s_eos(i) = 0.5*(tv%S(i,j,kmb) + tv%S(i,j,k))
845  else
846  t_eos(i) = 0.5*(tv%T(i,j,k-1) + tv%T(i,j,k))
847  s_eos(i) = 0.5*(tv%S(i,j,k-1) + tv%S(i,j,k))
848  endif
849  enddo
850  call calculate_density_derivs(t_eos, s_eos, pressure, drho_dt, drho_ds, &
851  tv%eqn_of_state, eosdom)
852  do i=is,ie
853  if ((k>kmb) .and. (k<kb(i))) then ; diff_work(i,j,k) = 0.0
854  else
855  if (k==kb(i)) then
856  drho = drho_dt(i) * (tv%T(i,j,k)-tv%T(i,j,kmb)) + &
857  drho_ds(i) * (tv%S(i,j,k)-tv%S(i,j,kmb))
858  else
859  drho = drho_dt(i) * (tv%T(i,j,k)-tv%T(i,j,k-1)) + &
860  drho_ds(i) * (tv%S(i,j,k)-tv%S(i,j,k-1))
861  endif
862  diff_work(i,j,k) = g_2dt * drho * &
863  (ea(i,j,k) * (h(i,j,k) + ea(i,j,k)) + &
864  eb(i,j,k-1)*(h(i,j,k-1) + eb(i,j,k-1)))
865  endif
866  enddo
867  enddo
868  else
869  do k=2,nz ; do i=is,ie
870  diff_work(i,j,k) = g_2dt * (gv%Rlay(k)-gv%Rlay(k-1)) * &
871  (ea(i,j,k) * (h(i,j,k) + ea(i,j,k)) + &
872  eb(i,j,k-1)*(h(i,j,k-1) + eb(i,j,k-1)))
873  enddo ; enddo
874  endif
875  endif
876 
877  if (present(kb_out)) then
878  do i=is,ie ; kb_out(i,j) = kb(i) ; enddo
879  endif
880 
881  enddo ! end of j loop
882 
883 ! Offer diagnostic fields for averaging.
884  if (cs%id_Kd > 0) call post_data(cs%id_Kd, kd_eff, cs%diag)
885  if (cs%id_Kd > 0) deallocate(kd_eff)
886  if (cs%id_diff_work > 0) call post_data(cs%id_diff_work, diff_work, cs%diag)
887  if (cs%id_diff_work > 0) deallocate(diff_work)
888 

◆ f_kb_to_ea_kb()

subroutine mom_entrain_diffusive::f_kb_to_ea_kb ( real, dimension(szi_(g),szk_(g)), intent(in)  h_bl,
real, dimension(szi_(g),szk_(g)), intent(in)  Sref,
real, dimension(szi_(g),szk_(g)), intent(in)  Ent_bl,
real, dimension(szi_(g)), intent(in)  I_dSkbp1,
real, dimension(szi_(g)), intent(in)  F_kb,
integer, intent(in)  kmb,
integer, intent(in)  i,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(entrain_diffusive_cs), pointer  CS,
real, dimension( g %isd: g %ied), intent(inout)  ea_kb,
real, intent(in), optional  tol_in 
)
private

Given an entrainment from below for layer kb, determine a consistent entrainment from above, such that dSkb * ea_kb = dSkbp1 * F_kb. The input value of ea_kb is both the maximum value that can be obtained and the first guess of the iterations. Ideally ea_kb should be an under-estimate.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]h_blLayer thickness, with the top interior
[in]srefThe coordinate reference potential density,
[in]ent_blThe average entrainment upward and downward
[in]i_dskbp1The inverse of the difference in reference potential density across the base of the uppermost interior layer [R-1 ~> m3 kg-1].
[in]f_kbThe entrainment from below by the uppermost interior layer [H ~> m or kg m-2]
[in]kmbThe number of mixed and buffer layers.
[in]iThe i-index to work on
csThis module's control structure.
[in,out]ea_kbThe entrainment from above by the layer below the buffer layer (i.e. layer kb) [H ~> m or kg m-2].
[in]tol_inA tolerance for the iterative determination of the entrainment [H ~> m or kg m-2].

Definition at line 1441 of file MOM_entrain_diffusive.F90.

1441  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1442  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
1443  real, dimension(SZI_(G),SZK_(G)), &
1444  intent(in) :: h_bl !< Layer thickness, with the top interior
1445  !! layer at k-index kmb+1 [H ~> m or kg m-2].
1446  real, dimension(SZI_(G),SZK_(G)), &
1447  intent(in) :: Sref !< The coordinate reference potential density,
1448  !! with the value of the topmost interior layer
1449  !! at index kmb+1 [R ~> kg m-3].
1450  real, dimension(SZI_(G),SZK_(G)), &
1451  intent(in) :: Ent_bl !< The average entrainment upward and downward
1452  !! across each interface around the buffer layers,
1453  !! [H ~> m or kg m-2].
1454  real, dimension(SZI_(G)), intent(in) :: I_dSkbp1 !< The inverse of the difference in reference
1455  !! potential density across the base of the
1456  !! uppermost interior layer [R-1 ~> m3 kg-1].
1457  real, dimension(SZI_(G)), intent(in) :: F_kb !< The entrainment from below by the
1458  !! uppermost interior layer [H ~> m or kg m-2]
1459  integer, intent(in) :: kmb !< The number of mixed and buffer layers.
1460  integer, intent(in) :: i !< The i-index to work on
1461  type(entrain_diffusive_CS), pointer :: CS !< This module's control structure.
1462  real, dimension(SZI_(G)), intent(inout) :: ea_kb !< The entrainment from above by the layer below
1463  !! the buffer layer (i.e. layer kb) [H ~> m or kg m-2].
1464  real, optional, intent(in) :: tol_in !< A tolerance for the iterative determination
1465  !! of the entrainment [H ~> m or kg m-2].
1466 
1467  real :: max_ea, min_ea
1468  real :: err, err_min, err_max
1469  real :: derr_dea
1470  real :: val, tolerance, tol1
1471  real :: ea_prev
1472  real :: dS_kbp1
1473  logical :: bisect_next, Newton
1474  real, dimension(SZI_(G)) :: dS_kb
1475  real, dimension(SZI_(G)) :: maxF, ent_maxF, zeros
1476  real, dimension(SZI_(G)) :: ddSkb_dE
1477  integer :: it
1478  integer, parameter :: MAXIT = 30
1479 
1480  ds_kbp1 = sref(i,kmb+2) - sref(i,kmb+1)
1481  max_ea = ea_kb(i) ; min_ea = 0.0
1482  val = ds_kbp1 * f_kb(i)
1483  err_min = -val
1484 
1485  tolerance = cs%Tolerance_Ent
1486  if (present(tol_in)) tolerance = tol_in
1487  bisect_next = .true.
1488 
1489  call determine_dskb(h_bl, sref, ent_bl, ea_kb, i, i, kmb, g, gv, .true., &
1490  ds_kb, ddskb_de)
1491 
1492  err = ds_kb(i) * ea_kb(i) - val
1493  derr_dea = ds_kb(i) + ddskb_de(i) * ea_kb(i)
1494  ! Return if Newton's method on the first guess would give a tolerably small
1495  ! change in the value of ea_kb.
1496  if ((err <= 0.0) .and. (abs(err) <= tolerance*abs(derr_dea))) return
1497 
1498  if (err == 0.0) then ; return ! The exact solution on the first guess...
1499  elseif (err > 0.0) then ! The root is properly bracketed.
1500  max_ea = ea_kb(i) ; err_max = err
1501  ! Use Newton's method (if it stays bounded) or the false position method
1502  ! to find the next value.
1503  if ((derr_dea > 0.0) .and. (derr_dea*(ea_kb(i) - min_ea) > err) .and. &
1504  (derr_dea*(max_ea - ea_kb(i)) > -1.0*err)) then
1505  ea_kb(i) = ea_kb(i) - err / derr_dea
1506  else ! Use the bisection for the next guess.
1507  ea_kb(i) = 0.5*(max_ea+min_ea)
1508  endif
1509  else
1510  ! Try to bracket the root first. If unable to bracket the root, return
1511  ! the maximum.
1512  zeros(i) = 0.0
1513  call find_maxf_kb(h_bl, sref, ent_bl, i_dskbp1, zeros, ea_kb, &
1514  kmb, i, i, g, gv, cs, maxf, ent_maxf, f_thresh = f_kb)
1515  err_max = ds_kbp1 * maxf(i) - val
1516  ! If err_max is negative, there is no good solution, so use the maximum
1517  ! value of F in the valid range.
1518  if (err_max <= 0.0) then
1519  ea_kb(i) = ent_maxf(i) ; return
1520  else
1521  max_ea = ent_maxf(i)
1522  ea_kb(i) = 0.5*(max_ea+min_ea) ! Use bisection for the next guess.
1523  endif
1524  endif
1525 
1526  ! Exit if the range between max_ea and min_ea already acceptable.
1527  ! if (abs(max_ea - min_ea) < 0.1*tolerance) return
1528 
1529  do it = 1, maxit
1530  call determine_dskb(h_bl, sref, ent_bl, ea_kb, i, i, kmb, g, gv, .true., &
1531  ds_kb, ddskb_de)
1532 
1533  err = ds_kb(i) * ea_kb(i) - val
1534  derr_dea = ds_kb(i) + ddskb_de(i) * ea_kb(i)
1535 
1536  ea_prev = ea_kb(i)
1537  ! Use Newton's method or the false position method to find the next value.
1538  newton = .false.
1539  if (err > 0.0) then
1540  max_ea = ea_kb(i) ; err_max = err
1541  if ((derr_dea > 0.0) .and. (derr_dea*(ea_kb(i)-min_ea) > err)) newton = .true.
1542  else
1543  min_ea = ea_kb(i) ; err_min = err
1544  if ((derr_dea > 0.0) .and. (derr_dea*(ea_kb(i)-max_ea) < err)) newton = .true.
1545  endif
1546 
1547  if (newton) then
1548  ea_kb(i) = ea_kb(i) - err / derr_dea
1549  elseif (bisect_next) then ! Use bisection to reduce the range.
1550  ea_kb(i) = 0.5*(max_ea+min_ea)
1551  bisect_next = .false.
1552  else ! Use the false-position method for the next guess.
1553  ea_kb(i) = min_ea + (max_ea-min_ea) * (err_min/(err_min - err_max))
1554  bisect_next = .true.
1555  endif
1556 
1557  tol1 = tolerance ; if (err > 0.0) tol1 = 0.099*tolerance
1558  if (ds_kb(i) <= ds_kbp1) then
1559  if (abs(ea_kb(i) - ea_prev) <= tol1) return
1560  else
1561  if (ds_kbp1*abs(ea_kb(i) - ea_prev) <= ds_kb(i)*tol1) return
1562  endif
1563  enddo
1564 

◆ f_to_ent()

subroutine mom_entrain_diffusive::f_to_ent ( real, dimension(szi_(g),szk_(g)), intent(in)  F,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
integer, dimension(szi_(g)), intent(in)  kb,
integer, intent(in)  kmb,
integer, intent(in)  j,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(entrain_diffusive_cs), intent(in)  CS,
real, dimension( g %isd: g %ied, g %ke), intent(in)  dsp1_ds,
real, dimension( g %isd: g %ied), intent(in)  eakb,
real, dimension( g %isd: g %ied, g %ke), intent(in)  Ent_bl,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  ea,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  eb,
logical, dimension( g %isd: g %ied), intent(in), optional  do_i_in 
)
private

This subroutine calculates the actual entrainments (ea and eb) and the amount of surface forcing that is applied to each layer if there is no bulk mixed layer.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]fThe density flux through a layer within a time step divided by the density difference across the interface below the layer [H ~> m or kg m-2].
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]kbThe index of the lightest layer denser than the deepest buffer layer.
[in]kmbThe number of mixed and buffer layers.
[in]jThe meridional index upon which to work.
[in]csThis module's control structure.
[in]dsp1_dsThe ratio of coordinate variable differences across the interfaces below a layer over the difference across the interface above the layer.
[in]eakbThe entrainment from above by the layer below the buffer layer [H ~> m or kg m-2].
[in]ent_blThe average entrainment upward and downward across each interface around the buffer layers [H ~> m or kg m-2].
[in,out]eaThe amount of fluid entrained from the layer
[in,out]ebThe amount of fluid entrained from the layer
[in]do_i_inIndicates which i-points to work on.

Definition at line 895 of file MOM_entrain_diffusive.F90.

895  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
896  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
897  real, dimension(SZI_(G),SZK_(G)), intent(in) :: F !< The density flux through a layer within
898  !! a time step divided by the density
899  !! difference across the interface below
900  !! the layer [H ~> m or kg m-2].
901  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
902  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
903  integer, dimension(SZI_(G)), intent(in) :: kb !< The index of the lightest layer denser than
904  !! the deepest buffer layer.
905  integer, intent(in) :: kmb !< The number of mixed and buffer layers.
906  integer, intent(in) :: j !< The meridional index upon which to work.
907  type(entrain_diffusive_CS), intent(in) :: CS !< This module's control structure.
908  real, dimension(SZI_(G),SZK_(G)), intent(in) :: dsp1_ds !< The ratio of coordinate variable
909  !! differences across the interfaces below
910  !! a layer over the difference across the
911  !! interface above the layer.
912  real, dimension(SZI_(G)), intent(in) :: eakb !< The entrainment from above by the layer
913  !! below the buffer layer [H ~> m or kg m-2].
914  real, dimension(SZI_(G),SZK_(G)), intent(in) :: Ent_bl !< The average entrainment upward and
915  !! downward across each interface around
916  !! the buffer layers [H ~> m or kg m-2].
917  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
918  intent(inout) :: ea !< The amount of fluid entrained from the layer
919  !! above within this time step [H ~> m or kg m-2].
920  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
921  intent(inout) :: eb !< The amount of fluid entrained from the layer
922  !! below within this time step [H ~> m or kg m-2].
923  logical, dimension(SZI_(G)), &
924  optional, intent(in) :: do_i_in !< Indicates which i-points to work on.
925 ! This subroutine calculates the actual entrainments (ea and eb) and the
926 ! amount of surface forcing that is applied to each layer if there is no bulk
927 ! mixed layer.
928 
929  real :: h1 ! The thickness in excess of the minimum that will remain
930  ! after exchange with the layer below [H ~> m or kg m-2].
931  logical :: do_i(SZI_(G))
932  integer :: i, k, is, ie, nz
933 
934  is = g%isc ; ie = g%iec ; nz = g%ke
935 
936  if (present(do_i_in)) then
937  do i=is,ie ; do_i(i) = do_i_in(i) ; enddo
938  do i=g%isc,g%iec ; if (do_i(i)) then
939  is = i ; exit
940  endif ; enddo
941  do i=g%iec,g%isc,-1 ; if (do_i(i)) then
942  ie = i ; exit
943  endif ; enddo
944  else
945  do i=is,ie ; do_i(i) = .true. ; enddo
946  endif
947 
948  do i=is,ie
949  ea(i,j,nz) = 0.0 ; eb(i,j,nz) = 0.0
950  enddo
951  if (cs%bulkmixedlayer) then
952  do i=is,ie
953  eb(i,j,kmb) = max(2.0*ent_bl(i,kmb+1) - eakb(i), 0.0)
954  enddo
955  do k=nz-1,kmb+1,-1 ; do i=is,ie ; if (do_i(i)) then
956  if (k > kb(i)) then
957  ! With a bulk mixed layer, surface buoyancy fluxes are applied
958  ! elsewhere, so F should always be nonnegative.
959  ea(i,j,k) = dsp1_ds(i,k)*f(i,k)
960  eb(i,j,k) = f(i,k)
961  elseif (k == kb(i)) then
962  ea(i,j,k) = eakb(i)
963  eb(i,j,k) = f(i,k)
964  elseif (k == kb(i)-1) then
965  ea(i,j,k) = ea(i,j,k+1)
966  eb(i,j,k) = eb(i,j,kmb)
967  else
968  ea(i,j,k) = ea(i,j,k+1)
969  ! Add the entrainment of the thin interior layers to eb going
970  ! up into the buffer layer.
971  eb(i,j,k) = eb(i,j,k+1) + max(0.0, h(i,j,k+1) - gv%Angstrom_H)
972  endif
973  endif ; enddo ; enddo
974  k = kmb
975  do i=is,ie ; if (do_i(i)) then
976  ! Adjust the previously calculated entrainment from below by the deepest
977  ! buffer layer to account for entrainment of thin interior layers .
978  if (kb(i) > kmb+1) &
979  eb(i,j,k) = eb(i,j,k+1) + max(0.0, h(i,j,k+1) - gv%Angstrom_H)
980 
981  ! Determine the entrainment from above for each buffer layer.
982  h1 = (h(i,j,k) - gv%Angstrom_H) + (eb(i,j,k) - ea(i,j,k+1))
983  ea(i,j,k) = max(ent_bl(i,k), ent_bl(i,k)-0.5*h1, -h1)
984  endif ; enddo
985  do k=kmb-1,2,-1 ; do i=is,ie ; if (do_i(i)) then
986  ! Determine the entrainment from below for each buffer layer.
987  eb(i,j,k) = max(2.0*ent_bl(i,k+1) - ea(i,j,k+1), 0.0)
988 
989  ! Determine the entrainment from above for each buffer layer.
990  h1 = (h(i,j,k) - gv%Angstrom_H) + (eb(i,j,k) - ea(i,j,k+1))
991  ea(i,j,k) = max(ent_bl(i,k), ent_bl(i,k)-0.5*h1, -h1)
992 ! if (h1 >= 0.0) then ; ea(i,j,k) = Ent_bl(i,K)
993 ! elseif (Ent_bl(i,K)+0.5*h1 >= 0.0) then ; ea(i,j,k) = Ent_bl(i,K)-0.5*h1
994 ! else ; ea(i,j,k) = -h1 ; endif
995  endif ; enddo ; enddo
996  do i=is,ie ; if (do_i(i)) then
997  eb(i,j,1) = max(2.0*ent_bl(i,2) - ea(i,j,2), 0.0)
998  ea(i,j,1) = 0.0
999  endif ; enddo
1000  else ! not BULKMIXEDLAYER
1001  ! Calculate the entrainment by each layer from above and below.
1002  ! Entrainment is always positive, but F may be negative due to
1003  ! surface buoyancy fluxes.
1004  do i=is,ie
1005  ea(i,j,1) = 0.0 ; eb(i,j,1) = max(f(i,1),0.0)
1006  ea(i,j,2) = dsp1_ds(i,2)*f(i,2) - min(f(i,1),0.0)
1007  enddo
1008 
1009  do k=2,nz-1 ; do i=is,ie
1010  eb(i,j,k) = max(f(i,k),0.0)
1011  ea(i,j,k+1) = dsp1_ds(i,k+1)*f(i,k+1) - (f(i,k)-eb(i,j,k))
1012  if (ea(i,j,k+1) < 0.0) then
1013  eb(i,j,k) = eb(i,j,k) - ea(i,j,k+1)
1014  ea(i,j,k+1) = 0.0
1015  endif
1016  enddo ; enddo
1017  endif ! end BULKMIXEDLAYER

◆ find_maxf_kb()

subroutine mom_entrain_diffusive::find_maxf_kb ( real, dimension(szi_(g),szk_(g)), intent(in)  h_bl,
real, dimension(szi_(g),szk_(g)), intent(in)  Sref,
real, dimension(szi_(g),szk_(g)), intent(in)  Ent_bl,
real, dimension(szi_(g)), intent(in)  I_dSkbp1,
real, dimension(szi_(g)), intent(in)  min_ent_in,
real, dimension(szi_(g)), intent(in)  max_ent_in,
integer, intent(in)  kmb,
integer, intent(in)  is,
integer, intent(in)  ie,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(entrain_diffusive_cs), pointer  CS,
real, dimension( g %isd: g %ied), intent(out)  maxF,
real, dimension( g %isd: g %ied), intent(out), optional  ent_maxF,
logical, dimension( g %isd: g %ied), intent(in), optional  do_i_in,
real, dimension( g %isd: g %ied), intent(out), optional  F_lim_maxent,
real, dimension( g %isd: g %ied), intent(in), optional  F_thresh 
)
private

Maximize F = ent*ds_kb*I_dSkbp1 in the range min_ent < ent < max_ent.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]h_blLayer thickness [H ~> m or kg m-2]
[in]srefReference potential density [R ~> kg m-3].
[in]ent_blThe average entrainment upward and
[in]i_dskbp1The inverse of the difference in reference potential density across the base of the uppermost interior layer [R-1 ~> m3 kg-1].
[in]min_ent_inThe minimum value of ent to search, [H ~> m or kg m-2].
[in]max_ent_inThe maximum value of ent to search, [H ~> m or kg m-2].
[in]kmbThe number of mixed and buffer layers.
[in]isThe start of the i-index range to work on.
[in]ieThe end of the i-index range to work on.
csThis module's control structure.
[out]maxfThe maximum value of F = ent*ds_kb*I_dSkbp1 found in the range min_ent < ent < max_ent [H ~> m or kg m-2].
[out]ent_maxfThe value of ent at that maximum [H ~> m or kg m-2].
[in]do_i_inA logical array indicating which columns
[out]f_lim_maxentIf present, do not apply the limit in
[in]f_threshIf F_thresh is present, return the first

Definition at line 1786 of file MOM_entrain_diffusive.F90.

1786  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1787  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1788  real, dimension(SZI_(G),SZK_(G)), &
1789  intent(in) :: h_bl !< Layer thickness [H ~> m or kg m-2]
1790  real, dimension(SZI_(G),SZK_(G)), &
1791  intent(in) :: Sref !< Reference potential density [R ~> kg m-3].
1792  real, dimension(SZI_(G),SZK_(G)), &
1793  intent(in) :: Ent_bl !< The average entrainment upward and
1794  !! downward across each interface around
1795  !! the buffer layers [H ~> m or kg m-2].
1796  real, dimension(SZI_(G)), intent(in) :: I_dSkbp1 !< The inverse of the difference in
1797  !! reference potential density across the
1798  !! base of the uppermost interior layer
1799  !! [R-1 ~> m3 kg-1].
1800  real, dimension(SZI_(G)), intent(in) :: min_ent_in !< The minimum value of ent to search,
1801  !! [H ~> m or kg m-2].
1802  real, dimension(SZI_(G)), intent(in) :: max_ent_in !< The maximum value of ent to search,
1803  !! [H ~> m or kg m-2].
1804  integer, intent(in) :: kmb !< The number of mixed and buffer layers.
1805  integer, intent(in) :: is !< The start of the i-index range to work on.
1806  integer, intent(in) :: ie !< The end of the i-index range to work on.
1807  type(entrain_diffusive_CS), pointer :: CS !< This module's control structure.
1808  real, dimension(SZI_(G)), intent(out) :: maxF !< The maximum value of F
1809  !! = ent*ds_kb*I_dSkbp1 found in the range
1810  !! min_ent < ent < max_ent [H ~> m or kg m-2].
1811  real, dimension(SZI_(G)), &
1812  optional, intent(out) :: ent_maxF !< The value of ent at that maximum [H ~> m or kg m-2].
1813  logical, dimension(SZI_(G)), &
1814  optional, intent(in) :: do_i_in !< A logical array indicating which columns
1815  !! to work on.
1816  real, dimension(SZI_(G)), &
1817  optional, intent(out) :: F_lim_maxent !< If present, do not apply the limit in
1818  !! finding the maximum value, but return the
1819  !! limited value at ent=max_ent_in in this
1820  !! array [H ~> m or kg m-2].
1821  real, dimension(SZI_(G)), &
1822  optional, intent(in) :: F_thresh !< If F_thresh is present, return the first
1823  !! value found that has F > F_thresh, or
1824  !! the maximum.
1825 
1826 ! Maximize F = ent*ds_kb*I_dSkbp1 in the range min_ent < ent < max_ent.
1827 ! ds_kb may itself be limited to positive values in determine_dSkb, which gives
1828 ! the prospect of two local maxima in the range - one at max_ent_in with that
1829 ! minimum value of ds_kb, and the other due to the unlimited (potentially
1830 ! negative) value. It is faster to find the true maximum by first finding the
1831 ! unlimited maximum and comparing it to the limited value at max_ent_in.
1832  real, dimension(SZI_(G)) :: &
1833  ent, &
1834  minent, maxent, ent_best, &
1835  F_max_ent_in, &
1836  F_maxent, F_minent, F, F_best, &
1837  dF_dent, dF_dE_max, dF_dE_min, dF_dE_best, &
1838  dS_kb, dS_kb_lim, ddSkb_dE, dS_anom_lim, &
1839  chg_prev, chg_pre_prev
1840  real :: dF_dE_mean, maxslope, minslope
1841  real :: tolerance
1842  real :: ratio_select_end
1843  real :: rat, max_chg, min_chg, chg1, chg2, chg
1844  logical, dimension(SZI_(G)) :: do_i, last_it, need_bracket, may_use_best
1845  logical :: doany, OK1, OK2, bisect, new_min_bound
1846  integer :: i, it, is1, ie1
1847  integer, parameter :: MAXIT = 20
1848 
1849  tolerance = cs%Tolerance_Ent
1850 
1851  if (present(do_i_in)) then
1852  do i=is,ie ; do_i(i) = do_i_in(i) ; enddo
1853  else
1854  do i=is,ie ; do_i(i) = .true. ; enddo
1855  endif
1856 
1857  ! The most likely value is at max_ent.
1858  call determine_dskb(h_bl, sref, ent_bl, max_ent_in, is, ie, kmb, g, gv, .false., &
1859  ds_kb, ddskb_de, ds_anom_lim=ds_anom_lim)
1860  ie1 = is-1 ; doany = .false.
1861  do i=is,ie
1862  ds_kb_lim(i) = ds_kb(i) + ds_anom_lim(i)
1863  f_max_ent_in(i) = max_ent_in(i)*ds_kb_lim(i)*i_dskbp1(i)
1864  maxent(i) = max_ent_in(i) ; minent(i) = min_ent_in(i)
1865  if ((abs(maxent(i) - minent(i)) < tolerance) .or. (.not.do_i(i))) then
1866  f_best(i) = max_ent_in(i)*ds_kb(i)*i_dskbp1(i)
1867  ent_best(i) = max_ent_in(i) ; ent(i) = max_ent_in(i)
1868  do_i(i) = .false.
1869  else
1870  f_maxent(i) = maxent(i) * ds_kb(i) * i_dskbp1(i)
1871  df_de_max(i) = (ds_kb(i) + maxent(i)*ddskb_de(i)) * i_dskbp1(i)
1872  doany = .true. ; last_it(i) = .false. ; need_bracket(i) = .true.
1873  endif
1874  enddo
1875 
1876  if (doany) then
1877  ie1 = is-1 ; do i=is,ie ; if (do_i(i)) ie1 = i ; enddo
1878  do i=ie1,is,-1 ; if (do_i(i)) is1 = i ; enddo
1879  ! Find the value of F and its derivative at min_ent.
1880  call determine_dskb(h_bl, sref, ent_bl, minent, is1, ie1, kmb, g, gv, .false., &
1881  ds_kb, ddskb_de, do_i_in = do_i)
1882  do i=is1,ie1 ; if (do_i(i)) then
1883  f_minent(i) = minent(i) * ds_kb(i) * i_dskbp1(i)
1884  df_de_min(i) = (ds_kb(i) + minent(i)*ddskb_de(i)) * i_dskbp1(i)
1885  endif ; enddo
1886 
1887  ratio_select_end = 0.9
1888  do it=1,maxit
1889  ratio_select_end = 0.5*ratio_select_end
1890  do i=is1,ie1 ; if (do_i(i)) then
1891  if (need_bracket(i)) then
1892  df_de_mean = (f_maxent(i) - f_minent(i)) / (maxent(i) - minent(i))
1893  maxslope = max(df_de_mean, df_de_min(i), df_de_max(i))
1894  minslope = min(df_de_mean, df_de_min(i), df_de_max(i))
1895  if (f_minent(i) >= f_maxent(i)) then
1896  if (df_de_min(i) > 0.0) then ; rat = 0.02 ! A small step should bracket the soln.
1897  elseif (maxslope < ratio_select_end*minslope) then
1898  ! The maximum of F is at minent.
1899  f_best(i) = f_minent(i) ; ent_best(i) = minent(i) ; rat = 0.0
1900  do_i(i) = .false.
1901  else ; rat = 0.382 ; endif ! Use the golden ratio
1902  else
1903  if (df_de_max(i) < 0.0) then ; rat = 0.98 ! A small step should bracket the soln.
1904  elseif (minslope > ratio_select_end*maxslope) then
1905  ! The maximum of F is at maxent.
1906  f_best(i) = f_maxent(i) ; ent_best(i) = maxent(i) ; rat = 1.0
1907  do_i(i) = .false.
1908  else ; rat = 0.618 ; endif ! Use the golden ratio
1909  endif
1910 
1911  if (rat >= 0.0) ent(i) = rat*maxent(i) + (1.0-rat)*minent(i)
1912  if (((maxent(i) - minent(i)) < tolerance) .or. (it==maxit)) &
1913  last_it(i) = .true.
1914  else ! The maximum is bracketed by minent, ent_best, and maxent.
1915  chg1 = 2.0*(maxent(i) - minent(i)) ; chg2 = chg1
1916  if (df_de_best(i) > 0) then
1917  max_chg = maxent(i) - ent_best(i) ; min_chg = 0.0
1918  else
1919  max_chg = 0.0 ; min_chg = minent(i) - ent_best(i) ! < 0
1920  endif
1921  if (max_chg - min_chg < 2.0*tolerance) last_it(i) = .true.
1922  if (df_de_max(i) /= df_de_best(i)) &
1923  chg1 = (maxent(i) - ent_best(i))*df_de_best(i) / &
1924  (df_de_best(i) - df_de_max(i))
1925  if (df_de_min(i) /= df_de_best(i)) &
1926  chg2 = (minent(i) - ent_best(i))*df_de_best(i) / &
1927  (df_de_best(i) - df_de_min(i))
1928  ok1 = ((chg1 < max_chg) .and. (chg1 > min_chg))
1929  ok2 = ((chg2 < max_chg) .and. (chg2 > min_chg))
1930  if (.not.(ok1 .or. ok2)) then ; bisect = .true. ; else
1931  if (ok1 .and. ok2) then ! Take the acceptable smaller change.
1932  chg = chg1 ; if (abs(chg2) < abs(chg1)) chg = chg2
1933  elseif (ok1) then ; chg = chg1
1934  else ; chg = chg2 ; endif
1935  if (abs(chg) > 0.5*abs(chg_pre_prev(i))) then ; bisect = .true.
1936  else ; bisect = .false. ; endif
1937  endif
1938  chg_pre_prev(i) = chg_prev(i)
1939  if (bisect) then
1940  if (df_de_best(i) > 0.0) then
1941  ent(i) = 0.5*(maxent(i) + ent_best(i))
1942  chg_prev(i) = 0.5*(maxent(i) - ent_best(i))
1943  else
1944  ent(i) = 0.5*(minent(i) + ent_best(i))
1945  chg_prev(i) = 0.5*(minent(i) - ent_best(i))
1946  endif
1947  else
1948  if (abs(chg) < tolerance) chg = sign(tolerance,chg)
1949  ent(i) = ent_best(i) + chg
1950  chg_prev(i) = chg
1951  endif
1952  endif
1953  endif ; enddo
1954 
1955  if (mod(it,3) == 0) then ! Re-determine the loop bounds.
1956  ie1 = is-1 ; do i=is1,ie ; if (do_i(i)) ie1 = i ; enddo
1957  do i=ie1,is,-1 ; if (do_i(i)) is1 = i ; enddo
1958  endif
1959 
1960  call determine_dskb(h_bl, sref, ent_bl, ent, is1, ie1, kmb, g, gv, .false., &
1961  ds_kb, ddskb_de, do_i_in = do_i)
1962  do i=is1,ie1 ; if (do_i(i)) then
1963  f(i) = ent(i)*ds_kb(i)*i_dskbp1(i)
1964  df_dent(i) = (ds_kb(i) + ent(i)*ddskb_de(i)) * i_dskbp1(i)
1965  endif ; enddo
1966 
1967  if (present(f_thresh)) then ; do i=is1,ie1 ; if (do_i(i)) then
1968  if (f(i) >= f_thresh(i)) then
1969  f_best(i) = f(i) ; ent_best(i) = ent(i) ; do_i(i) = .false.
1970  endif
1971  endif ; enddo ; endif
1972 
1973  doany = .false.
1974  do i=is1,ie1 ; if (do_i(i)) then
1975  if (.not.last_it(i)) doany = .true.
1976  if (last_it(i)) then
1977  if (need_bracket(i)) then
1978  if ((f(i) > f_maxent(i)) .and. (f(i) > f_minent(i))) then
1979  f_best(i) = f(i) ; ent_best(i) = ent(i)
1980  elseif (f_maxent(i) > f_minent(i)) then
1981  f_best(i) = f_maxent(i) ; ent_best(i) = maxent(i)
1982  else
1983  f_best(i) = f_minent(i) ; ent_best(i) = minent(i)
1984  endif
1985  elseif (f(i) > f_best(i)) then
1986  f_best(i) = f(i) ; ent_best(i) = ent(i)
1987  endif
1988  do_i(i) = .false.
1989  elseif (need_bracket(i)) then
1990  if ((f(i) > f_maxent(i)) .and. (f(i) > f_minent(i))) then
1991  need_bracket(i) = .false. ! The maximum is now bracketed.
1992  chg_prev(i) = (maxent(i) - minent(i))
1993  chg_pre_prev(i) = 2.0*chg_prev(i)
1994  ent_best(i) = ent(i) ; f_best(i) = f(i) ; df_de_best(i) = df_dent(i)
1995  elseif ((f(i) <= f_maxent(i)) .and. (f(i) > f_minent(i))) then
1996  new_min_bound = .true. ! We have a new minimum bound.
1997  elseif ((f(i) <= f_maxent(i)) .and. (f(i) > f_minent(i))) then
1998  new_min_bound = .false. ! We have a new maximum bound.
1999  else ! This case would bracket a minimum. Wierd.
2000  ! Unless the derivative indicates that there is a maximum near the
2001  ! lower bound, try keeping the end with the larger value of F
2002  ! in a tie keep the minimum as the answer here will be compared
2003  ! with the maximum input value later.
2004  new_min_bound = .true.
2005  if (df_de_min(i) > 0.0 .or. (f_minent(i) >= f_maxent(i))) &
2006  new_min_bound = .false.
2007  endif
2008  if (need_bracket(i)) then ! Still not bracketed.
2009  if (new_min_bound) then
2010  minent(i) = ent(i) ; f_minent(i) = f(i) ; df_de_min(i) = df_dent(i)
2011  else
2012  maxent(i) = ent(i) ; f_maxent(i) = f(i) ; df_de_max(i) = df_dent(i)
2013  endif
2014  endif
2015  else ! The root was previously bracketed.
2016  if (f(i) >= f_best(i)) then ! There is a new maximum.
2017  if (ent(i) > ent_best(i)) then ! Replace minent with ent_prev.
2018  minent(i) = ent_best(i) ; f_minent(i) = f_best(i) ; df_de_min(i) = df_de_best(i)
2019  else ! Replace maxent with ent_best.
2020  maxent(i) = ent_best(i) ; f_maxent(i) = f_best(i) ; df_de_max(i) = df_de_best(i)
2021  endif
2022  ent_best(i) = ent(i) ; f_best(i) = f(i) ; df_de_best(i) = df_dent(i)
2023  else
2024  if (ent(i) < ent_best(i)) then ! Replace the minent with ent.
2025  minent(i) = ent(i) ; f_minent(i) = f(i) ; df_de_min(i) = df_dent(i)
2026  else ! Replace maxent with ent_prev.
2027  maxent(i) = ent(i) ; f_maxent(i) = f(i) ; df_de_max(i) = df_dent(i)
2028  endif
2029  endif
2030  if ((maxent(i) - minent(i)) <= tolerance) do_i(i) = .false. ! Done.
2031  endif ! need_bracket.
2032  endif ; enddo
2033  if (.not.doany) exit
2034  enddo
2035  endif
2036 
2037  if (present(f_lim_maxent)) then
2038  ! Return the unlimited maximum in maxF, and the limited value of F at maxent.
2039  do i=is,ie
2040  maxf(i) = f_best(i)
2041  f_lim_maxent(i) = f_max_ent_in(i)
2042  if (present(ent_maxf)) ent_maxf(i) = ent_best(i)
2043  enddo
2044  else
2045  ! Now compare the two? potential maxima using the limited value of dF_kb.
2046  doany = .false.
2047  do i=is,ie
2048  may_use_best(i) = (ent_best(i) /= max_ent_in(i))
2049  if (may_use_best(i)) doany = .true.
2050  enddo
2051  if (doany) then
2052  ! For efficiency, could save previous value of dS_anom_lim_best?
2053  call determine_dskb(h_bl, sref, ent_bl, ent_best, is, ie, kmb, g, gv, .true., &
2054  ds_kb_lim)
2055  do i=is,ie
2056  f_best(i) = ent_best(i)*ds_kb_lim(i)*i_dskbp1(i)
2057  ! The second test seems necessary because of roundoff differences that
2058  ! can arise during compilation.
2059  if ((f_best(i) > f_max_ent_in(i)) .and. (may_use_best(i))) then
2060  maxf(i) = f_best(i)
2061  if (present(ent_maxf)) ent_maxf(i) = ent_best(i)
2062  else
2063  maxf(i) = f_max_ent_in(i)
2064  if (present(ent_maxf)) ent_maxf(i) = max_ent_in(i)
2065  endif
2066  enddo
2067  else
2068  ! All of the maxima are at the maximum entrainment.
2069  do i=is,ie ; maxf(i) = f_max_ent_in(i) ; enddo
2070  if (present(ent_maxf)) then
2071  do i=is,ie ; ent_maxf(i) = max_ent_in(i) ; enddo
2072  endif
2073  endif
2074  endif
2075 

◆ set_ent_bl()

subroutine mom_entrain_diffusive::set_ent_bl ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
real, dimension( g %isd: g %ied, g %ke+1), intent(in)  dtKd_int,
type(thermo_var_ptrs), intent(in)  tv,
integer, dimension( g %isd: g %ied), intent(inout)  kb,
integer, intent(in)  kmb,
logical, dimension( g %isd: g %ied), intent(in)  do_i,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(entrain_diffusive_cs), pointer  CS,
integer, intent(in)  j,
real, dimension(szi_(g),szk_(g)+1), intent(out)  Ent_bl,
real, dimension(szi_(g),szk_(g)), intent(out)  Sref,
real, dimension(szi_(g),szk_(g)), intent(out)  h_bl 
)
private

This subroutine sets the average entrainment across each of the interfaces between buffer layers within a timestep. It also causes thin and relatively light interior layers to be entrained by the deepest buffer layer. Also find the initial coordinate potential densities (Sref) of each layer.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]dtkd_intThe diapycnal diffusivity across
[in]tvA structure containing pointers to any available thermodynamic fields. Absent fields have NULL ptrs.
[in,out]kbThe index of the lightest layer denser than the buffer layer or 1 if there is no buffer layer.
[in]kmbThe number of mixed and buffer layers.
[in]do_iA logical variable indicating which i-points to work on.
[in]usA dimensional unit scaling type
csThis module's control structure.
[in]jThe meridional index upon which to work.
[out]ent_blThe average entrainment upward and
[out]srefThe coordinate potential density minus 1000 for each layer [R ~> kg m-3].
[out]h_blThe thickness of each layer [H ~> m or kg m-2].

Definition at line 1025 of file MOM_entrain_diffusive.F90.

1025  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1026  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1027  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1028  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1029  real, dimension(SZI_(G),SZK_(G)+1), &
1030  intent(in) :: dtKd_int !< The diapycnal diffusivity across
1031  !! each interface times the time step
1032  !! [H2 ~> m2 or kg2 m-4].
1033  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers to any
1034  !! available thermodynamic fields. Absent
1035  !! fields have NULL ptrs.
1036  integer, dimension(SZI_(G)), intent(inout) :: kb !< The index of the lightest layer denser
1037  !! than the buffer layer or 1 if there is
1038  !! no buffer layer.
1039  integer, intent(in) :: kmb !< The number of mixed and buffer layers.
1040  logical, dimension(SZI_(G)), intent(in) :: do_i !< A logical variable indicating which
1041  !! i-points to work on.
1042  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1043  type(entrain_diffusive_CS), pointer :: CS !< This module's control structure.
1044  integer, intent(in) :: j !< The meridional index upon which to work.
1045  real, dimension(SZI_(G),SZK_(G)+1), &
1046  intent(out) :: Ent_bl !< The average entrainment upward and
1047  !! downward across each interface around
1048  !! the buffer layers [H ~> m or kg m-2].
1049  real, dimension(SZI_(G),SZK_(G)), intent(out) :: Sref !< The coordinate potential density minus
1050  !! 1000 for each layer [R ~> kg m-3].
1051  real, dimension(SZI_(G),SZK_(G)), intent(out) :: h_bl !< The thickness of each layer [H ~> m or kg m-2].
1052 
1053 ! This subroutine sets the average entrainment across each of the interfaces
1054 ! between buffer layers within a timestep. It also causes thin and relatively
1055 ! light interior layers to be entrained by the deepest buffer layer.
1056 ! Also find the initial coordinate potential densities (Sref) of each layer.
1057 ! Does there need to be limiting when the layers below are all thin?
1058 
1059  ! Local variables
1060  real, dimension(SZI_(G)) :: &
1061  b1, d1, & ! Variables used by the tridiagonal solver [H-1 ~> m-1 or m2 kg-1] and [nondim].
1062  Rcv, & ! Value of the coordinate variable (potential density)
1063  ! based on the simulated T and S and P_Ref [R ~> kg m-3].
1064  pres, & ! Reference pressure (P_Ref) [R L2 T-2 ~> Pa].
1065  frac_rem, & ! The fraction of the diffusion remaining [nondim].
1066  h_interior ! The interior thickness available for entrainment [H ~> m or kg m-2].
1067  real, dimension(SZI_(G), SZK_(G)) :: &
1068  S_est ! An estimate of the coordinate potential density - 1000 after
1069  ! entrainment for each layer [R ~> kg m-3].
1070  real :: max_ent ! The maximum possible entrainment [H ~> m or kg m-2].
1071  real :: dh ! An available thickness [H ~> m or kg m-2].
1072  real :: Kd_x_dt ! The diffusion that remains after thin layers are
1073  ! entrained [H2 ~> m2 or kg2 m-4].
1074  real :: h_neglect ! A thickness that is so small it is usually lost
1075  ! in roundoff and can be neglected [H ~> m or kg m-2].
1076  integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
1077  integer :: i, k, is, ie, nz
1078  is = g%isc ; ie = g%iec ; nz = g%ke
1079 
1080 ! max_ent = 1.0e14*GV%Angstrom_H ! This is set to avoid roundoff problems.
1081  max_ent = 1.0e4*gv%m_to_H
1082  h_neglect = gv%H_subroundoff
1083 
1084  do i=is,ie ; pres(i) = tv%P_Ref ; enddo
1085  eosdom(:) = eos_domain(g%HI)
1086  do k=1,kmb
1087  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), pres, rcv, tv%eqn_of_state, eosdom)
1088  do i=is,ie
1089  h_bl(i,k) = h(i,j,k) + h_neglect
1090  sref(i,k) = rcv(i) - cs%Rho_sig_off
1091  enddo
1092  enddo
1093 
1094  do i=is,ie
1095  h_interior(i) = 0.0 ; ent_bl(i,1) = 0.0
1096 ! if (kb(i) > nz) Ent_bl(i,Kmb+1) = 0.0
1097  enddo
1098 
1099  do k=2,kmb ; do i=is,ie
1100  if (do_i(i)) then
1101  ent_bl(i,k) = min(2.0 * dtkd_int(i,k) / (h(i,j,k-1) + h(i,j,k) + h_neglect), &
1102  max_ent)
1103  else ; ent_bl(i,k) = 0.0 ; endif
1104  enddo ; enddo
1105 
1106  ! Determine the coordinate density of the bottommost buffer layer if there
1107  ! is no entrainment from the layers below. This is a partial solver, based
1108  ! on the first pass of a tridiagonal solver, as the values in the upper buffer
1109  ! layers are not needed.
1110 
1111  do i=is,ie
1112  b1(i) = 1.0 / (h_bl(i,1) + ent_bl(i,2))
1113  d1(i) = h_bl(i,1) * b1(i) ! = 1.0 - Ent_bl(i,2)*b1(i)
1114  s_est(i,1) = (h_bl(i,1)*sref(i,1)) * b1(i)
1115  enddo
1116  do k=2,kmb-1 ; do i=is,ie
1117  b1(i) = 1.0 / ((h_bl(i,k) + ent_bl(i,k+1)) + d1(i)*ent_bl(i,k))
1118  d1(i) = (h_bl(i,k) + d1(i)*ent_bl(i,k)) * b1(i) ! = 1.0 - Ent_bl(i,K+1)*b1(i)
1119  s_est(i,k) = (h_bl(i,k)*sref(i,k) + ent_bl(i,k)*s_est(i,k-1)) * b1(i)
1120  enddo ; enddo
1121  do i=is,ie
1122  s_est(i,kmb) = (h_bl(i,kmb)*sref(i,kmb) + ent_bl(i,kmb)*s_est(i,kmb-1)) / &
1123  (h_bl(i,kmb) + d1(i)*ent_bl(i,kmb))
1124  frac_rem(i) = 1.0
1125  enddo
1126 
1127  ! Entrain any thin interior layers that are lighter (in the coordinate
1128  ! potential density) than the deepest buffer layer will be, and adjust kb.
1129  do i=is,ie ; kb(i) = nz+1 ; if (do_i(i)) kb(i) = kmb+1 ; enddo
1130 
1131  do k=kmb+1,nz ; do i=is,ie ; if (do_i(i)) then
1132  if ((k == kb(i)) .and. (s_est(i,kmb) > (gv%Rlay(k) - cs%Rho_sig_off))) then
1133  if (4.0*dtkd_int(i,kmb+1)*frac_rem(i) > &
1134  (h_bl(i,kmb) + h(i,j,k)) * (h(i,j,k) - gv%Angstrom_H)) then
1135  ! Entrain this layer into the buffer layer and move kb down.
1136  dh = max((h(i,j,k) - gv%Angstrom_H), 0.0)
1137  if (dh > 0.0) then
1138  frac_rem(i) = frac_rem(i) - ((h_bl(i,kmb) + h(i,j,k)) * dh) / &
1139  (4.0*dtkd_int(i,kmb+1))
1140  sref(i,kmb) = (h_bl(i,kmb)*sref(i,kmb) + dh*(gv%Rlay(k)-cs%Rho_sig_off)) / &
1141  (h_bl(i,kmb) + dh)
1142  h_bl(i,kmb) = h_bl(i,kmb) + dh
1143  s_est(i,kmb) = (h_bl(i,kmb)*sref(i,kmb) + ent_bl(i,kmb)*s_est(i,kmb-1)) / &
1144  (h_bl(i,kmb) + d1(i)*ent_bl(i,kmb))
1145  endif
1146  kb(i) = kb(i) + 1
1147  endif
1148  endif
1149  endif ; enddo ; enddo
1150 
1151  ! This is where variables are be set up with a different vertical grid
1152  ! in which the (newly?) massless layers are taken out.
1153  do k=nz,kmb+1,-1 ; do i=is,ie
1154  if (k >= kb(i)) h_interior(i) = h_interior(i) + (h(i,j,k)-gv%Angstrom_H)
1155  if (k==kb(i)) then
1156  h_bl(i,kmb+1) = h(i,j,k) ; sref(i,kmb+1) = gv%Rlay(k) - cs%Rho_sig_off
1157  elseif (k==kb(i)+1) then
1158  h_bl(i,kmb+2) = h(i,j,k) ; sref(i,kmb+2) = gv%Rlay(k) - cs%Rho_sig_off
1159  endif
1160  enddo ; enddo
1161  do i=is,ie ; if (kb(i) >= nz) then
1162  h_bl(i,kmb+1) = h(i,j,nz)
1163  sref(i,kmb+1) = gv%Rlay(nz) - cs%Rho_sig_off
1164  h_bl(i,kmb+2) = gv%Angstrom_H
1165  sref(i,kmb+2) = sref(i,kmb+1) + (gv%Rlay(nz) - gv%Rlay(nz-1))
1166  endif ; enddo
1167 
1168  ! Perhaps we should revisit the way that the average entrainment between the
1169  ! buffer layer and the interior is calculated so that it is not unduly
1170  ! limited when the layers are less than sqrt(Kd * dt) thick?
1171  do i=is,ie ; if (do_i(i)) then
1172  kd_x_dt = frac_rem(i) * dtkd_int(i,kmb+1)
1173  if ((kd_x_dt <= 0.0) .or. (h_interior(i) <= 0.0)) then
1174  ent_bl(i,kmb+1) = 0.0
1175  else
1176  ! If the combined layers are exceptionally thin, use sqrt(Kd*dt) as the
1177  ! estimate of the thickness in the denominator of the thickness diffusion.
1178  ent_bl(i,kmb+1) = min(0.5*h_interior(i), sqrt(kd_x_dt), &
1179  kd_x_dt / (0.5*(h_bl(i,kmb) + h_bl(i,kmb+1))))
1180  endif
1181  else
1182  ent_bl(i,kmb+1) = 0.0
1183  endif ; enddo
1184