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. More... | |
type | forcing_diags |
Structure that defines the id handles for the forcing type. More... | |
type | mech_forcing |
Structure that contains pointers to the mechanical forcing at the surface 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 versions: solo, coupled, and ice-shelf. More... | |
Functions/Subroutines | |
subroutine, public | extractfluxes1d (G, GV, US, fluxes, optics, nsw, j, dt_in_T, FluxRescaleDepth, useRiverHeatContent, useCalvingHeatContent, h, T, netMassInOut, netMassOut, net_heat, net_salt, pen_SW_bnd, tv, aggregate_FW, 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, US, fluxes, optics, nsw, dt_in_T, FluxRescaleDepth, useRiverHeatContent, useCalvingHeatContent, h, T, netMassInOut, netMassOut, net_heat, Net_salt, Pen_SW_bnd, tv, aggregate_FW) |
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, US, fluxes, optics, nsw, 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, US, 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, US, haloshift) |
Write out chksums for thermodynamic fluxes. More... | |
subroutine, public | mom_mech_forcing_chksum (mesg, forces, G, US, haloshift) |
Write out chksums for the driving mechanical forces. More... | |
subroutine | mech_forcing_singlepointprint (forces, G, i, j, mesg) |
Write out values of the mechanical forcing arrays at the i,j location. This is a debugging tool. 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, US, use_temperature, handles, use_berg_fluxes) |
Register members of the forcing type for diagnostics. More... | |
subroutine, public | forcing_accumulate (flux_tmp, forces, fluxes, G, wt2) |
Accumulate the forcing over time steps, taking input from a mechanical forcing type and a temporary forcing-flux type. More... | |
subroutine, public | fluxes_accumulate (flux_tmp, fluxes, G, wt2, forces) |
Accumulate the thermodynamic fluxes over time steps. More... | |
subroutine, public | copy_common_forcing_fields (forces, fluxes, G, skip_pres) |
This subroutine copies the computational domains of common forcing fields from a mech_forcing type to a (thermodynamic) forcing type. More... | |
subroutine, public | set_derived_forcing_fields (forces, fluxes, G, US, Rho0) |
This subroutine calculates certain derived forcing fields based on information from a mech_forcing type and stores them in a (thermodynamic) forcing type. More... | |
subroutine, public | set_net_mass_forcing (fluxes, forces, G, US) |
This subroutine determines the net mass source to the ocean from a (thermodynamic) forcing type and stores it in a mech_forcing type. More... | |
subroutine, public | get_net_mass_forcing (fluxes, G, US, net_mass_src) |
This subroutine calculates determines the net mass source to the ocean from a (thermodynamic) forcing type and stores it in a provided array. More... | |
subroutine, public | copy_back_forcing_fields (fluxes, forces, G) |
This subroutine copies the computational domains of common forcing fields from a mech_forcing type to a (thermodynamic) forcing type. More... | |
subroutine, public | mech_forcing_diags (forces, dt, G, time_end, 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, sfc_state, G, US, time_end, 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, water, heat, ustar, press, shelf, iceberg, salt) |
Conditionally allocate fields within the forcing type. More... | |
subroutine, public | allocate_mech_forcing (G, forces, stress, ustar, shelf, press, iceberg) |
Conditionally allocate fields within the mechanical forcing type. More... | |
subroutine | myalloc (array, is, ie, js, je, flag) |
Allocates and zeroes-out array. More... | |
subroutine, public | deallocate_forcing_type (fluxes) |
Deallocate the forcing type. More... | |
subroutine, public | deallocate_mech_forcing (forces) |
Deallocate the mechanical 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 | water, | ||
logical, intent(in), optional | heat, | ||
logical, intent(in), optional | ustar, | ||
logical, intent(in), optional | press, | ||
logical, intent(in), optional | shelf, | ||
logical, intent(in), optional | iceberg, | ||
logical, intent(in), optional | salt | ||
) |
Conditionally allocate fields within the forcing type.
[in] | g | Ocean grid structure |
[in,out] | fluxes | A structure containing thermodynamic forcing fields |
[in] | water | If present and true, allocate water fluxes |
[in] | heat | If present and true, allocate heat fluxes |
[in] | ustar | If present and true, allocate ustar and related fields |
[in] | press | If present and true, allocate p_surf and related fields |
[in] | shelf | If present and true, allocate fluxes for ice-shelf |
[in] | iceberg | If present and true, allocate fluxes for icebergs |
[in] | salt | If present and true, allocate salt fluxes |
Definition at line 2811 of file MOM_forcing_type.F90.
References myalloc().
subroutine, public mom_forcing_type::allocate_mech_forcing | ( | type(ocean_grid_type), intent(in) | G, |
type(mech_forcing), intent(inout) | forces, | ||
logical, intent(in), optional | stress, | ||
logical, intent(in), optional | ustar, | ||
logical, intent(in), optional | shelf, | ||
logical, intent(in), optional | press, | ||
logical, intent(in), optional | iceberg | ||
) |
Conditionally allocate fields within the mechanical forcing type.
[in] | g | Ocean grid structure |
[in,out] | forces | Forcing fields structure |
[in] | stress | If present and true, allocate taux, tauy |
[in] | ustar | If present and true, allocate ustar and related fields |
[in] | shelf | If present and true, allocate forces for ice-shelf |
[in] | press | If present and true, allocate p_surf and related fields |
[in] | iceberg | If present and true, allocate forces for icebergs |
Definition at line 2879 of file MOM_forcing_type.F90.
References myalloc().
Referenced by mom_surface_forcing_nuopc::convert_iob_to_forces(), mom_surface_forcing_mct::convert_iob_to_forces(), idealized_hurricane::idealized_hurricane_wind_forcing(), neverland_surface_forcing::neverland_wind_forcing(), 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(unit_scale_type), intent(in) | US, | ||
type(forcing), intent(inout) | fluxes, | ||
type(optics_type), pointer | optics, | ||
integer, intent(in) | nsw, | ||
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] | us | A dimensional unit scaling type |
[in,out] | fluxes | surface fluxes |
optics | penetrating SW optics | |
[in] | nsw | The number of frequency bands of penetrating shortwave radiation |
[in] | h | layer thickness [H ~> m or kg m-2] |
[in] | temp | prognostic temp [degC] |
[in] | salt | salinity [ppt] |
[in,out] | tv | thermodynamics type |
[in] | j | j-row to work on |
[in,out] | buoyancyflux | buoyancy fluxes [L2 T-3 ~> m2 s-3] |
[in,out] | netheatminussw | surf Heat flux [degC H s-1 ~> degC m s-1 or degC kg m-2 s-1] |
[in,out] | netsalt | surf salt flux [ppt H s-1 ~> ppt m s-1 or ppt kg m-2 s-1] |
[in] | skip_diags | If present and true, skip calculating diagnostics inside extractFluxes1d() |
Definition at line 879 of file MOM_forcing_type.F90.
References extractfluxes1d(), mom_opacity::optics_nbands(), and mom_opacity::sumswoverbands().
Referenced by calculatebuoyancyflux2d().
subroutine, public mom_forcing_type::calculatebuoyancyflux2d | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
type(unit_scale_type), intent(in) | US, | ||
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] | us | A dimensional unit scaling type |
[in,out] | fluxes | surface fluxes |
optics | SW ocean optics | |
[in] | h | layer thickness [H ~> m or kg m-2] |
[in] | temp | temperature [degC] |
[in] | salt | salinity [ppt] |
[in,out] | tv | thermodynamics type |
[in,out] | buoyancyflux | buoyancy fluxes [L2 T-3 ~> m2 s-3] |
[in,out] | netheatminussw | surf temp flux [degC H ~> degC m or degC kg m-2] |
[in,out] | netsalt | surf salt flux [ppt H ~> ppt m or ppt kg m-2] |
[in] | skip_diags | If present and true, skip calculating diagnostics inside extractFluxes1d() |
Definition at line 976 of file MOM_forcing_type.F90.
References calculatebuoyancyflux1d(), and mom_opacity::optics_nbands().
Referenced by mom_diabatic_driver::diabatic_ale(), mom_diabatic_driver::diabatic_ale_legacy(), and mom_diabatic_driver::layered_diabatic().
subroutine, public mom_forcing_type::copy_back_forcing_fields | ( | type(forcing), intent(in) | fluxes, |
type(mech_forcing), intent(inout) | forces, | ||
type(ocean_grid_type), intent(in) | G | ||
) |
This subroutine copies the computational domains of common forcing fields from a mech_forcing type to a (thermodynamic) forcing type.
[in] | fluxes | A structure containing thermodynamic forcing fields |
[in,out] | forces | A structure with the driving mechanical forces |
[in] | g | grid type |
Definition at line 2182 of file MOM_forcing_type.F90.
subroutine, public mom_forcing_type::copy_common_forcing_fields | ( | type(mech_forcing), intent(in) | forces, |
type(forcing), intent(inout) | fluxes, | ||
type(ocean_grid_type), intent(in) | G, | ||
logical, intent(in), optional | skip_pres | ||
) |
This subroutine copies the computational domains of common forcing fields from a mech_forcing type to a (thermodynamic) forcing type.
[in] | forces | A structure with the driving mechanical forces |
[in,out] | fluxes | A structure containing thermodynamic forcing fields |
[in] | g | grid type |
[in] | skip_pres | If present and true, do not copy pressure fields. |
Definition at line 2046 of file MOM_forcing_type.F90.
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 2931 of file MOM_forcing_type.F90.
subroutine, public mom_forcing_type::deallocate_mech_forcing | ( | type(mech_forcing), intent(inout) | forces | ) |
Deallocate the mechanical forcing type.
[in,out] | forces | Forcing fields structure |
Definition at line 2983 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(unit_scale_type), intent(in) | US, | ||
type(forcing), intent(inout) | fluxes, | ||
type(optics_type), pointer | optics, | ||
integer, intent(in) | nsw, | ||
integer, intent(in) | j, | ||
real, intent(in) | dt_in_T, | ||
real, intent(in) | FluxRescaleDepth, | ||
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(max(1,nsw),g%isd:g%ied), intent(out) | pen_SW_bnd, | ||
type(thermo_var_ptrs), intent(inout) | tv, | ||
logical, intent(in) | aggregate_FW, | ||
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(max(1,nsw),g%isd:g%ied), 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] | us | A dimensional unit scaling type |
[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_in_t | The time step for these fluxes [T ~> s] |
[in] | fluxrescaledepth | min ocean depth before fluxes are scaled away [H ~> m or kg m-2] |
[in] | useriverheatcontent | logical for river heat content |
[in] | usecalvingheatcontent | logical for calving heat content |
[in] | h | layer thickness [H ~> m or kg m-2] |
[in] | t | layer temperatures [degC] |
[out] | netmassinout | net mass flux (non-Bouss) or volume flux (if Bouss) of water in/out of ocean over a time step [H ~> m or kg m-2] |
[out] | netmassout | net mass flux (non-Bouss) or volume flux (if Bouss) of water leaving ocean surface over a time step [H ~> m or kg m-2]. 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). [degC H ~> degC m or degC kg m-2]. |
[out] | net_salt | surface salt flux into the ocean accumulated over a time step [ppt H ~> ppt m or ppt kg m-2]. |
[out] | pen_sw_bnd | penetrating SW flux, split into bands. [degC H ~> degC m or degC kg m-2] 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 | For determining how to aggregate forcing. |
[out] | nonpensw | Non-penetrating SW used in net_heat |
[out] | net_heat_rate | Rate of net surface heating |
[out] | net_salt_rate | Surface salt flux into the ocean |
[out] | netmassinout_rate | Rate of net mass flux into the ocean |
[out] | pen_sw_bnd_rate | Rate of penetrative shortwave heating |
[in] | skip_diags | If present and true, skip calculating diagnostics |
Definition at line 346 of file MOM_forcing_type.F90.
References mom_opacity::extract_optics_slice(), mom_error_handler::mom_error(), and mom_opacity::optics_nbands().
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(unit_scale_type), intent(in) | US, | ||
type(forcing), intent(inout) | fluxes, | ||
type(optics_type), pointer | optics, | ||
integer, intent(in) | nsw, | ||
real, intent(in) | dt_in_T, | ||
real, intent(in) | FluxRescaleDepth, | ||
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(max(1,nsw),g%isd:g%ied,g%jsd:g%jed), intent(out) | Pen_SW_bnd, | ||
type(thermo_var_ptrs), intent(inout) | tv, | ||
logical, intent(in) | aggregate_FW | ||
) |
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] | us | A dimensional unit scaling type |
[in,out] | fluxes | structure containing pointers to forcing. |
optics | pointer to optics | |
[in] | nsw | number of bands of penetrating SW |
[in] | dt_in_t | The time step for these fluxes [T ~> s] |
[in] | fluxrescaledepth | min ocean depth before fluxes are scaled away [H ~> m or kg m-2] |
[in] | useriverheatcontent | logical for river heat content |
[in] | usecalvingheatcontent | logical for calving heat content |
[in] | h | layer thickness [H ~> m or kg m-2] |
[in] | t | layer temperatures [degC] |
[out] | netmassinout | net mass flux (non-Bouss) or volume flux (if Bouss) of water in/out of ocean over a time step [H ~> m or kg m-2] |
[out] | netmassout | net mass flux (non-Bouss) or volume flux (if Bouss) of water leaving ocean surface over a time step [H ~> m or kg m-2]. |
[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). [degC H ~> degC m or degC kg m-2] |
[out] | net_salt | surface salt flux into the ocean accumulated over a time step [ppt H ~> ppt m or ppt kg m-2] |
[out] | pen_sw_bnd | penetrating SW flux, by frequency band [degC H ~> degC m or degC kg m-2] with 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 | For determining how to aggregate the forcing. |
Definition at line 817 of file MOM_forcing_type.F90.
References extractfluxes1d().
subroutine, public mom_forcing_type::fluxes_accumulate | ( | type(forcing), intent(in) | flux_tmp, |
type(forcing), intent(inout) | fluxes, | ||
type(ocean_grid_type), intent(inout) | G, | ||
real, intent(out) | wt2, | ||
type(mech_forcing), intent(in), optional | forces | ||
) |
Accumulate the thermodynamic fluxes over time steps.
[in] | flux_tmp | A temporary structure with current thermodynamic forcing fields |
[in,out] | fluxes | A structure containing time-averaged thermodynamic forcing fields |
[in,out] | g | The ocean's grid structure |
[out] | wt2 | The relative weight of the new fluxes |
[in] | forces | A structure with the driving mechanical forces |
Definition at line 1901 of file MOM_forcing_type.F90.
References mom_error_handler::mom_error().
Referenced by forcing_accumulate().
subroutine, public mom_forcing_type::forcing_accumulate | ( | type(forcing), intent(in) | flux_tmp, |
type(mech_forcing), intent(in) | forces, | ||
type(forcing), intent(inout) | fluxes, | ||
type(ocean_grid_type), intent(inout) | G, | ||
real, intent(out) | wt2 | ||
) |
Accumulate the forcing over time steps, taking input from a mechanical forcing type and a temporary forcing-flux type.
[in] | flux_tmp | A temporary structure with current thermodynamic forcing fields |
[in] | forces | A structure with the driving mechanical forces |
[in,out] | fluxes | A structure containing time-averaged thermodynamic forcing fields |
[in,out] | g | The ocean's grid structure |
[out] | wt2 | The relative weight of the new fluxes |
Definition at line 1882 of file MOM_forcing_type.F90.
References fluxes_accumulate().
subroutine, public mom_forcing_type::forcing_diagnostics | ( | type(forcing), intent(in) | fluxes, |
type(surface), intent(in) | sfc_state, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(unit_scale_type), intent(in) | US, | ||
type(time_type), intent(in) | time_end, | ||
type(diag_ctrl), intent(inout) | 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 | A structure containing thermodynamic forcing fields |
[in] | sfc_state | A structure containing fields that describe the surface state of the ocean. |
[in] | g | grid type |
[in] | us | A dimensional unit scaling type |
[in] | time_end | The end time of the diagnostic interval. |
[in,out] | diag | diagnostic regulator |
[in,out] | handles | diagnostic ids |
Definition at line 2238 of file MOM_forcing_type.F90.
References mom_diag_mediator::disable_averaging(), mom_diag_mediator::enable_averages(), mom_spatial_means::global_area_integral(), and mom_spatial_means::global_area_mean().
Referenced by mom_main(), mom_ocean_model_mct::update_ocean_model(), and mom_ocean_model_nuopc::update_ocean_model().
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 | A structure containing thermodynamic forcing fields |
[in] | g | Grid type |
[in] | mesg | Message |
[in] | i | i-index |
[in] | j | j-index |
Definition at line 1161 of file MOM_forcing_type.F90.
References locmsg().
Referenced by mom_diabatic_aux::applyboundaryfluxesinout().
subroutine, public mom_forcing_type::get_net_mass_forcing | ( | type(forcing), intent(in) | fluxes, |
type(ocean_grid_type), intent(in) | G, | ||
type(unit_scale_type), intent(in) | US, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out) | net_mass_src | ||
) |
This subroutine calculates determines the net mass source to the ocean from a (thermodynamic) forcing type and stores it in a provided array.
[in] | fluxes | A structure containing thermodynamic forcing fields |
[in] | g | The ocean grid type |
[in] | us | A dimensional unit scaling type |
[out] | net_mass_src | The net mass flux of water into the ocean [kg m-2 s-1]. |
Definition at line 2142 of file MOM_forcing_type.F90.
Referenced by set_net_mass_forcing().
subroutine, public mom_forcing_type::mech_forcing_diags | ( | type(mech_forcing), intent(in) | forces, |
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(time_type), intent(in) | time_end, | ||
type(diag_ctrl), intent(inout) | 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] | forces | A structure with the driving mechanical forces |
[in] | dt | time step for the forcing [s] |
[in] | g | grid type |
[in] | time_end | The end time of the diagnostic interval. |
[in,out] | diag | diagnostic type |
[in,out] | handles | diagnostic id for diag_manager |
Definition at line 2201 of file MOM_forcing_type.F90.
References mom_diag_mediator::disable_averaging(), and mom_diag_mediator::enable_averaging().
Referenced by mom_main(), mom_ocean_model_mct::update_ocean_model(), and mom_ocean_model_nuopc::update_ocean_model().
|
private |
Write out values of the mechanical forcing arrays at the i,j location. This is a debugging tool.
[in] | forces | A structure with the driving mechanical forces |
[in] | g | Grid type |
[in] | mesg | Message |
[in] | i | i-index |
[in] | j | j-index |
Definition at line 1133 of file MOM_forcing_type.F90.
References locmsg().
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, | ||
type(unit_scale_type), intent(in) | US, | ||
integer, intent(in), optional | haloshift | ||
) |
Write out chksums for thermodynamic fluxes.
[in] | mesg | message |
[in] | fluxes | A structure containing thermodynamic forcing fields |
[in] | g | grid type |
[in] | us | A dimensional unit scaling type |
[in] | haloshift | shift in halo |
Definition at line 1012 of file MOM_forcing_type.F90.
Referenced by mom_main(), and mom::step_mom_thermo().
subroutine, public mom_forcing_type::mom_mech_forcing_chksum | ( | character(len=*), intent(in) | mesg, |
type(mech_forcing), intent(in) | forces, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(unit_scale_type), intent(in) | US, | ||
integer, intent(in), optional | haloshift | ||
) |
Write out chksums for the driving mechanical forces.
[in] | mesg | message |
[in] | forces | A structure with the driving mechanical forces |
[in] | g | grid type |
[in] | us | A dimensional unit scaling type |
[in] | haloshift | shift in halo |
Definition at line 1104 of file MOM_forcing_type.F90.
Referenced by mom_main(), and mom::step_mom().
|
private |
Allocates and zeroes-out array.
array | Array to be allocated | |
[in] | is | Start i-index |
[in] | ie | End i-index |
[in] | js | Start j-index |
[in] | je | End j-index |
[in] | flag | Flag to indicate to allocate |
Definition at line 2917 of file MOM_forcing_type.F90.
Referenced by allocate_forcing_type(), and allocate_mech_forcing().
subroutine, public mom_forcing_type::register_forcing_type_diags | ( | type(time_type), intent(in) | Time, |
type(diag_ctrl), intent(inout) | diag, | ||
type(unit_scale_type), intent(in) | US, | ||
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] | us | A dimensional unit scaling 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 1221 of file MOM_forcing_type.F90.
subroutine, public mom_forcing_type::set_derived_forcing_fields | ( | type(mech_forcing), intent(in) | forces, |
type(forcing), intent(inout) | fluxes, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(unit_scale_type), intent(in) | US, | ||
real, intent(in) | Rho0 | ||
) |
This subroutine calculates certain derived forcing fields based on information from a mech_forcing type and stores them in a (thermodynamic) forcing type.
[in] | forces | A structure with the driving mechanical forces |
[in,out] | fluxes | A structure containing thermodynamic forcing fields |
[in] | g | grid type |
[in] | us | A dimensional unit scaling type |
[in] | rho0 | A reference density of seawater [R ~> kg m-3], as used to calculate ustar. |
Definition at line 2089 of file MOM_forcing_type.F90.
subroutine, public mom_forcing_type::set_net_mass_forcing | ( | type(forcing), intent(in) | fluxes, |
type(mech_forcing), intent(inout) | forces, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(unit_scale_type), intent(in) | US | ||
) |
This subroutine determines the net mass source to the ocean from a (thermodynamic) forcing type and stores it in a mech_forcing type.
[in] | fluxes | A structure containing thermodynamic forcing fields |
[in,out] | forces | A structure with the driving mechanical forces |
[in] | us | A dimensional unit scaling type |
[in] | g | The ocean grid type |
Definition at line 2129 of file MOM_forcing_type.F90.
References get_net_mass_forcing().