rossby_front_2d_initialization module reference

Initial conditions for the 2D Rossby front test.

More…

Functions/Subroutines

rossby_front_initialize_thickness()

Initialization of thicknesses in 2D Rossby front test.

rossby_front_initialize_temperature_salinity()

Initialization of temperature and salinity in the Rossby front test.

rossby_front_initialize_velocity()

Initialization of u and v in the Rossby front test.

ypseudo()

Pseudo coordinate across domain used by Hml() and dTdy() returns a coordinate from -PI/2 ..

hml()

Analytic prescription of mixed layer depth in 2d Rossby front test, in the same units as max_depth (usually [Z ~> m] or [H ~> m or kg m-2])

dtdy()

Analytic prescription of mixed layer temperature gradient in [C L-1 ~> degC m-1] in 2d Rossby front test.

Detailed Description

Description of the 2d Rossby front initial conditions

Consistent with a linear equation of state, the system has a constant stratification below the mixed layer, stratified in temperature only. Isotherms are flat below the mixed layer and vertical within. Salinity is constant. The mixed layer has a half sine form so that there are no mixed layer or temperature gradients at the side walls.

Below the mixed layer the potential temperature, \(\theta(z)\), is given by

\[\theta(z) = \theta_0 - \Delta \theta \left( z + h_{ML} \right)\]

where \(\theta_0\) and \(\Delta \theta\) are external model parameters.

The depth of the mixed layer, \(H_{ML}\) is

\[h_{ML}(y) = h_{min} + \left( h_{max} - h_{min} \right) \cos{\pi y/L}\]

. The temperature in mixed layer is given by the reference temperature at \(z=h_{ML}\) so that

\[\begin{eqnarray} \theta(y,z) = \theta_0 - \Delta \theta \left( z + h_{ML} \right) & \forall & z < h_{ML}(y) T.B.D. \end{eqnarray}\]

Function/Subroutine Documentation

subroutine rossby_front_2d_initialization/rossby_front_initialize_thickness(h, G, GV, US, param_file, just_read)

Initialization of thicknesses in 2D Rossby front test.

Parameters:
  • g :: [in] Grid structure

  • gv :: [in] Vertical grid structure

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

  • h :: [out] The thickness that is being initialized [H ~> m or kg m-2]

  • param_file :: [in] A structure indicating the open file to parse for model parameter values.

  • just_read :: [in] If true, this call will only read parameters without changing h.

Call to:

hml mdl mom_error_handler::mom_error mom_error_handler::mom_mesg regrid_consts::regridding_rho regrid_consts::regridding_sigma

subroutine rossby_front_2d_initialization/rossby_front_initialize_temperature_salinity(T, S, h, G, GV, US, param_file, just_read)

Initialization of temperature and salinity in the Rossby front test.

Parameters:
  • g :: [in] Grid structure

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

  • t :: [out] Potential temperature [C ~> degC]

  • s :: [out] Salinity [S ~> ppt]

  • h :: [in] Thickness [H ~> m or kg m-2]

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

  • param_file :: [in] Parameter file handle

  • just_read :: [in] If true, this call will only read parameters without changing T & S.

Call to:

hml mdl mom_error_handler::mom_error

subroutine rossby_front_2d_initialization/rossby_front_initialize_velocity(u, v, h, G, GV, US, param_file, just_read)

Initialization of u and v in the Rossby front test.

Parameters:
  • g :: [in] Grid structure

  • gv :: [in] Vertical grid structure

  • u :: [out] i-component of velocity [L T-1 ~> m s-1]

  • v :: [out] j-component of velocity [L T-1 ~> m s-1]

  • h :: [in] Thickness [H ~> m or kg m-2]

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

  • param_file :: [in] A structure indicating the open file to parse for model parameter values.

  • just_read :: [in] If present and true, this call will only read parameters without setting u & v.

Call to:

dtdy hml mdl mom_error_handler::mom_error

Called from:

mom_state_initialization::mom_initialize_state

function rossby_front_2d_initialization/ypseudo(G, lat) [real]

Pseudo coordinate across domain used by Hml() and dTdy() returns a coordinate from -PI/2 .. PI/2 squashed towards the center of the domain [radians].

Parameters:
  • g :: [in] Grid structure

  • lat :: [in] Latitude in arbitrary units, often [km]

Call to:

frontfractionalwidth

Called from:

dtdy hml

function rossby_front_2d_initialization/hml(G, lat, max_depth) [real]

Analytic prescription of mixed layer depth in 2d Rossby front test, in the same units as max_depth (usually [Z ~> m] or [H ~> m or kg m-2])

Parameters:
  • g :: [in] Grid structure

  • lat :: [in] Latitude in arbitrary units, often [km]

  • max_depth :: [in] The maximum depth of the ocean [Z ~> m] or [H ~> m or kg m-2]

Call to:

hmlmax hmlmin ypseudo

Called from:

rossby_front_initialize_temperature_salinity rossby_front_initialize_thickness rossby_front_initialize_velocity

function rossby_front_2d_initialization/dtdy(G, dT, lat, US) [real]

Analytic prescription of mixed layer temperature gradient in [C L-1 ~> degC m-1] in 2d Rossby front test.

Parameters:
  • g :: [in] Grid structure

  • dt :: [in] Top to bottom temperature difference [C ~> degC]

  • lat :: [in] Latitude in [km]

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

Call to:

frontfractionalwidth hmlmax hmlmin ypseudo

Called from:

rossby_front_initialize_velocity