namespace mom_forcing_type

Overview

This module implements boundary forcing for MOM6. More…

namespace mom_forcing_type {

// global functions

subroutine, public extractfluxes1d(
    G G,
    GV GV,
    fluxes fluxes,
    optics optics,
    nsw nsw,
    j j,
    dt dt,
    FluxRescaleDepth FluxRescaleDepth,
    useRiverHeatContent useRiverHeatContent,
    useCalvingHeatContent useCalvingHeatContent,
    h h,
    T T,
    netMassInOut netMassInOut,
    netMassOut netMassOut,
    net_heat net_heat,
    net_salt net_salt,
    pen_SW_bnd pen_SW_bnd,
    tv tv,
    aggregate_FW aggregate_FW,
    nonpenSW nonpenSW,
    netmassInOut_rate netmassInOut_rate,
    net_Heat_Rate net_Heat_Rate,
    net_salt_rate net_salt_rate,
    pen_sw_bnd_Rate pen_sw_bnd_Rate,
    skip_diags skip_diags
    );

subroutine, public extractfluxes2d(
    G G,
    GV GV,
    fluxes fluxes,
    optics optics,
    nsw nsw,
    dt dt,
    FluxRescaleDepth FluxRescaleDepth,
    useRiverHeatContent useRiverHeatContent,
    useCalvingHeatContent useCalvingHeatContent,
    h h,
    T T,
    netMassInOut netMassInOut,
    netMassOut netMassOut,
    net_heat net_heat,
    Net_salt Net_salt,
    Pen_SW_bnd Pen_SW_bnd,
    tv tv,
    aggregate_FW aggregate_FW
    );

subroutine, public calculatebuoyancyflux1d(
    G G,
    GV GV,
    US US,
    fluxes fluxes,
    optics optics,
    nsw nsw,
    h h,
    Temp Temp,
    Salt Salt,
    tv tv,
    j j,
    buoyancyFlux buoyancyFlux,
    netHeatMinusSW netHeatMinusSW,
    netSalt netSalt,
    skip_diags skip_diags
    );

subroutine, public calculatebuoyancyflux2d(
    G G,
    GV GV,
    US US,
    fluxes fluxes,
    optics optics,
    h h,
    Temp Temp,
    Salt Salt,
    tv tv,
    buoyancyFlux buoyancyFlux,
    netHeatMinusSW netHeatMinusSW,
    netSalt netSalt,
    skip_diags skip_diags
    );

subroutine, public mom_forcing_chksum(mesg mesg, fluxes fluxes, G G, US US, haloshift haloshift);
subroutine, public mom_mech_forcing_chksum(mesg mesg, forces forces, G G, US US, haloshift haloshift);
subroutine, public forcing_singlepointprint(fluxes fluxes, G G, i i, j j, mesg mesg);

subroutine, public register_forcing_type_diags(
    Time Time,
    diag diag,
    US US,
    use_temperature use_temperature,
    handles handles,
    use_berg_fluxes use_berg_fluxes
    );

subroutine, public forcing_accumulate(
    flux_tmp flux_tmp,
    forces forces,
    fluxes fluxes,
    dt dt,
    G G,
    wt2 wt2
    );

subroutine, public fluxes_accumulate(
    flux_tmp flux_tmp,
    fluxes fluxes,
    dt dt,
    G G,
    wt2 wt2,
    forces forces
    );

subroutine, public copy_common_forcing_fields(forces forces, fluxes fluxes, G G, skip_pres skip_pres);
subroutine, public set_derived_forcing_fields(forces forces, fluxes fluxes, G G, US US, Rho0 Rho0);
subroutine, public set_net_mass_forcing(fluxes fluxes, forces forces, G G);
subroutine, public get_net_mass_forcing(fluxes fluxes, G G, net_mass_src net_mass_src);
subroutine, public copy_back_forcing_fields(fluxes fluxes, forces forces, G G);
subroutine, public mech_forcing_diags(forces forces, dt dt, G G, diag diag, handles handles);

subroutine, public forcing_diagnostics(
    fluxes fluxes,
    sfc_state sfc_state,
    dt dt,
    G G,
    diag diag,
    handles handles
    );

subroutine, public allocate_forcing_type(
    G G,
    fluxes fluxes,
    water water,
    heat heat,
    ustar ustar,
    press press,
    shelf shelf,
    iceberg iceberg,
    salt salt
    );

subroutine, public allocate_mech_forcing(
    G G,
    forces forces,
    stress stress,
    ustar ustar,
    shelf shelf,
    press press,
    iceberg iceberg
    );

subroutine, public deallocate_forcing_type(fluxes fluxes);
subroutine, public deallocate_mech_forcing(forces forces);

} // namespace mom_forcing_type

