mom module reference

The central module of the MOM6 ocean model.

More…

Data Types

mom_diag_ids

A structure with diagnostic IDs of the state variables.

mom_control_struct

Control structure for the MOM module, including the variables that describe the state of the ocean.

Functions/Subroutines

step_mom()

This subroutine orchestrates the time stepping of MOM.

step_mom_dynamics()

Time step the ocean dynamics, including the momentum and continuity equations.

step_mom_tracer_dyn()

step_MOM_tracer_dyn does tracer advection and lateral diffusion, bringing the tracers up to date with the changes in state due to the dynamics. Surface sources and sinks and remapping are handled via step_MOM_thermo.

step_mom_thermo()

MOM_step_thermo orchestrates the thermodynamic time stepping and vertical remapping, via calls to diabatic (or adiabatic) and ALE_regrid.

step_offline()

step_offline is the main driver for running tracers offline in MOM6. This has been primarily developed with ALE configurations in mind. Some work has been done in isopycnal configuration, but the work is very preliminary. Some more detail about this capability along with some of the subroutines called here can be found in tracers/MOM_offline_control.F90

initialize_mom()

Initialize MOM, including memory allocation, setting up parameters and diagnostics, initializing the ocean state variables, and initializing subsidiary modules.

finish_mom_initialization()

Finishes initializing MOM and writes out the initial conditions.

register_diags()

Register certain diagnostics.

mom_timing_init()

Set up CPU clock IDs for timing various subroutines.

set_restart_fields()

Set the fields that are needed for bitwise identical restarting the time stepping scheme.

adjust_ssh_for_p_atm()

Apply a correction to the sea surface height to compensate for the atmospheric pressure (the inverse barometer).

extract_surface_state()

Set the surface (return) properties of the ocean model by setting the appropriate fields in sfc_state.

rotate_initial_state()

Rotate initialization fields from input to rotated arrays.

mom_state_is_synchronized()

Return true if all phases of step_MOM are at the same point in time.

get_mom_state_elements()

This subroutine offers access to values or pointers to other types from within the MOM_control_struct, allowing the MOM_control_struct to be opaque.

get_ocean_stocks()

Find the global integrals of various quantities.

save_mom_restart()

Save restart/pickup files required to initialize the MOM6 internal state.

mom_end()

End of ocean model, including memory deallocation.

Detailed Description

Modular Ocean Model (MOM) Version 6.0 (MOM6)

Additional contributions from: * Whit Anderson

  • Brian Arbic

  • Will Cooke

  • Anand Gnanadesikan

  • Matthew Harrison

  • Mehmet Ilicak

  • Laura Jackson

  • Jasmine John

  • John Krasting

  • Zhi Liang

  • Bonnie Samuels

  • Harper Simmons

  • Laurent White

  • Niki Zadeh

MOM ice-shelf code was developed by * Daniel Goldberg

  • Robert Hallberg

  • Chris Little

  • Olga Sergienko

Overview of MOM

This program (MOM) simulates the ocean by numerically solving the hydrostatic primitive equations in generalized Lagrangian vertical coordinates, typically tracking stretched pressure (p*) surfaces or following isopycnals in the ocean’s interior, and general orthogonal horizontal coordinates. Unlike earlier versions of MOM, in MOM6 these equations are horizontally discretized on an Arakawa C-grid. (It remains to be seen whether a B-grid dynamic core will be revived in MOM6 at a later date; for now applications requiring a B-grid discretization should use MOM5.1.) MOM6 offers a range of options for the physical parameterizations, from those most appropriate to highly idealized models for geophysical fluid dynamics studies to a rich suite of processes appropriate for realistic ocean simulations. The thermodynamic options typically use conservative temperature and preformed salinity as conservative state variables and a full nonlinear equation of state, but there are also idealized adiabatic configurations of the model that use fixed density layers. Version 6.0 of MOM continues in the long tradition of a commitment to climate-quality ocean simulations embodied in previous versions of MOM, even as it draws extensively on the lessons learned in the development of the Generalized Ocean Layered Dynamics (GOLD) ocean model, which was also primarily developed at NOAA/GFDL. MOM has also benefited tremendously from the FMS infrastructure, which it utilizes and shares with other component models developed at NOAA/GFDL.

When run is isopycnal-coordinate mode, the uppermost few layers are often used to describe a bulk mixed layer, including the effects of penetrating shortwave radiation. Either a split- explicit time stepping scheme or a non-split scheme may be used for the dynamics, while the time stepping may be split (and use different numbers of steps to cover the same interval) for the forcing, the thermodynamics, and for the dynamics. Most of the numerics are second order accurate in space. MOM can run with an absurdly thin minimum layer thickness. A variety of non-isopycnal vertical coordinate options are under development, but all exploit the advantages of a Lagrangian vertical coordinate, as discussed in detail by Adcroft and Hallberg (Ocean Modelling, 2006).

Details of the numerics and physical parameterizations are provided in the appropriate source files. All of the available options are selected at run-time by parsing the input files, usually MOM_input and MOM_override, and the options choices are then documented for each run in MOM_param_docs.

MOM6 integrates the equations forward in time in three distinct phases. In one phase, the dynamic equations for the velocities and layer thicknesses are advanced, capturing the propagation of external and internal inertia-gravity waves, Rossby waves, and other strictly adiabatic processes, including lateral stresses, vertical viscosity and momentum forcing, and interface height diffusion (commonly called Gent-McWilliams diffusion in depth- coordinate models). In the second phase, all tracers are advected and diffused along the layers. The third phase applies diabatic processes, vertical mixing of water properties, and perhaps vertical remapping to cause the layers to track the desired vertical coordinate.

The present file ( MOM.F90) orchestrates the main time stepping loops. One time integration option for the dynamics uses a split explicit time stepping scheme to rapidly step the barotropic pressure and velocity fields. The barotropic velocities are averaged over the baroclinic time step before they are used to advect thickness and determine the baroclinic accelerations. As described in Hallberg and Adcroft (2009), a barotropic correction is applied to the time-mean layer velocities to ensure that the sum of the layer transports agrees with the time-mean barotropic transport, thereby ensuring that the estimates of the free surface from the sum of the layer thicknesses agrees with the final free surface height as calculated by the barotropic solver. The barotropic and baroclinic velocities are kept consistent by recalculating the barotropic velocities from the baroclinic transports each time step. This scheme is described in Hallberg, 1997, J. Comp. Phys. 135, 54-65 and in Hallberg and Adcroft, 2009, Ocean Modelling, 29, 15-26.

The other time integration options use non-split time stepping schemes based on the 3-step third order Runge-Kutta scheme described in Matsuno, 1966, J. Met. Soc. Japan, 44, 85-88, or on a two-step quasi-2nd order Runge-Kutta scheme. These are much slower than the split time-stepping scheme, but they are useful for providing a more robust solution for debugging cases where the more complicated split time-stepping scheme may be giving suspect solutions.

