MOM6
|
This module implements boundary forcing for MOM6.
The ocean is a forced-dissipative system. Forcing occurs at the boundaries, and this module mediates the various forcing terms from momentum, heat, salt, and mass. Boundary fluxes from other tracers are treated by coupling to biogeochemical models. We here present elements of how MOM6 assumes boundary fluxes are passed into the ocean.
Note that all fluxes are positive into the ocean. For surface boundary fluxes, that means fluxes are positive downward. For example, a positive shortwave flux warms the ocean.
The ocean surface exchanges momentum with the overlying atmosphere, sea ice, and land ice. The momentum is exchanged as a horizontal stress (Newtons per squared meter: N/m2) imposed on the upper ocean grid cell.
The ocean gains or loses mass through evaporation, precipitation, sea ice melt/form, and and river runoff. Positive mass fluxes add mass to the liquid ocean. The boundary mass flux units are (kilogram per square meter per sec: kg/(m2/sec)).
Over most of the ocean, there is no exchange of salt with the atmosphere. However, the liquid ocean exchanges salt with sea ice. When ice forms, it extracts salt from ice pockets and discharges the salt into the liquid ocean. The salt concentration of sea ice is therefore much lower (around 5ppt) than liquid seawater (around 30-35ppt in high latitudes).
For ocean-ice models run with a prescribed atmosphere, such as in the CORE/OMMIP simulations, it is necessary to employ a surface restoring term to the k=1 salinity equation, thus imposing a salt flux onto the ocean even outside of sea ice regimes. This salt flux is non-physical, and represents a limitation of the ocean-ice models run without an interactive atmosphere. Sometimes this salt flux is converted to an implied fresh water flux. However, doing so generally leads to changes in the sea level, unless a global normalization is provided to zero-out the net water flux. As a complement, for models with a restoring salt flux, one may choose to zero-out the net salt entering the ocean. There are pros/cons of each approach.
There are many terms that contribute to boundary-related heating of the k=1 surface model grid cell. We here outline details of this heat, with each term having units W/m2.
The net flux of heat crossing ocean surface is stored in the diagnostic array "hfds". This array is computed as
\[ \mbox{hfds = shortwave + longwave + latent + sensible + mass transfer + frazil + restore + flux adjustments} \]
The shortwave field itself is split into two pieces:
The convergence of boundary-related heat into surface grid cell is given by the difference in the net heat entering the top of the k=1 cell and the penetrative SW leaving the bottom of the cell.
\begin{eqnarray*} Q(k=1) &=& \mbox{hfds} - \mbox{pen_SW(leaving bottom of k=1)} \\ &=& \mbox{nonpen_SW} + (\mbox{pen_SW(enter k=1)}-\mbox{pen_SW(leave k=1)}) + \mbox{LW+LAT+SENS+MASS+FRAZ+RES} \\ &=& \mbox{nonpen_SW}+ \mbox{LW+LAT+SENS+MASS+FRAZ+RES} + [\mbox{pen_SW(enter k=1)} - \mbox{pen_SW(leave k=1)}] \end{eqnarray*}
The convergence of the penetrative shortwave flux is given by \( \mbox{pen_SW (enter k)}-\mbox{pen_SW (leave k)}\). This term appears for all cells k=1,nz. It is diagnosed as "rsdoabsorb" inside module MOM6/src/parameterizations/vertical/MOM_diabatic_aux.F90
Data Types | |
type | forcing |
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by MOM. Data in this type is allocated in the module MOM_surface_forcing.F90, of which there are three: solo, coupled, and ice-shelf. Alternatively, they are allocated in MESO_surface_forcing.F90, which is a special case of solo_driver/MOM_surface_forcing.F90. More... | |
type | forcing_diags |
Structure that defines the id handles for the forcing type. More... | |
Functions/Subroutines | |
subroutine, public | extractfluxes1d (G, GV, fluxes, optics, nsw, j, dt, DepthBeforeScalingFluxes, useRiverHeatContent, useCalvingHeatContent, h, T, netMassInOut, netMassOut, net_heat, net_salt, pen_SW_bnd, tv, aggregate_FW_forcing, nonpenSW, netmassInOut_rate, net_Heat_Rate, net_salt_rate, pen_sw_bnd_Rate, skip_diags) |
This subroutine extracts fluxes from the surface fluxes type. It works on a j-row for optimization purposes. The 2d (i,j) wrapper is the next subroutine below. This routine multiplies fluxes by dt, so that the result is an accumulation of fluxes over a time step. More... | |
subroutine, public | extractfluxes2d (G, GV, fluxes, optics, nsw, dt, DepthBeforeScalingFluxes, useRiverHeatContent, useCalvingHeatContent, h, T, netMassInOut, netMassOut, net_heat, Net_salt, Pen_SW_bnd, tv, aggregate_FW_forcing) |
2d wrapper for 1d extract fluxes from surface fluxes type. This subroutine extracts fluxes from the surface fluxes type. It multiplies the fluxes by dt, so that the result is an accumulation of the fluxes over a time step. More... | |
subroutine, public | calculatebuoyancyflux1d (G, GV, fluxes, optics, h, Temp, Salt, tv, j, buoyancyFlux, netHeatMinusSW, netSalt, skip_diags) |
This routine calculates surface buoyancy flux by adding up the heat, FW & salt fluxes. These are actual fluxes, with units of stuff per time. Setting dt=1 in the call to extractFluxes routine allows us to get "stuf per time" rather than the time integrated fluxes needed in other routines that call extractFluxes. More... | |
subroutine, public | calculatebuoyancyflux2d (G, GV, fluxes, optics, h, Temp, Salt, tv, buoyancyFlux, netHeatMinusSW, netSalt, skip_diags) |
Calculates surface buoyancy flux by adding up the heat, FW and salt fluxes, for 2d arrays. This is a wrapper for calculateBuoyancyFlux1d. More... | |
subroutine, public | mom_forcing_chksum (mesg, fluxes, G, haloshift) |
Write out chksums for basic state variables. More... | |
subroutine, public | forcing_singlepointprint (fluxes, G, i, j, mesg) |
Write out values of the fluxes arrays at the i,j location. This is a debugging tool. More... | |
subroutine, public | register_forcing_type_diags (Time, diag, use_temperature, handles, use_berg_fluxes) |
Register members of the forcing type for diagnostics. More... | |
subroutine, public | forcing_accumulate (flux_tmp, fluxes, dt, G, wt2) |
Accumulate the forcing over time steps. More... | |
subroutine, public | mech_forcing_diags (fluxes, dt, G, diag, handles) |
Offer mechanical forcing fields for diagnostics for those fields registered as part of register_forcing_type_diags. More... | |
subroutine, public | forcing_diagnostics (fluxes, state, dt, G, diag, handles) |
Offer buoyancy forcing fields for diagnostics for those fields registered as part of register_forcing_type_diags. More... | |
subroutine, public | allocate_forcing_type (G, fluxes, stress, ustar, water, heat, shelf, press, iceberg) |
Conditionally allocate fields within the forcing type. More... | |
subroutine, public | deallocate_forcing_type (fluxes) |
Deallocate the forcing type. More... | |
subroutine, public mom_forcing_type::allocate_forcing_type | ( | type(ocean_grid_type), intent(in) | G, |
type(forcing), intent(inout) | fluxes, | ||
logical, intent(in), optional | stress, | ||
logical, intent(in), optional | ustar, | ||
logical, intent(in), optional | water, | ||
logical, intent(in), optional | heat, | ||
logical, intent(in), optional | shelf, | ||
logical, intent(in), optional | press, | ||
logical, intent(in), optional | iceberg | ||
) |
Conditionally allocate fields within the forcing type.
[in] | g | Ocean grid structure |
[in,out] | fluxes | Forcing fields structure |
[in] | stress | If present and true, allocate taux, tauy |
[in] | ustar | If present and true, allocate ustar |
[in] | water | If present and true, allocate water fluxes |
[in] | heat | If present and true, allocate heat fluxes |
[in] | shelf | If present and true, allocate fluxes for ice-shelf |
[in] | press | If present and true, allocate p_surf |
[in] | iceberg | If present and true, allocate fluxes for icebergs |
Definition at line 2334 of file MOM_forcing_type.F90.
References myalloc().
Referenced by mom_surface_forcing::buoyancy_forcing_allocate(), meso_surface_forcing::meso_wind_forcing(), ocean_model_mod::ocean_model_init(), scm_cvmix_tests::scm_cvmix_tests_buoyancy_forcing(), scm_cvmix_tests::scm_cvmix_tests_wind_forcing(), scm_idealized_hurricane::scm_idealized_hurricane_wind_forcing(), mom_surface_forcing::set_forcing(), user_surface_forcing::user_wind_forcing(), and mom_surface_forcing::wind_forcing_by_data_override().
subroutine, public mom_forcing_type::calculatebuoyancyflux1d | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(forcing), intent(inout) | fluxes, | ||
type(optics_type), pointer | optics, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | h, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | Temp, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | Salt, | ||
type(thermo_var_ptrs), intent(inout) | tv, | ||
integer, intent(in) | j, | ||
real, dimension( g %isd: g %ied, g %ke+1), intent(inout) | buoyancyFlux, | ||
real, dimension( g %isd: g %ied), intent(inout) | netHeatMinusSW, | ||
real, dimension( g %isd: g %ied), intent(inout) | netSalt, | ||
logical, intent(in), optional | skip_diags | ||
) |
This routine calculates surface buoyancy flux by adding up the heat, FW & salt fluxes. These are actual fluxes, with units of stuff per time. Setting dt=1 in the call to extractFluxes routine allows us to get "stuf per time" rather than the time integrated fluxes needed in other routines that call extractFluxes.
[in] | g | ocean grid |
[in] | gv | ocean vertical grid structure |
[in,out] | fluxes | surface fluxes |
optics | penetrating SW optics | |
[in] | h | layer thickness (H) |
[in] | temp | prognostic temp(deg C) |
[in] | salt | salinity (ppt) |
[in,out] | tv | thermodynamics type |
[in] | j | j-row to work on |
[in,out] | buoyancyflux | buoyancy flux (m^2/s^3) |
[in,out] | netheatminussw | surf Heat flux (K H/s) |
[in,out] | netsalt | surf salt flux (ppt H/s) |
[in] | skip_diags | If present and true, skip calculating diagnostics inside extractFluxes1d() |
Definition at line 775 of file MOM_forcing_type.F90.
References mom_eos::calculate_density_derivs(), extractfluxes1d(), and mom_shortwave_abs::sumswoverbands().
Referenced by calculatebuoyancyflux2d().
subroutine, public mom_forcing_type::calculatebuoyancyflux2d | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(forcing), intent(inout) | fluxes, | ||
type(optics_type), pointer | optics, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | h, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | Temp, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | Salt, | ||
type(thermo_var_ptrs), intent(inout) | tv, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke+1), intent(inout) | buoyancyFlux, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(inout), optional | netHeatMinusSW, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(inout), optional | netSalt, | ||
logical, intent(in), optional | skip_diags | ||
) |
Calculates surface buoyancy flux by adding up the heat, FW and salt fluxes, for 2d arrays. This is a wrapper for calculateBuoyancyFlux1d.
[in] | g | ocean grid |
[in] | gv | ocean vertical grid structure |
[in,out] | fluxes | surface fluxes |
optics | SW ocean optics | |
[in] | h | layer thickness (H) |
[in] | temp | temperature (deg C) |
[in] | salt | salinity (ppt) |
[in,out] | tv | thermodynamics type |
[in,out] | buoyancyflux | buoy flux (m^2/s^3) |
[in,out] | netheatminussw | surf temp flux (K H) |
[in,out] | netsalt | surf salt flux (ppt H) |
[in] | skip_diags | If present and true, skip calculating diagnostics inside extractFluxes1d() |
Definition at line 864 of file MOM_forcing_type.F90.
References calculatebuoyancyflux1d().
Referenced by mom_diabatic_driver::diabatic().
subroutine, public mom_forcing_type::deallocate_forcing_type | ( | type(forcing), intent(inout) | fluxes | ) |
Deallocate the forcing type.
[in,out] | fluxes | Forcing fields structure |
Definition at line 2425 of file MOM_forcing_type.F90.
subroutine, public mom_forcing_type::extractfluxes1d | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(forcing), intent(inout) | fluxes, | ||
type(optics_type), pointer | optics, | ||
integer, intent(in) | nsw, | ||
integer, intent(in) | j, | ||
real, intent(in) | dt, | ||
real, intent(in) | DepthBeforeScalingFluxes, | ||
logical, intent(in) | useRiverHeatContent, | ||
logical, intent(in) | useCalvingHeatContent, | ||
real, dimension( g %isd: g %ied, g %ke), intent(in) | h, | ||
real, dimension( g %isd: g %ied, g %ke), intent(in) | T, | ||
real, dimension( g %isd: g %ied), intent(out) | netMassInOut, | ||
real, dimension( g %isd: g %ied), intent(out) | netMassOut, | ||
real, dimension( g %isd: g %ied), intent(out) | net_heat, | ||
real, dimension( g %isd: g %ied), intent(out) | net_salt, | ||
real, dimension(:,:), intent(out) | pen_SW_bnd, | ||
type(thermo_var_ptrs), intent(inout) | tv, | ||
logical, intent(in) | aggregate_FW_forcing, | ||
real, dimension( g %isd: g %ied), intent(out), optional | nonpenSW, | ||
real, dimension( g %isd: g %ied), intent(out), optional | netmassInOut_rate, | ||
real, dimension( g %isd: g %ied), intent(out), optional | net_Heat_Rate, | ||
real, dimension( g %isd: g %ied), intent(out), optional | net_salt_rate, | ||
real, dimension(:,:), intent(out), optional | pen_sw_bnd_Rate, | ||
logical, intent(in), optional | skip_diags | ||
) |
This subroutine extracts fluxes from the surface fluxes type. It works on a j-row for optimization purposes. The 2d (i,j) wrapper is the next subroutine below. This routine multiplies fluxes by dt, so that the result is an accumulation of fluxes over a time step.
[in] | g | ocean grid structure |
[in] | gv | ocean vertical grid structure |
[in,out] | fluxes | structure containing pointers to possible forcing fields. NULL unused fields. |
optics | pointer to optics | |
[in] | nsw | number of bands of penetrating SW |
[in] | j | j-index to work on |
[in] | dt | time step in seconds |
[in] | depthbeforescalingfluxes | min ocean depth before scale away fluxes (H) |
[in] | useriverheatcontent | logical for river heat content |
[in] | usecalvingheatcontent | logical for calving heat content |
[in] | h | layer thickness (in H units) |
[in] | t | layer temperatures (deg C) |
[out] | netmassinout | net mass flux (non-Bouss) or volume flux (if Bouss) of water in/out of ocean over a time step (H units) |
[out] | netmassout | net mass flux (non-Bouss) or volume flux (if Bouss) of water leaving ocean surface over a time step (H units). netMassOut < 0 means mass leaves ocean. |
[out] | net_heat | net heat at the surface accumulated over a time step for coupler + restoring. Exclude two terms from net_heat: (1) downwelling (penetrative) SW, (2) evaporation heat content, (since do not yet know evap temperature). Units of net_heat are (K * H). |
[out] | net_salt | surface salt flux into the ocean accumulated over a time step (ppt * H) |
[out] | pen_sw_bnd | penetrating SW flux, split into bands. Units are (deg K * H) and array size nsw x G isd: G ied, where nsw=number of SW bands in pen_SW_bnd. This heat flux is not part of net_heat. |
[in,out] | tv | structure containing pointers to available thermodynamic fields. Used to keep track of the heat flux associated with net mass fluxes into the ocean. |
[in] | aggregate_fw_forcing | For determining how to aggregate forcing. |
[out] | nonpensw | non-downwelling SW; use in net_heat. Sum over SW bands when diagnosing nonpenSW. Units are (K * H). |
[out] | net_heat_rate | Optional outputs of contributions to surface |
[out] | net_salt_rate | buoyancy flux which do not include dt |
[out] | netmassinout_rate | and therefore are used to compute the rate. |
[out] | pen_sw_bnd_rate | Perhaps just a temporary fix. |
[in] | skip_diags | If present and true, skip calculating diagnostics |
Definition at line 278 of file MOM_forcing_type.F90.
References mom_error_handler::mom_error().
Referenced by mom_diabatic_aux::applyboundaryfluxesinout(), mom_bulk_mixed_layer::bulkmixedlayer(), calculatebuoyancyflux1d(), and extractfluxes2d().
subroutine, public mom_forcing_type::extractfluxes2d | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(forcing), intent(inout) | fluxes, | ||
type(optics_type), pointer | optics, | ||
integer, intent(in) | nsw, | ||
real, intent(in) | dt, | ||
real, intent(in) | DepthBeforeScalingFluxes, | ||
logical, intent(in) | useRiverHeatContent, | ||
logical, intent(in) | useCalvingHeatContent, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | h, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | T, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out) | netMassInOut, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out) | netMassOut, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out) | net_heat, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out) | Net_salt, | ||
real, dimension(:,:,:), intent(out) | Pen_SW_bnd, | ||
type(thermo_var_ptrs), intent(inout) | tv, | ||
logical, intent(in) | aggregate_FW_forcing | ||
) |
2d wrapper for 1d extract fluxes from surface fluxes type. This subroutine extracts fluxes from the surface fluxes type. It multiplies the fluxes by dt, so that the result is an accumulation of the fluxes over a time step.
[in] | g | ocean grid structure |
[in] | gv | ocean vertical grid structure |
[in,out] | fluxes | structure containing pointers to forcing. |
optics | pointer to optics | |
[in] | nsw | number of bands of penetrating SW |
[in] | dt | time step in seconds |
[in] | depthbeforescalingfluxes | min ocean depth before scale away fluxes (H) |
[in] | useriverheatcontent | logical for river heat content |
[in] | usecalvingheatcontent | logical for calving heat content |
[in] | h | layer thickness (in H units) |
[in] | t | layer temperatures (deg C) |
[out] | netmassinout | net mass flux (non-Bouss) or volume flux (if Bouss) of water in/out of ocean over a time step (H units) |
[out] | netmassout | net mass flux (non-Bouss) or volume flux (if Bouss) of water leaving ocean surface over a time step (H units). |
[out] | net_heat | net heat at the surface accumulated over a time step associated with coupler + restore. Exclude two terms from net_heat: (1) downwelling (penetrative) SW, (2) evaporation heat content, (since do not yet know temperature of evap). Units of net_heat are (K * H). |
[out] | net_salt | surface salt flux into the ocean accumulated over a time step (ppt * H) |
[out] | pen_sw_bnd | penetrating shortwave flux, split into bands. Units (deg K * H) & array size nsw x G isd: G ied, where nsw=number of SW bands in pen_SW_bnd. This heat flux is not in net_heat. |
[in,out] | tv | structure containing pointers to available thermodynamic fields. Here it is used to keep track of the heat flux associated with net mass fluxes into the ocean. |
[in] | aggregate_fw_forcing | For determining how to aggregate the forcing. |
Definition at line 716 of file MOM_forcing_type.F90.
References extractfluxes1d().
subroutine, public mom_forcing_type::forcing_accumulate | ( | type(forcing), intent(in) | flux_tmp, |
type(forcing), intent(inout) | fluxes, | ||
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(inout) | G, | ||
real, intent(out) | wt2 | ||
) |
Accumulate the forcing over time steps.
[in] | dt | The elapsed time since the last call to this subroutine, in s |
[in,out] | g | The ocean's grid structure |
Definition at line 1657 of file MOM_forcing_type.F90.
References mom_error_handler::mom_error().
subroutine, public mom_forcing_type::forcing_diagnostics | ( | type(forcing), intent(in) | fluxes, |
type(surface), intent(in) | state, | ||
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(diag_ctrl), intent(in) | diag, | ||
type(forcing_diags), intent(inout) | handles | ||
) |
Offer buoyancy forcing fields for diagnostics for those fields registered as part of register_forcing_type_diags.
[in] | fluxes | flux type |
[in] | state | ocean state |
[in] | dt | time step |
[in] | g | grid type |
[in] | diag | diagnostic regulator |
[in,out] | handles | diagnostic ids |
Definition at line 1848 of file MOM_forcing_type.F90.
References mom_spatial_means::global_area_integral(), mom_spatial_means::global_area_mean(), and mom_diag_mediator::query_averaging_enabled().
Referenced by mom_main().
subroutine, public mom_forcing_type::forcing_singlepointprint | ( | type(forcing), intent(in) | fluxes, |
type(ocean_grid_type), intent(in) | G, | ||
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
character(len=*), intent(in) | mesg | ||
) |
Write out values of the fluxes arrays at the i,j location. This is a debugging tool.
[in] | fluxes | Fluxes type |
[in] | g | Grid type |
[in] | mesg | Message |
[in] | i | i-index |
[in] | j | j-index |
Definition at line 980 of file MOM_forcing_type.F90.
References locmsg().
Referenced by mom_diabatic_aux::applyboundaryfluxesinout().
subroutine, public mom_forcing_type::mech_forcing_diags | ( | type(forcing), intent(in) | fluxes, |
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(diag_ctrl), intent(in) | diag, | ||
type(forcing_diags), intent(inout) | handles | ||
) |
Offer mechanical forcing fields for diagnostics for those fields registered as part of register_forcing_type_diags.
[in] | fluxes | fluxes type |
[in] | dt | time step |
[in] | g | grid type |
[in] | diag | diagnostic type |
[in,out] | handles | diagnostic id for diag_manager |
Definition at line 1809 of file MOM_forcing_type.F90.
References mom_diag_mediator::query_averaging_enabled().
Referenced by mom_main().
subroutine, public mom_forcing_type::mom_forcing_chksum | ( | character(len=*), intent(in) | mesg, |
type(forcing), intent(in) | fluxes, | ||
type(ocean_grid_type), intent(in) | G, | ||
integer, intent(in), optional | haloshift | ||
) |
Write out chksums for basic state variables.
[in] | mesg | message |
[in] | fluxes | fluxes type |
[in] | g | grid type |
[in] | haloshift | shift in halo |
Definition at line 899 of file MOM_forcing_type.F90.
Referenced by mom_main(), mom_ice_shelf::shelf_calc_flux(), mom::step_mom(), and mom::step_mom_thermo().
subroutine, public mom_forcing_type::register_forcing_type_diags | ( | type(time_type), intent(in) | Time, |
type(diag_ctrl), intent(inout) | diag, | ||
logical, intent(in) | use_temperature, | ||
type(forcing_diags), intent(inout) | handles, | ||
logical, intent(in), optional | use_berg_fluxes | ||
) |
Register members of the forcing type for diagnostics.
[in] | time | time type |
[in,out] | diag | diagnostic control type |
[in] | use_temperature | True if T/S are in use |
[in,out] | handles | handles for diagnostics |
[in] | use_berg_fluxes | If true, allow iceberg flux diagnostics |
Definition at line 1041 of file MOM_forcing_type.F90.