Detailed Documentation

This module implements boundary forcing for MOM6.

Boundary fluxes

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.

Surface boundary momentum fluxes

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.

Surface boundary mass fluxes

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)).

  • Evaporation field can in fact represent a mass loss (evaporation) or mass gain (condensation in foggy areas).

  • sea ice formation leads to mass moving from the liquid ocean to the ice model, and melt adds liquid to the ocean.

  • Precipitation can be liquid or frozen (snow). Furthermore, in some versions of the GFDL coupler, precipitation can be negative. The reason is that the ice model combines precipitation with ice melt and ice formation. This limitation of the ice model diagnostics should be overcome future versions.

  • River runoff can be liquid or frozen. Frozen runoff is often associated with calving land-ice and/or ice bergs.

Surface boundary salt fluxes

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.

Surface boundary heat fluxes

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}\]
  • shortwave (SW) = shortwave radiation (always warms ocean)

  • longwave (LW) = longwave radiation (generally cools ocean)

  • latent (LAT) = turbulent latent heat loss due to evaporation (liquid to vapor) or melt (snow to liquid); generally cools the ocean

  • sensible (SENS) = turbulent heat transfer due to differences in air-sea or ice-sea temperature

  • mass transfer (MASS) = heat transfer due to heat content of mass (e.g., E-P+R) transferred across ocean surface; computed relative to 0 Celsius

  • frazil (FRAZ) = heat transferred to form frazil sea ice (positive heating of liquid ocean)

  • restore (RES) = heat from surface damping sometimes imposed in non-coupled model simulations .

  • restore (flux adjustments) = heat from surface flux adjustment.

Treatment of shortwave

The shortwave field itself is split into two pieces:

  • shortwave = penetrative SW + non-penetrative SW

  • non-penetrative = non-downwelling shortwave; portion of SW totally absorbed in the k=1 cell. The non-penetrative SW is combined with LW+LAT+SENS+seaice_melt_heat in net_heat inside routine extractFluxes1d. Notably, for many cases, non-penetrative SW = 0.

  • penetrative = that portion of shortwave penetrating below a tiny surface layer. This is the downwelling shortwave. Penetrative SW participates in the penetrative SW heating of k=1,nz cells, with the amount of penetration dependent on optical properties.

Convergence of heat into the k=1 cell

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{split}\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*}\end{split}\]

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

Global Functions

subroutine, public extractfluxes1d(
    G G,
    GV GV,
    fluxes fluxes,
    optics optics,
    nsw nsw,
    j j,
    dt dt,
    FluxRescaleDepth FluxRescaleDepth,
    useRiverHeatContent useRiverHeatContent,
    useCalvingHeatContent useCalvingHeatContent,
    h h,
    T T,
    netMassInOut netMassInOut,
    netMassOut netMassOut,
    net_heat net_heat,
    net_salt net_salt,
    pen_SW_bnd pen_SW_bnd,
    tv tv,
    aggregate_FW aggregate_FW,
    nonpenSW nonpenSW,
    netmassInOut_rate netmassInOut_rate,
    net_Heat_Rate net_Heat_Rate,
    net_salt_rate net_salt_rate,
    pen_sw_bnd_Rate pen_sw_bnd_Rate,
    skip_diags 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.

Parameters:

g

ocean grid structure

gv

ocean vertical grid structure

fluxes

structure containing pointers to possible forcing fields. NULL unused fields.

optics

pointer to optics

nsw

number of bands of penetrating SW

j

j-index to work on

dt

time step [s]

fluxrescaledepth

min ocean depth before fluxes are scaled away [H ~> m or kg m-2]

useriverheatcontent

logical for river heat content

usecalvingheatcontent

logical for calving heat content

h

layer thickness [H ~> m or kg m-2]

t

layer temperatures [degC]

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]

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.

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].