There are a range of closure options available. Horizontal velocities are subject to a combination of horizontal biharmonic and Laplacian friction (based on a stress tensor formalism) and a vertical Fickian viscosity (perhaps using the kinematic viscosity of water). The horizontal viscosities may be constant, spatially varying or may be dynamically calculated using Smagorinsky’s approach. A diapycnal diffusion of density and thermodynamic quantities is also allowed, but not required, as is horizontal diffusion of interface heights (akin to the Gent-McWilliams closure of geopotential coordinate models). The diapycnal mixing may use a fixed diffusivity or it may use the shear Richardson number dependent closure, like that described in Jackson et al. (JPO, 2008). When there is diapycnal diffusion, it applies to momentum as well. As this is in addition to the vertical viscosity, the vertical Prandtl always exceeds 1. A refined bulk-mixed layer is often used to describe the planetary boundary layer in realistic ocean simulations.

MOM has a number of noteworthy debugging capabilities. Excessively large velocities are truncated and MOM will stop itself after a number of such instances to keep the model from crashing altogether. This is useful in diagnosing failures, or (by accepting some truncations) it may be useful for getting the model past the adjustment from an ill-balanced initial condition. In addition, all of the accelerations in the columns with excessively large velocities may be directed to a text file. Parallelization errors may be diagnosed using the DEBUG option, which causes extensive checksums to be written out along with comments indicating where in the algorithm the sums originate and what variable is being summed. The point where these checksums differ between runs is usually a good indication of where in the code the problem lies. All of the test cases provided with MOM are routinely tested to ensure that they give bitwise identical results regardless of the domain decomposition, or whether they use static or dynamic memory allocation.

Structure of MOM

About 115 other files of source code and 4 header files comprise the MOM code, although there are several hundred more files that make up the FMS infrastructure upon which MOM is built. Each of the MOM files contains comments documenting what it does, and most of the file names are fairly self-evident. In addition, all subroutines and data types are referenced via a module use, only statement, and the module names are consistent with the file names, so it is not too hard to find the source file for a subroutine.

The typical MOM directory tree is as follows:

../MOM
|-- ac
|-- config_src
|   |-- drivers
|   !   |-- FMS_cap
|   !   |-- ice_solo_driver
|   !   |-- mct_cap
|   !   |-- nuopc_cap
|   !   |-- solo_driver
|   !   `-- unit_drivers
|   |-- external
|   !   |-- drifters
|   !   |-- GFDL_ocean_BGC
|   !   `-- ODA_hooks
|   |-- infra
|   !   |-- FMS1
|   !   `-- FMS2
|   `-- memory
|   !   |-- dynamic_nonsymmetric
|   !   `-- dynamic_symmetric
|-- docs
|-- pkg
|   |-- CVMix-src
|   |-- ...
|   `-- MOM6_DA_hooks
`-- src
    |-- ALE
    |-- core
    |-- diagnostics
    |-- equation_of_state
    |-- framework
    |-- ice_shelf
    |-- initialization
    |-- ocean_data_assim
    |-- parameterizations
    |   |-- CVMix
    |   |-- lateral
    |   `-- vertical
    |-- tracer
    `-- user

Rather than describing each file here, selected directory contents will be described to give a broad overview of the MOM code structure.

The directories under config_src contain files that are used for configuring the code, for instance for coupled or ocean-only runs. Only one or two of these directories are used in compiling any, particular run.

  • config_src/drivers/FMS-cap: The files here are used to couple MOM as a component in a larger run driven by the FMS coupler. This includes code that converts various forcing fields into the code structures and flux and unit conventions used by MOM, and converts the MOM surface fields back to the forms used by other FMS components.

  • config_src/drivers/nuopc-cap: The files here are used to couple MOM as a component in a larger run driven by the NUOPC coupler. This includes code that converts various forcing fields into the code structures and flux and unit conventions used by MOM, and converts the MOM surface fields back to the forms used by other NUOPC components.

  • config_src/drivers/solo_driver: The files here are include the _main driver that is used when MOM is configured as an ocean-only model, as well as the files that specify the surface forcing in this configuration.

  • config_src/external: The files here are mostly just stubs, so that MOM6 can compile with calls to the public interfaces external packages, but without actually requiring those packages themselves. In more elaborate configurations, would be linked to the actual code for those external packages rather than these simple stubs.

  • config_src/memory/dynamic-symmetric: The only file here is the version of MOM_memory.h that is used for dynamic memory configurations of MOM.

The directories under src contain most of the MOM files. These files are used in every configuration using MOM.

  • src/core: The files here constitute the MOM dynamic core. This directory also includes files with the types that describe the model’s lateral grid and have defined types that are shared across various MOM modules to allow for more succinct and flexible subroutine argument lists.

  • src/diagnostics: The files here calculate various diagnostics that are ancilliary to the model itself. While most of these diagnostics do not directly affect the model’s solution, there are some, like the calculation of the deformation radius, that are used in some of the process parameterizations.

  • src/equation_of_state: These files describe the physical properties of sea-water, including both the equation of state and when it freezes.

  • src/framework: These files provide infrastructure utilities for MOM. Many are simply wrappers for capabilities provided by FMS, although others provide capabilities (like the file_parser) that are unique to MOM. When MOM is adapted to use a modeling infrastructure distinct from FMS, most of the required changes are in this directory.

  • src/initialization: These are the files that are used to initialize the MOM grid or provide the initial physical state for MOM. These files are not intended to be modified, but provide a means for calling user-specific initialization code like the examples in src/user.

  • src/parameterizations/lateral: These files implement a number of quasi-lateral (along-layer) process parameterizations, including lateral viscosities, parameterizations of eddy effects, and the calculation of tidal forcing.

  • src/parameterizations/vertical: These files implement a number of vertical mixing or diabatic processes, including the effects of vertical viscosity and code to parameterize the planetary boundary layer. There is a separate driver that orchestrates this portion of the algorithm, and there is a diversity of parameterizations to be found here.

  • src/tracer: These files handle the lateral transport and diffusion of tracers, or are the code to implement various passive tracer packages. Additional tracer packages are readily accommodated.

  • src/user: These are either stub routines that a user could use to change the model’s initial conditions or forcing, or are examples that implement specific test cases. These files can easily be hand edited to create new analytically specified configurations.

Most simulations can be set up by modifying only the files MOM_input, and possibly one or two of the files in src/user. In addition, the diag_table (MOM_diag_table) will commonly be modified to tailor the output to the needs of the question at hand. The FMS utility mkmf works with a file called path_names to build an appropriate makefile, and path_names should be edited to reflect the actual location of the desired source code.

The separate MOM-examples git repository provides a large number of working configurations of MOM, along with reference solutions for several different compilers on GFDL’s latest large computer. The versions of MOM_memory.h in these directories need not be used if dynamic memory allocation is desired, and the answers should be unchanged.

There are 3 publicly visible subroutines in this file ( MOM.F90). * step_MOM steps MOM over a specified interval of time.

  • MOM_initialize calls initialize and does other initialization that does not warrant user modification.

  • extract_surface_state determines the surface (bulk mixed layer if traditional isopycnal vertical coordinate) properties of the current model state and packages pointers to these fields into an exported structure.

The remaining subroutines in this file ( src/core/MOM.F90) are:

  • find_total_transport determines the barotropic mass transport.

  • register_diags registers many diagnostic fields for the dynamic solver, or of the main model variables.

  • MOM_timing_init initializes various CPU time clocks.

  • write_static_fields writes out various time-invariant fields.

  • set_restart_fields is used to specify those fields that are written to and read from the restart file.

Diagnosing MOM heat budget

Here are some example heat budgets for the ALE version of MOM6.

Depth integrated heat budget

Depth integrated heat budget diagnostic for MOM.

  • OPOTTEMPTEND_2d = T_ADVECTION_XY_2d + OPOTTEMPPMDIFF_2d + HFDS + HFGEOU

  • T_ADVECTION_XY_2d = horizontal advection

  • OPOTTEMPPMDIFF_2d = neutral diffusion

  • HFDS = net surface boundary heat flux

  • HFGEOU = geothermal heat flux

  • HFDS = net surface boundary heat flux entering the ocean = rsntds + rlntds + hfls + hfss + heat_pme + hfsifrazil

  • More heat flux cross-checks * hfds = net_heat_coupler + hfsifrazil + heat_pme

    • heat_pme = heat_content_surfwater = heat_content_massin + heat_content_massout = heat_content_fprec + heat_content_cond + heat_content_vprec * hfrunoffds + hfevapds + hfrainds

Depth integrated heat budget

Here is an example 3d heat budget diagnostic for MOM.

  • OPOTTEMPTEND = T_ADVECTION_XY + TH_TENDENCY_VERT_REMAP + OPOTTEMPDIFF + OPOTTEMPPMDIFF * BOUNDARY_FORCING_HEAT_TENDENCY + FRAZIL_HEAT_TENDENCY

  • OPOTTEMPTEND = net tendency of heat as diagnosed in MOM.F90

  • T_ADVECTION_XY = heating of a cell from lateral advection

  • TH_TENDENCY_VERT_REMAP = heating of a cell from vertical remapping

  • OPOTTEMPDIFF = heating of a cell from diabatic diffusion

  • OPOTTEMPPMDIFF = heating of a cell from neutral diffusion

  • BOUNDARY_FORCING_HEAT_TENDENCY = heating of cell from boundary fluxes

  • FRAZIL_HEAT_TENDENCY = heating of cell from frazil

  • TH_TENDENCY_VERT_REMAP has zero vertical sum, as it redistributes heat in vertical.

  • OPOTTEMPDIFF has zero vertical sum, as it redistributes heat in the vertical.

  • BOUNDARY_FORCING_HEAT_TENDENCY generally has 3d structure, with k > 1 contributions from penetrative shortwave, and from other fluxes for the case when layers are tiny, in which case MOM6 partitions tendencies into k > 1 layers.

  • FRAZIL_HEAT_TENDENCY generally has 3d structure, since MOM6 frazil calculation checks the full ocean column.

  • FRAZIL_HEAT_TENDENCY[k=@sum] = HFSIFRAZIL = column integrated frazil heating.

  • HFDS = FRAZIL_HEAT_TENDENCY[k=@sum] + BOUNDARY_FORCING_HEAT_TENDENCY[k=@sum]

Here is an example 2d heat budget (depth summed) diagnostic for MOM.

  • OPOTTEMPTEND_2d = T_ADVECTION_XY_2d + OPOTTEMPPMDIFF_2d + HFDS

Here is an example 3d salt budget diagnostic for MOM.

  • OSALTTEND = S_ADVECTION_XY + SH_TENDENCY_VERT_REMAP + OSALTDIFF + OSALTPMDIFF * BOUNDARY_FORCING_SALT_TENDENCY

  • OSALTTEND = net tendency of salt as diagnosed in MOM.F90

  • S_ADVECTION_XY = salt convergence to cell from lateral advection

  • SH_TENDENCY_VERT_REMAP = salt convergence to cell from vertical remapping

  • OSALTDIFF = salt convergence to cell from diabatic diffusion

  • OSALTPMDIFF = salt convergence to cell from neutral diffusion

  • BOUNDARY_FORCING_SALT_TENDENCY = salt convergence to cell from boundary fluxes

  • SH_TENDENCY_VERT_REMAP has zero vertical sum, as it redistributes salt in vertical.

  • OSALTDIFF has zero vertical sum, as it redistributes salt in the vertical.

  • BOUNDARY_FORCING_SALT_TENDENCY generally has 3d structure, with k > 1 contributions from the case when layers are tiny, in which case MOM6 partitions tendencies into k > 1 layers.

  • SFDSI = BOUNDARY_FORCING_SALT_TENDENCY[k=@sum]

Here is an example 2d salt budget (depth summed) diagnostic for MOM.

  • OSALTTEND_2d = S_ADVECTION_XY_2d + OSALTPMDIFF_2d + SFDSI (+ SALT_FLUX_RESTORE)

Type Documentation

type mom/mom_diag_ids

A structure with diagnostic IDs of the state variables.

Type fields:
  • % id_u [integer,private] :: 3-d state field diagnostic IDs

  • % id_v [integer,private] :: 3-d state field diagnostic IDs

  • % id_h [integer,private] :: 3-d state field diagnostic IDs

  • % id_ssh_inst [integer,private] :: 2-d state field diagnostic ID

type mom/mom_control_struct

Control structure for the MOM module, including the variables that describe the state of the ocean.

Type fields:
  • % real (* eta_av_bc [*) :: layer thickness [H ~> m or kg m-2]

  • % real :: potential temperature [C ~> degC]

  • % real :: salinity [S ~> ppt]

  • % real :: zonal velocity component [L T-1 ~> m s-1]

  • % real :: uh = u * h * dy at u grid points [H L2 T-1 ~> m3 s-1 or kg s-1]

  • % real :: accumulated zonal thickness fluxes to advect tracers [H L2 ~> m3 or kg]

  • % real :: meridional velocity [L T-1 ~> m s-1]

  • % real :: vh = v * h * dx at v grid points [H L2 T-1 ~> m3 s-1 or kg s-1]

  • % real :: accumulated meridional thickness fluxes to advect tracers [H L2 ~> m3 or kg]

  • % real :: A running time integral of the sea surface height [T Z ~> s m].

  • % real :: time-averaged (over a forcing time step) sea surface height with a correction for the inverse barometer [Z ~> m]

  • % real :: free surface height or column mass time averaged over the last baroclinic dynamics time step [H ~> m or kg m-2]

  • % hml [real(:,:),pointer] :: active mixed layer depth [Z ~> m]

  • % time_in_cycle [real] :: The running time of the current time-stepping cycle in calls that step the dynamics, and also the length of the time integral of ssh_rint [T ~> s].

  • % time_in_thermo_cycle [real] :: The running time of the current time-stepping cycle in calls that step the thermodynamics [T ~> s].

  • % g_in [type( ocean_grid_type )] :: Input grid metric.

  • % g [type( ocean_grid_type ),pointer] :: Model grid metric.

  • % rotate_index [logical] :: True if index map is rotated.

  • % homogenize_forcings [logical] :: True if all inputs are homogenized.

  • % update_ustar [logical] :: True to update ustar from homogenized tau.

  • % gv [type( verticalgrid_type ),pointer] :: structure containing vertical grid info

  • % us [type( unit_scale_type ),pointer] :: structure containing various unit conversion factors

  • % tv [type( thermo_var_ptrs )] :: structure containing pointers to available thermodynamic fields

  • % t_dyn_rel_adv [real] :: The time of the dynamics relative to tracer advection and lateral mixing [T ~> s], or equivalently the elapsed time since advectively updating the tracers. t_dyn_rel_adv is invariably positive and may span multiple coupling timesteps.

  • % n_dyn_steps_in_adv [integer] :: The number of dynamics time steps that contributed to uhtr and vhtr since the last time tracer advection occured.

  • % t_dyn_rel_thermo [real] :: The time of the dynamics relative to diabatic processes and remapping [T ~> s]. t_dyn_rel_thermo can be negative or positive depending on whether the diabatic processes are applied before or after the dynamics and may span multiple coupling timesteps.

  • % t_dyn_rel_diag [real] :: The time of the diagnostics relative to diabatic processes and remapping [T ~> s]. t_dyn_rel_diag is always positive, since the diagnostics must lag.

  • % preadv_h_stored [logical] :: If true, the thicknesses from before the advective cycle have been stored for use in diagnostics.

  • % diag [type( diag_ctrl )] :: structure to regulate diagnostic output timing

  • % visc [type( vertvisc_type )] :: structure containing vertical viscosities, bottom drag viscosities, and related fields

  • % meke [type( meke_type )] :: Fields related to the Mesoscale Eddy Kinetic Energy.

  • % adiabatic [logical] :: If true, there are no diapycnal mass fluxes, and no calls to routines to calculate or apply diapycnal fluxes.

  • % diabatic_first [logical] :: If true, apply diabatic and thermodynamic processes before time stepping the dynamics.

  • % use_ale_algorithm [logical] :: If true, use the ALE algorithm rather than layered isopycnal/stacked shallow water mode. This logical is set by calling the function useRegridding() from the MOM_regridding module.

  • % remap_aux_vars [logical] :: If true, apply ALE remapping to all of the auxiliary 3-D variables that are needed to reproduce across restarts, similarly to what is done with the primary state variables.

  • % remap_uv_using_old_alg [logical] :: If true, use the old “remapping via a delta z” method for velocities. If false, remap between two grids described by thicknesses.

  • % stoch_eos_cs [type( mom_stoch_eos_cs )] :: structure containing random pattern for stoch EOS

  • % alternate_first_direction [logical] :: If true, alternate whether the x- or y-direction updates occur first in directionally split parts of the calculation.

  • % first_dir_restart [real] :: A real copy of Gfirst_direction for use in restart files [nondim].

  • % offline_tracer_mode [logical] :: If true,

  • % meke_in_dynamics [logical] :: If .true. (default), MEKE is called in the dynamics routine otherwise it is called during the tracer dynamics.

  • % time [type(time_type),pointer] :: pointer to the ocean clock

  • % dt [real] :: (baroclinic) dynamics time step [T ~> s]

  • % dt_therm [real] :: thermodynamics time step [T ~> s]

  • % thermo_spans_coupling [logical] :: If true, thermodynamic and tracer time steps can span multiple coupled time steps.

  • % nstep_tot [integer] :: The total number of dynamic timesteps taken so far in this run segment.

  • % count_calls [logical] :: If true, count the calls to step_MOM, rather than the number of dynamics steps in nstep_tot.

  • % debug [logical] :: If true, write verbose checksums for debugging purposes.

  • % ntrunc [integer] :: number u,v truncations since last call to write_energy

  • % cont_stencil [integer] :: The stencil for thickness from the continuity solver.

  • % do_dynamics [logical] :: If false, does not call step_MOM_dyn_*. This is an undocumented run-time flag that is fragile.

  • % split [logical] :: If true, use the split time stepping scheme.

  • % use_alt_split [logical] :: If true, use a version of the split explicit time stepping scheme that exchanges velocities with step_MOM that have the average barotropic phase over a baroclinic timestep rather than the instantaneous barotropic phase.

  • % use_rk2 [logical] :: If true, use RK2 instead of RK3 in unsplit mode (i.e., no split between barotropic and baroclinic).

  • % interface_filter [logical] :: If true, apply an interface height filter immediately after any calls to thickness_diffuse.

  • % thickness_diffuse [logical] :: If true, diffuse interface height w/ a diffusivity KHTH.

  • % thickness_diffuse_first [logical] :: If true, diffuse thickness before dynamics.

  • % mixedlayer_restrat [logical] :: If true, use submesoscale mixed layer restratifying scheme.

  • % usemeke [logical] :: If true, call the MEKE parameterization.

  • % use_stochastic_eos [logical] :: If true, use the stochastic EOS parameterizations.

  • % usewaves [logical] :: If true, update Stokes drift.

  • % use_diabatic_time_bug [logical] :: If true, uses the wrong calendar time for diabatic processes, as was done in MOM6 versions prior to February 2018.

  • % dtbt_reset_period [real] :: The time interval between dynamic recalculation of the barotropic time step [T ~> s]. If this is negative dtbt is never calculated, and if it is 0, dtbt is calculated every step.

  • % dtbt_reset_interval [type(time_type)] :: A time_time representation of dtbt_reset_period.

  • % dtbt_reset_time [type(time_type)] :: The next time DTBT should be calculated.

  • % dt_obc_seg_period [real] :: The time interval between OBC segment updates for OBGC tracers [T ~> s], or a negative value if the segment data are time-invarant, or zero to update the OBGC segment data with every call to update_OBC_segment_data.

  • % dt_obc_seg_interval [type(time_type)] :: A time_time representation of dt_obc_seg_period.

  • % dt_obc_seg_time [type(time_type)] :: The next time OBC segment update is applied to OBGC tracers.

  • % frac_shelf_h [real(:,:),pointer] :: fraction of total area occupied by ice shelf [nondim]

  • % mass_shelf [real(:,:),pointer] :: Mass of ice shelf [R Z ~> kg m-2].

  • % adp [type( accel_diag_ptrs )] :: structure containing pointers to accelerations, for derived diagnostics (e.g., energy budgets)

  • % cdp [type( cont_diag_ptrs )] :: structure containing pointers to continuity equation terms, for derived diagnostics (e.g., energy budgets)

  • % u_prev [real(:,:,:),pointer] :: previous value of u stored for diagnostics [L T-1 ~> m s-1]

  • % v_prev [real(:,:,:),pointer] :: previous value of v stored for diagnostics [L T-1 ~> m s-1]

  • % interp_p_surf [logical] :: If true, linearly interpolate surface pressure over the coupling time step, using specified value at the end of the coupling step. False by default.

  • % p_surf_prev_set [logical] :: If true, p_surf_prev has been properly set from a previous time-step or the ocean restart file. This is only valid when interp_p_surf is true.

  • % p_surf_prev [real(:,:),pointer] :: surface pressure [R L2 T-2 ~> Pa] at end previous call to step_MOM

  • % p_surf_begin [real(:,:),pointer] :: surface pressure [R L2 T-2 ~> Pa] at start of step_MOM_dyn_

  • % p_surf_end [real(:,:),pointer] :: surface pressure [R L2 T-2 ~> Pa] at end of step_MOM_dyn_

  • % write_ic [logical] :: If true, then the initial conditions will be written to file.

  • % ic_file [character (len=120)] :: A file into which the initial conditions are written in a new run if SAVE_INITIAL_CONDS is true.

  • % calc_rho_for_sea_lev [logical] :: If true, calculate rho to convert pressure to sea level.

  • % hmix [real] :: Diagnostic mixed layer thickness over which to average surface tracer properties when a bulk mixed layer is not used [H ~> m or kg m-2], or a negative value if a bulk mixed layer is being used.

  • % hfrz [real] :: If HFrz > 0, the nominal depth over which melt potential is computed [H ~> m or kg m-2]. The actual depth over which melt potential is computed is min(HFrz, OBLD), where OBLD is the boundary layer depth. If HFrz <= 0 (default), melt potential will not be computed.

  • % hmix_uv [real] :: Depth scale over which to average surface flow to feedback to the coupler/driver [H ~> m or kg m-2] when bulk mixed layer is not used, or a negative value if a bulk mixed layer is being used.

  • % check_bad_sfc_vals [logical] :: If true, scan surface state for ridiculous values.

  • % bad_val_ssh_max [real] :: Maximum SSH before triggering bad value message [Z ~> m].

  • % bad_val_sst_max [real] :: Maximum SST before triggering bad value message [C ~> degC].

  • % bad_val_sst_min [real] :: Minimum SST before triggering bad value message [C ~> degC].

  • % bad_val_sss_max [real] :: Maximum SSS before triggering bad value message [S ~> ppt].

  • % bad_val_col_thick [real] :: Minimum column thickness before triggering bad value message [Z ~> m].

  • % answer_date [integer] :: The vintage of the expressions for the surface properties. Values below 20190101 recover the answers from the end of 2018, while higher values use more appropriate expressions that differ at roundoff for non-Boussinesq cases.

  • % use_particles [logical] :: Turns on the particles package.

  • % use_uh_particles [logical] :: particles are advected by uh/h

  • % use_dbclient [logical] :: Turns on the database client used for ML inference/analysis.

  • % particle_type [character (len=10)] :: Particle types include: surface(default), profiling and sail drone.

  • % ids [type( mom_diag_ids )] :: Handles used for diagnostics.

  • % transport_ids [type( transport_diag_ids )] :: Handles used for transport diagnostics.

  • % sfc_ids [type( surface_diag_ids )] :: Handles used for surface diagnostics.

  • % diag_pre_sync [type( diag_grid_storage )] :: The grid (thicknesses) before remapping.

  • % diag_pre_dyn [type( diag_grid_storage )] :: The grid (thicknesses) before dynamics.

  • % dyn_unsplit_csp [type( mom_dyn_unsplit_cs ),pointer] :: Pointer to the control structure used for the unsplit dynamics.

  • % dyn_unsplit_rk2_csp [type( mom_dyn_unsplit_rk2_cs ),pointer] :: Pointer to the control structure used for the unsplit RK2 dynamics.

  • % dyn_split_rk2_csp [type( mom_dyn_split_rk2_cs ),pointer] :: Pointer to the control structure used for the mode-split RK2 dynamics.

  • % dyn_split_rk2b_csp [type( mom_dyn_split_rk2b_cs ),pointer] :: Pointer to the control structure used for an alternate version of the mode-split RK2 dynamics.

  • % thickness_diffuse_csp [type( thickness_diffuse_cs )] :: Pointer to the control structure used for the isopycnal height diffusive transport. This is also common referred to as Gent-McWilliams diffusion.

  • % interface_filter_csp [type( interface_filter_cs )] :: Control structure used for the interface height smoothing operator.

  • % mixedlayer_restrat_csp [type( mixedlayer_restrat_cs )] :: Pointer to the control structure used for the mixed layer restratification.

  • % set_visc_csp [type( set_visc_cs )] :: Pointer to the control structure used to set viscosities.

  • % diabatic_csp [type( diabatic_cs ),pointer] :: Pointer to the control structure for the diabatic driver.

  • % meke_csp [type( meke_cs )] :: Pointer to the control structure for the MEKE updates.

  • % varmix [type( varmix_cs )] :: Control structure for the variable mixing module.

  • % tracer_reg [type(tracer_registry_type),pointer] :: Pointer to the MOM tracer registry.

  • % tracer_adv_csp [type( tracer_advect_cs ),pointer] :: Pointer to the MOM tracer advection control structure.

  • % tracer_diff_csp [type( tracer_hor_diff_cs ),pointer] :: Pointer to the MOM along-isopycnal tracer diffusion control structure.

  • % tracer_flow_csp [type( tracer_flow_control_cs ),pointer] :: Pointer to the control structure that orchestrates the calling of tracer packages.

  • % update_obc_csp [type( update_obc_cs ),pointer] :: Pointer to the control structure for updating open boundary condition properties.

  • % obc [type( ocean_obc_type ),pointer] :: Pointer to the MOM open boundary condition type.

  • % sponge_csp [type( sponge_cs ),pointer] :: Pointer to the layered-mode sponge control structure.

  • % ale_sponge_csp [type( ale_sponge_cs ),pointer] :: Pointer to the ALE-mode sponge control structure.

  • % oda_incupd_csp [type( oda_incupd_cs ),pointer] :: Pointer to the oda incremental update control structure.

  • % int_tide_csp [type( int_tide_cs ),pointer] :: Pointer to the internal tides control structure.

  • % ale_csp [type( ale_cs ),pointer] :: Pointer to the Arbitrary Lagrangian Eulerian (ALE) vertical coordinate control structure.

  • % sum_output_csp [type( sum_output_cs ),pointer] :: Pointer to the globally summed output control structure.

  • % diagnostics_csp [type( diagnostics_cs )] :: Pointer to the MOM diagnostics control structure.

  • % offline_csp [type( offline_transport_cs ),pointer] :: Pointer to the offline tracer transport control structure.

  • % por_bar_cs [type( porous_barrier_cs )] :: Control structure for porous barrier.

  • % ensemble_ocean [logical] :: if true, this run is part of a larger ensemble for the purpose of data assimilation or statistical analysis.

  • % odacs [type( oda_cs ),pointer] :: a pointer to the control structure for handling ensemble model state vectors and data assimilation increments and priors

  • % dbcomms_cs [type( dbcomms_cs_type )] :: Control structure for database client used for online ML/AI.

  • % use_porbar [logical] :: If true, use porous barrier to constrain the widths and face areas at the edges of the grid cells.

  • % pbv [type( porous_barrier_type )] :: porous barrier fractional cell metrics

  • % particles [type(particles),pointer] :: Lagrangian particles.

  • % stoch_cs [type( stochastic_cs ),pointer] :: a pointer to the stochastics control structure

  • % restart_cs [type( mom_restart_cs ),pointer] :: Pointer to MOM’s restart control structure.

Function/Subroutine Documentation

subroutine mom/step_mom(forces_in, fluxes_in, sfc_state, Time_start, time_int_in, CS, Waves, do_dynamics, do_thermodynamics, start_cycle, end_cycle, cycle_length, reset_therm)

This subroutine orchestrates the time stepping of MOM. The adiabatic dynamics are stepped by calls to one of the step_MOM_dyn_…routines. The action of lateral processes on tracers occur in calls to advect_tracer and tracer_hordiff. Vertical mixing and possibly remapping occur inside of diabatic.

Parameters:
  • forces_in :: [inout] A structure with the driving mechanical forces

  • fluxes_in :: [inout] A structure with pointers to themodynamic, tracer and mass exchange forcing fields

  • sfc_state :: [inout] surface ocean state

  • time_start :: [in] starting time of a segment, as a time type

  • time_int_in :: [in] time interval covered by this run segment [T ~> s].

  • cs :: [inout] control structure from initialize_MOM

  • waves :: An optional pointer to a wave property CS

  • do_dynamics :: [in] Present and false, do not do updates due to the dynamics.

  • do_thermodynamics :: [in] Present and false, do not do updates due to the thermodynamics or remapping.

  • start_cycle :: [in] This indicates whether this call is to be treated as the first call to step_MOM in a time-stepping cycle; missing is like true.

  • end_cycle :: [in] This indicates whether this call is to be treated as the last call to step_MOM in a time-stepping cycle; missing is like true.

  • cycle_length :: [in] The amount of time in a coupled time stepping cycle [T ~> s].

  • reset_therm :: [in] This indicates whether the running sums of thermodynamic quantities should be reset. If missing, this is like start_cycle.

Call to:

adjust_ssh_for_p_atm mom_lateral_mixing_coeffs::calc_depth_function mom_interface_heights::calc_derived_thermo mom_lateral_mixing_coeffs::calc_resoln_function mom_error_handler::calltree_enter mom_error_handler::calltree_leave mom_error_handler::calltree_waypoint mom_diag_mediator::diag_copy_diag_to_storage extract_surface_state mom_forcing_type::homogenize_forcing mom_forcing_type::homogenize_mech_forcing id_clock_diagnostics id_clock_dynamics id_clock_ocean id_clock_other id_clock_pass mom_state_is_synchronized mom_variables::rotate_surface_state mom_oda_driver_mod::set_analysis_time mom_oda_driver_mod::set_prior_tracer step_mom_dynamics step_mom_thermo step_mom_tracer_dyn mom_stochastics::update_stochastics mom_wave_interface::update_stokes_drift

Called from:

mom6

subroutine mom/step_mom_dynamics(forces, p_surf_begin, p_surf_end, dt, dt_thermo, bbl_time_int, CS, Time_local, Waves)

Time step the ocean dynamics, including the momentum and continuity equations.

Parameters:
  • forces :: [in] A structure with the driving mechanical forces

  • p_surf_begin :: A pointer (perhaps NULL) to the surface pressure at the beginning of this dynamic step, intent in [R L2 T-2 ~> Pa].

  • p_surf_end :: A pointer (perhaps NULL) to the surface pressure at the end of this dynamic step, intent in [R L2 T-2 ~> Pa].

  • dt :: [in] time interval covered by this call [T ~> s].

  • dt_thermo :: [in] time interval covered by any updates that may span multiple dynamics steps [T ~> s].

  • bbl_time_int :: [in] time interval over which updates to the bottom boundary layer properties will apply [T ~> s], or zero not to update the properties.

  • cs :: [inout] control structure from initialize_MOM

  • time_local :: [in] End time of a segment, as a time type

  • waves :: Container for wave related parameters; the

Call to:

mom_error_handler::calltree_waypoint id_clock_bbl_visc id_clock_diagnostics id_clock_dynamics id_clock_int_filter id_clock_ml_restrat id_clock_other id_clock_pass id_clock_stoch id_clock_thick_diff id_clock_vart mom_stoch_eos::mom_calc_vart mom_stoch_eos::post_stoch_eos_diags mom_grid::set_first_direction

Called from:

step_mom

subroutine mom/step_mom_tracer_dyn(CS, G, GV, US, h, Time_local)

step_MOM_tracer_dyn does tracer advection and lateral diffusion, bringing the tracers up to date with the changes in state due to the dynamics. Surface sources and sinks and remapping are handled via step_MOM_thermo.

Parameters:
  • cs :: [inout] control structure

  • g :: [inout] ocean grid structure

  • gv :: [in] ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • h :: [in] layer thicknesses after the transports [H ~> m or kg m-2]

  • time_local :: [in] The model time at the end of the time step.

Call to:

mom_interface_heights::calc_derived_thermo mom_error_handler::calltree_waypoint id_clock_diagnostics id_clock_other id_clock_pass id_clock_thermo id_clock_tracer

Called from:

step_mom

subroutine mom/step_mom_thermo(CS, G, GV, US, u, v, h, tv, fluxes, dtdia, Time_end_thermo, update_BBL, Waves)

MOM_step_thermo orchestrates the thermodynamic time stepping and vertical remapping, via calls to diabatic (or adiabatic) and ALE_regrid.

Parameters:
  • cs :: [inout] Master MOM control structure

  • g :: [inout] ocean grid structure

  • gv :: [inout] ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • u :: [inout] zonal velocity [L T-1 ~> m s-1]

  • v :: [inout] meridional velocity [L T-1 ~> m s-1]

  • h :: [inout] layer thickness [H ~> m or kg m-2]

  • tv :: [inout] A structure pointing to various thermodynamic variables

  • fluxes :: [inout] pointers to forcing fields

  • dtdia :: [in] The time interval over which to advance [T ~> s]

  • time_end_thermo :: [in] End of averaging interval for thermo diags

  • update_bbl :: [in] If true, calculate the bottom boundary layer properties.

  • waves :: Container for wave related parameters

Call to:

mom_ale::ale_update_regrid_weights mom_oda_driver_mod::apply_oda_tracer_increments mom_interface_heights::calc_derived_thermo mom_error_handler::calltree_enter mom_error_handler::calltree_leave mom_error_handler::calltree_waypoint id_clock_adiabatic id_clock_ale id_clock_bbl_visc id_clock_diabatic id_clock_pass id_clock_thermo mom_particles_mod::particles_to_k_space mom_particles_mod::particles_to_z_space mom_ale::pre_ale_diagnostics mom_dynamics_split_rk2::remap_dyn_split_rk2_aux_vars mom_dynamics_split_rk2b::remap_dyn_split_rk2b_aux_vars

Called from:

step_mom

subroutine mom/step_offline(forces, fluxes, sfc_state, Time_start, time_interval, CS)

step_offline is the main driver for running tracers offline in MOM6. This has been primarily developed with ALE configurations in mind. Some work has been done in isopycnal configuration, but the work is very preliminary. Some more detail about this capability along with some of the subroutines called here can be found in tracers/MOM_offline_control.F90

Parameters:
  • forces :: [in] A structure with the driving mechanical forces

  • fluxes :: [inout] pointers to forcing fields

  • sfc_state :: [inout] surface ocean state

  • time_start :: [in] starting time of a segment, as a time type

  • time_interval :: [in] time interval [T ~> s]

  • cs :: [inout] control structure from initialize_MOM

Call to:

adjust_ssh_for_p_atm mom_lateral_mixing_coeffs::calc_depth_function mom_interface_heights::calc_derived_thermo mom_lateral_mixing_coeffs::calc_resoln_function extract_surface_state id_clock_ale id_clock_offline_tracer mom_offline_main::offline_advection_layer

Called from:

mom6

subroutine mom/initialize_mom(Time, Time_init, param_file, dirs, CS, Time_in, offline_tracer_mode, input_restart_file, diag_ptr, count_calls, tracer_flow_CSp, ice_shelf_CSp, waves_CSp, ensemble_num)

Initialize MOM, including memory allocation, setting up parameters and diagnostics, initializing the ocean state variables, and initializing subsidiary modules.

Parameters:
  • time :: [inout] model time, set in this routine

  • time_init :: [in] The start time for the coupled model’s calendar

  • param_file :: [out] structure indicating parameter file to parse

  • dirs :: [out] structure with directory paths

  • cs :: [inout] pointer set in this routine to MOM control structure

  • time_in :: [in] time passed to MOM_initialize_state when model is not being started from a restart file

  • offline_tracer_mode :: [out] True is returned if tracers are being run offline

  • input_restart_file :: [in] If present, name of restart file to read

  • diag_ptr :: A pointer set in this routine to the diagnostic regulatory structure

  • tracer_flow_csp :: A pointer set in this routine to

  • count_calls :: [in] If true, nstep_tot counts the number of calls to step_MOM instead of the number of dynamics timesteps.

  • ice_shelf_csp :: A pointer to an ice shelf control structure

  • waves_csp :: An optional pointer to a wave property CS

  • ensemble_num :: Ensemble index provided by the cap (instead of FMS ensemble manager)

Call to:

mom_ale::ale_register_diags mom_interface_heights::calc_derived_thermo mom_boundary_update::call_obc_register mom_tracer_flow_control::call_tracer_register_obc_segments mom_error_handler::calltree_enter mom_error_handler::calltree_leave mom_error_handler::calltree_waypoint mom_check_scaling::check_mom6_scaling_factors mom_transcribe_grid::copy_dyngrid_to_mom_grid mom_transcribe_grid::copy_mom_grid_to_dyngrid mom_database_comms::database_comms_init mom_restart::determine_is_new_run mom_diag_mediator::diag_copy_diag_to_storage mom_eos::eos_init mom_obsolete_params::find_obsolete_params mom_get_input::get_mom_input mom_verticalgrid::get_tr_flux_units mom_ice_shelf::ice_shelf_query id_clock_init id_clock_mom_init id_clock_pass_init id_clock_unit_tests mom_ale_sponge::init_ale_sponge_diags mom_oda_incupd::init_oda_incupd_diags mom_sponge::init_sponge_diags mom_ice_shelf::initialize_ice_shelf mom_restart::is_new_run mom_tracer_registry::lock_tracer_registry mom_meke::meke_init mom_mixed_layer_restrat::mixedlayer_restrat_register_restarts mom_coord_initialization::mom_initialize_coord mom_fixed_initialization::mom_initialize_fixed mom_state_initialization::mom_initialize_state mom_timing_init mom_diabatic_driver::register_diabatic_restarts register_diags mom_obsolete_diagnostics::register_obsolete_diagnostics mom_restart::restart_init mom_dyn_horgrid::rotate_dyn_horgrid mom_hor_index::rotate_hor_index rotate_initial_state mom_open_boundary::rotate_obc_config mom_open_boundary::rotate_obc_init mom_grid::set_first_direction set_restart_fields mom_set_visc::set_visc_init mom_stoch_eos::stoch_eos_register_restarts mom_stochastics::stochastics_init mom_unit_scaling::unit_scaling_init mom_unit_tests::unit_tests mom_shared_initialization::write_ocean_geometry_file mom_coord_initialization::write_vertgrid_file

Called from:

mom6

subroutine mom/finish_mom_initialization(Time, dirs, CS)

Finishes initializing MOM and writes out the initial conditions.

Parameters:
  • time :: [in] model time, used in this routine

  • dirs :: [in] structure with directory paths

  • cs :: [inout] MOM control structure

Call to:

mom_error_handler::calltree_enter mom_error_handler::calltree_leave id_clock_init

Called from:

mom6

subroutine mom/register_diags(Time, G, GV, US, IDs, diag)

Register certain diagnostics.

Parameters:
  • time :: [in] current model time

  • g :: [in] ocean grid structure

  • gv :: [in] ocean vertical grid structure

  • us :: [inout] A dimensional unit scaling type

  • ids :: [inout] A structure with the diagnostic IDs.

  • diag :: [inout] regulates diagnostic output

Call to:

mom_verticalgrid::get_thickness_units

Called from:

initialize_mom

subroutine mom/mom_timing_init(CS)

Set up CPU clock IDs for timing various subroutines.

Parameters:

cs :: [in] control structure set up by initialize_MOM.

Call to:

id_clock_adiabatic id_clock_ale id_clock_bbl_visc id_clock_continuity id_clock_diabatic id_clock_diagnostics id_clock_dynamics id_clock_int_filter id_clock_ml_restrat id_clock_mom_init id_clock_ocean id_clock_offline_tracer id_clock_other id_clock_pass id_clock_pass_init id_clock_stoch id_clock_thermo id_clock_thick_diff id_clock_tracer id_clock_vart id_clock_z_diag

Called from:

initialize_mom

subroutine mom/set_restart_fields(GV, US, param_file, CS, restart_CSp)

Set the fields that are needed for bitwise identical restarting the time stepping scheme. In addition to those specified here directly, there may be fields related to the forcing or to the barotropic solver that are needed; these are specified in sub- routines that are called from this one.

This routine should be altered if there are any changes to the time stepping scheme. The CHECK_RESTART facility may be used to confirm that all needed restart fields have been included.

Parameters:
  • gv :: [inout] ocean vertical grid structure

  • us :: [inout] A dimensional unit scaling type

  • param_file :: [in] opened file for parsing to get parameters

  • cs :: [in] control structure set up by initialize_MOM

  • restart_csp :: pointer to the restart control structure that will be used for MOM.

Call to:

mom_verticalgrid::get_flux_units mom_verticalgrid::get_thickness_units

Called from:

initialize_mom

subroutine mom/adjust_ssh_for_p_atm(tv, G, GV, US, ssh, p_atm, use_EOS)

Apply a correction to the sea surface height to compensate for the atmospheric pressure (the inverse barometer).

Parameters:
  • tv :: [in] A structure pointing to various thermodynamic variables

  • g :: [in] ocean grid structure

  • gv :: [in] ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • ssh :: [inout] time mean surface height [Z ~> m]

  • p_atm :: Ocean surface pressure [R L2 T-2 ~> Pa]

  • use_eos :: [in] If true, calculate the density for the SSH correction using the equation of state.

Call to:

mom_eos::eos_domain

Called from:

step_mom step_offline

subroutine mom/extract_surface_state(CS, sfc_state_in)

Set the surface (return) properties of the ocean model by setting the appropriate fields in sfc_state. Unused fields are set to NULL or are unallocated.

Parameters:
  • cs :: [inout] Master MOM control structure

  • sfc_state_in :: [inout] transparent ocean surface state structure shared with the calling routine data in this structure is intent out.

Call to:

mom_error_handler::calltree_enter mom_error_handler::calltree_leave mom_checksum_packages::mom_surface_chksum mom_variables::rotate_surface_state

Called from:

mom6 step_mom step_offline

subroutine mom/rotate_initial_state(u_in, v_in, h_in, T_in, S_in, use_temperature, turns, u, v, h, T, S)

Rotate initialization fields from input to rotated arrays.

Parameters:
  • u_in :: [in] Zonal velocity on the initial grid [L T-1 ~> m s-1]

  • v_in :: [in] Meridional velocity on the initial grid [L T-1 ~> m s-1]

  • h_in :: [in] Layer thickness on the initial grid [H ~> m or kg m-2]

  • t_in :: [in] Temperature on the initial grid [C ~> degC]

  • s_in :: [in] Salinity on the initial grid [S ~> ppt]

  • use_temperature :: [in] If true, temperature and salinity are active

  • turns :: [in] The number quarter-turns to apply

  • u :: [out] Zonal velocity on the rotated grid [L T-1 ~> m s-1]

  • v :: [out] Meridional velocity on the rotated grid [L T-1 ~> m s-1]

  • h :: [out] Layer thickness on the rotated grid [H ~> m or kg m-2]

  • t :: [out] Temperature on the rotated grid [C ~> degC]

  • s :: [out] Salinity on the rotated grid [S ~> ppt]

Called from:

initialize_mom

function mom/mom_state_is_synchronized(CS, adv_dyn) [logical]

Return true if all phases of step_MOM are at the same point in time.

Parameters:
  • cs :: [inout] MOM control structure

  • adv_dyn :: [in] If present and true, only check whether the advection is up-to-date with the dynamics.

Return:

undefined :: True if all phases of the update are synchronized.

Called from:

mom6 step_mom

subroutine mom/get_mom_state_elements(CS, G, GV, US, C_p, C_p_scaled, use_temp)

This subroutine offers access to values or pointers to other types from within the MOM_control_struct, allowing the MOM_control_struct to be opaque.

Parameters:
  • cs :: [inout] MOM control structure

  • g :: structure containing metrics and grid info

  • gv :: structure containing vertical grid info

  • us :: A dimensional unit scaling type

  • c_p :: [out] The heat capacity [J kg degC-1]

  • c_p_scaled :: [out] The heat capacity in scaled units [Q C-1 ~> J kg-1 degC-1]

  • use_temp :: [out] True if temperature is a state variable

Called from:

mom6

subroutine mom/get_ocean_stocks(CS, mass, heat, salt, on_PE_only)

Find the global integrals of various quantities.

Parameters:
  • cs :: [inout] MOM control structure

  • heat :: [out] The globally integrated integrated ocean heat [J].

  • salt :: [out] The globally integrated integrated ocean salt [kg].

  • mass :: [out] The globally integrated integrated ocean mass [kg].

  • on_pe_only :: [in] If present and true, only sum on the local PE.

Call to:

mom_spatial_means::global_mass_integral

subroutine mom/save_mom_restart(CS, directory, time, G, time_stamped, filename, GV, num_rest_files, write_IC)

Save restart/pickup files required to initialize the MOM6 internal state.

Parameters:
  • cs :: [inout] MOM control structure

  • directory :: [in] The directory where the restart files are to be written

  • time :: [in] The current model time

  • g :: [inout] The ocean’s grid structure

  • time_stamped :: [in] If present and true, add time-stamp to the restart file names

  • filename :: [in] A filename that overrides the name in CSrestartfile

  • gv :: [in] The ocean’s vertical grid structure

  • num_rest_files :: [out] number of restart files written

  • write_ic :: [in] If present and true, initial conditions are being written

Called from:

mom6 ocean_model_mod::ocean_model_restart ocean_model_mod::ocean_model_save_restart

subroutine mom/mom_end(CS)

End of ocean model, including memory deallocation.

Parameters:

cs :: [inout] MOM control structure

Call to:

mom_meke::meke_end mom_diagnostics::mom_diagnostics_end mom_boundary_update::obc_register_end mom_offline_main::offline_transport_end mom_set_visc::set_visc_end mom_thickness_diffuse::thickness_diffuse_end mom_tracer_advect::tracer_advect_end mom_tracer_flow_control::tracer_flow_control_end mom_tracer_hor_diff::tracer_hor_diff_end mom_tracer_registry::tracer_registry_end mom_unit_scaling::unit_scaling_end

Called from:

mom6