net_salt

surface salt flux into the ocean accumulated over a time step [ppt H ~> ppt m or ppt kg m-2].

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.

tv

structure containing pointers to available thermodynamic fields. Used to keep track of the heat flux associated with net mass fluxes into the ocean.

aggregate_fw

For determining how to aggregate forcing.

nonpensw

Non-penetrating SW used in net_heat

net_heat_rate

Rate of net surface heating

net_salt_rate

Surface salt flux into the ocean

netmassinout_rate

Rate of net mass flux into the ocean

pen_sw_bnd_rate

Rate of penetrative shortwave heating

skip_diags

If present and true, skip calculating diagnostics

subroutine, public extractfluxes2d(
    G G,
    GV GV,
    fluxes fluxes,
    optics optics,
    nsw nsw,
    dt dt,
    FluxRescaleDepth FluxRescaleDepth,
    useRiverHeatContent useRiverHeatContent,
    useCalvingHeatContent useCalvingHeatContent,
    h h,
    T T,
    netMassInOut netMassInOut,
    netMassOut netMassOut,
    net_heat net_heat,
    Net_salt Net_salt,
    Pen_SW_bnd Pen_SW_bnd,
    tv tv,
    aggregate_FW 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.

Parameters:

g

ocean grid structure

gv

ocean vertical grid structure

fluxes

structure containing pointers to forcing.

optics

pointer to optics

nsw

number of bands of penetrating SW

dt

time step [s]

fluxrescaledepth

min ocean depth before fluxes are scaled away [H ~> m or kg m-2]

useriverheatcontent

logical for river heat content

usecalvingheatcontent

logical for calving heat content

h

layer thickness [H ~> m or kg m-2]

t

layer temperatures [degC]

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]

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].

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]

net_salt

surface salt flux into the ocean accumulated over a time step [ppt H ~> ppt m or ppt kg m-2]

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.

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.

aggregate_fw

For determining how to aggregate the forcing.

subroutine, public calculatebuoyancyflux1d(
    G G,
    GV GV,
    US US,
    fluxes fluxes,
    optics optics,
    nsw nsw,
    h h,
    Temp Temp,
    Salt Salt,
    tv tv,
    j j,
    buoyancyFlux buoyancyFlux,
    netHeatMinusSW netHeatMinusSW,
    netSalt netSalt,
    skip_diags 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.

Parameters:

g

ocean grid

gv

ocean vertical grid structure

us

A dimensional unit scaling type

fluxes

surface fluxes

optics

penetrating SW optics

nsw

The number of frequency bands of penetrating shortwave radiation

h

layer thickness [H ~> m or kg m-2]

temp

prognostic temp [degC]

salt

salinity [ppt]

tv

thermodynamics type

j

j-row to work on

buoyancyflux

buoyancy fluxes [L2 T-3 ~> m2 s-3]

netheatminussw

surf Heat flux [degC H s-1 ~> degC m s-1 or degC kg m-2 s-1]

netsalt

surf salt flux [ppt H s-1 ~> ppt m s-1 or ppt kg m-2 s-1]

skip_diags

If present and true, skip calculating diagnostics inside extractFluxes1d()

subroutine, public calculatebuoyancyflux2d(
    G G,
    GV GV,
    US US,
    fluxes fluxes,
    optics optics,
    h h,
    Temp Temp,
    Salt Salt,
    tv tv,
    buoyancyFlux buoyancyFlux,
    netHeatMinusSW netHeatMinusSW,
    netSalt netSalt,
    skip_diags skip_diags
    )

Calculates surface buoyancy flux by adding up the heat, FW and salt fluxes, for 2d arrays. This is a wrapper for calculateBuoyancyFlux1d.

Parameters:

g

ocean grid

gv

ocean vertical grid structure

us

A dimensional unit scaling type

fluxes

surface fluxes

optics

SW ocean optics

h

layer thickness [H ~> m or kg m-2]

temp

temperature [degC]

salt

salinity [ppt]

tv

thermodynamics type

buoyancyflux

buoyancy fluxes [L2 T-3 ~> m2 s-3]

netheatminussw

surf temp flux [degC H ~> degC m or degC kg m-2]

netsalt

surf salt flux [ppt H ~> ppt m or ppt kg m-2]

skip_diags

If present and true, skip calculating diagnostics inside extractFluxes1d()

subroutine, public mom_forcing_chksum(
    mesg mesg,
    fluxes fluxes,
    G G,
    US US,
    haloshift haloshift
    )

Write out chksums for thermodynamic fluxes.

Parameters:

mesg

message

fluxes

A structure containing thermodynamic forcing fields

g

grid type

us

A dimensional unit scaling type

haloshift

shift in halo

subroutine, public mom_mech_forcing_chksum(
    mesg mesg,
    forces forces,
    G G,
    US US,
    haloshift haloshift
    )

Write out chksums for the driving mechanical forces.

Parameters:

mesg

message

forces

A structure with the driving mechanical forces

g

grid type

us

A dimensional unit scaling type

haloshift

shift in halo

subroutine, public forcing_singlepointprint(
    fluxes fluxes,
    G G,
    i i,
    j j,
    mesg mesg
    )

Write out values of the fluxes arrays at the i,j location. This is a debugging tool.

Parameters:

fluxes

A structure containing thermodynamic forcing fields

g

Grid type

mesg

Message

i

i-index

j

j-index

subroutine, public register_forcing_type_diags(
    Time Time,
    diag diag,
    US US,
    use_temperature use_temperature,
    handles handles,
    use_berg_fluxes use_berg_fluxes
    )

Register members of the forcing type for diagnostics.

Parameters:

time

time type

diag

diagnostic control type

us

A dimensional unit scaling type

use_temperature

True if T/S are in use

handles

handles for diagnostics

use_berg_fluxes

If true, allow iceberg flux diagnostics

subroutine, public forcing_accumulate(
    flux_tmp flux_tmp,
    forces forces,
    fluxes fluxes,
    dt dt,
    G G,
    wt2 wt2
    )

Accumulate the forcing over time steps, taking input from a mechanical forcing type and a temporary forcing-flux type.

Parameters:

flux_tmp

A temporary structure with current thermodynamic forcing fields

forces

A structure with the driving mechanical forces

fluxes

A structure containing time-averaged thermodynamic forcing fields

dt

The elapsed time since the last call to this subroutine [s]

g

The ocean’s grid structure

wt2

The relative weight of the new fluxes

subroutine, public fluxes_accumulate(
    flux_tmp flux_tmp,
    fluxes fluxes,
    dt dt,
    G G,
    wt2 wt2,
    forces forces
    )

Accumulate the thermodynamic fluxes over time steps.

Parameters:

flux_tmp

A temporary structure with current thermodynamic forcing fields

fluxes

A structure containing time-averaged thermodynamic forcing fields

dt

The elapsed time since the last call to this subroutine [s]

g

The ocean’s grid structure

wt2

The relative weight of the new fluxes

forces

A structure with the driving mechanical forces

subroutine, public copy_common_forcing_fields(
    forces forces,
    fluxes fluxes,
    G G,
    skip_pres skip_pres
    )

This subroutine copies the computational domains of common forcing fields from a mech_forcing type to a (thermodynamic) forcing type.

Parameters:

forces

A structure with the driving mechanical forces

fluxes

A structure containing thermodynamic forcing fields

g

grid type

skip_pres

If present and true, do not copy pressure fields.

subroutine, public set_derived_forcing_fields(
    forces forces,
    fluxes fluxes,
    G G,
    US US,
    Rho0 Rho0
    )

This subroutine calculates certain derived forcing fields based on information from a mech_forcing type and stores them in a (thermodynamic) forcing type.

Parameters:

forces

A structure with the driving mechanical forces

fluxes

A structure containing thermodynamic forcing fields

g

grid type

us

A dimensional unit scaling type

rho0

A reference density of seawater [kg m-3], as used to calculate ustar.

subroutine, public set_net_mass_forcing(fluxes fluxes, forces forces, G G)

This subroutine determines the net mass source to the ocean from a (thermodynamic) forcing type and stores it in a mech_forcing type.

Parameters:

fluxes

A structure containing thermodynamic forcing fields

forces

A structure with the driving mechanical forces

g

The ocean grid type

subroutine, public get_net_mass_forcing(
    fluxes fluxes,
    G G,
    net_mass_src 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.

Parameters:

fluxes

A structure containing thermodynamic forcing fields

g

The ocean grid type

net_mass_src

The net mass flux of water into the ocean [kg m-2 s-1].

subroutine, public copy_back_forcing_fields(fluxes fluxes, forces forces, G G)

This subroutine copies the computational domains of common forcing fields from a mech_forcing type to a (thermodynamic) forcing type.

Parameters:

fluxes

A structure containing thermodynamic forcing fields

forces

A structure with the driving mechanical forces

g

grid type

subroutine, public mech_forcing_diags(
    forces forces,
    dt dt,
    G G,
    diag diag,
    handles handles
    )

Offer mechanical forcing fields for diagnostics for those fields registered as part of register_forcing_type_diags.

Parameters:

forces

A structure with the driving mechanical forces

dt

time step

g

grid type

diag

diagnostic type

handles

diagnostic id for diag_manager

subroutine, public forcing_diagnostics(
    fluxes fluxes,
    sfc_state sfc_state,
    dt dt,
    G G,
    diag diag,
    handles handles
    )

Offer buoyancy forcing fields for diagnostics for those fields registered as part of register_forcing_type_diags.

Parameters:

fluxes

A structure containing thermodynamic forcing fields

sfc_state

A structure containing fields that describe the surface state of the ocean.

dt

time step

g

grid type

diag

diagnostic regulator

handles

diagnostic ids

subroutine, public allocate_forcing_type(
    G G,
    fluxes fluxes,
    water water,
    heat heat,
    ustar ustar,
    press press,
    shelf shelf,
    iceberg iceberg,
    salt salt
    )

Conditionally allocate fields within the forcing type.

Parameters:

g

Ocean grid structure

fluxes

A structure containing thermodynamic forcing fields

water

If present and true, allocate water fluxes

heat

If present and true, allocate heat fluxes

ustar

If present and true, allocate ustar and related fields

press

If present and true, allocate p_surf and related fields

shelf

If present and true, allocate fluxes for ice-shelf

iceberg

If present and true, allocate fluxes for icebergs

salt

If present and true, allocate salt fluxes

subroutine, public allocate_mech_forcing(
    G G,
    forces forces,
    stress stress,
    ustar ustar,
    shelf shelf,
    press press,
    iceberg iceberg
    )

Conditionally allocate fields within the mechanical forcing type.

Parameters:

g

Ocean grid structure

forces

Forcing fields structure

stress

If present and true, allocate taux, tauy

ustar

If present and true, allocate ustar and related fields

shelf

If present and true, allocate forces for ice-shelf

press

If present and true, allocate p_surf and related fields

iceberg

If present and true, allocate forces for icebergs

subroutine, public deallocate_forcing_type(fluxes fluxes)

Deallocate the forcing type.

Parameters:

fluxes

Forcing fields structure

subroutine, public deallocate_mech_forcing(forces forces)

Deallocate the mechanical forcing type.

Parameters:

forces

Forcing fields structure