MOM6
|
Initial conditions for the 2D Rossby front test.
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
\[ \theta(y,z) = \theta_0 - \Delta \theta \left( z + h_{ML} \right) & \forall & z < h_{ML}(y) T.B.D. \]
Functions/Subroutines | |
subroutine, public | rossby_front_initialize_thickness (h, G, GV, param_file, just_read_params) |
Initialization of thicknesses in 2D Rossby front test. More... | |
subroutine, public | rossby_front_initialize_temperature_salinity (T, S, h, G, param_file, eqn_of_state, just_read_params) |
Initialization of temperature and salinity in the Rossby front test. More... | |
subroutine, public | rossby_front_initialize_velocity (u, v, h, G, GV, param_file, just_read_params) |
Initialization of u and v in the Rossby front test. More... | |
real function | ypseudo (G, lat) |
Pseudo coordinate across domain used by Hml() and dTdy() returns a coordinate from -PI/2 .. PI/2 squashed towards the center of the domain. More... | |
real function | hml (G, lat) |
Analytic prescription of mixed layer depth in 2d Rossby front test. More... | |
real function | dtdy (G, dT, lat) |
Analytic prescription of mixed layer temperature gradient in 2d Rossby front test. More... | |
Variables | |
character(len=40) | mdl = "Rossby_front_2d_initialization" |
This module's name. More... | |
real, parameter | frontfractionalwidth = 0.5 |
Width of front as fraction of domain. More... | |
real, parameter | hmlmin = 0.25 |
Shallowest ML as fractional depth of ocean. More... | |
real, parameter | hmlmax = 0.75 |
Deepest ML as fractional depth of ocean. More... | |
|
private |
Analytic prescription of mixed layer temperature gradient in 2d Rossby front test.
[in] | g | Grid structure |
[in] | dt | Top to bottom temperature difference |
[in] | lat | Latitude |
Definition at line 247 of file Rossby_front_2d_initialization.F90.
References frontfractionalwidth, hmlmax, hmlmin, and ypseudo().
Referenced by rossby_front_initialize_velocity().
|
private |
Analytic prescription of mixed layer depth in 2d Rossby front test.
[in] | g | Grid structure |
[in] | lat | Latitude |
Definition at line 234 of file Rossby_front_2d_initialization.F90.
References hmlmax, hmlmin, and ypseudo().
Referenced by rossby_front_initialize_temperature_salinity(), rossby_front_initialize_thickness(), and rossby_front_initialize_velocity().
subroutine, public rossby_front_2d_initialization::rossby_front_initialize_temperature_salinity | ( | real, dimension(szi_(g),szj_(g), szk_(g)), intent(out) | T, |
real, dimension(szi_(g),szj_(g), szk_(g)), intent(out) | S, | ||
real, dimension(szi_(g),szj_(g), szk_(g)), intent(in) | h, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(param_file_type), intent(in) | param_file, | ||
type(eos_type), pointer | eqn_of_state, | ||
logical, intent(in), optional | just_read_params | ||
) |
Initialization of temperature and salinity in the Rossby front test.
[in] | g | Grid structure |
[out] | t | Potential temperature [deg C] |
[out] | s | Salinity [ppt] |
[in] | h | Thickness |
[in] | param_file | Parameter file handle |
eqn_of_state | Equation of state structure | |
[in] | just_read_params | If present and true, this call will only read parameters without changing h. |
Definition at line 114 of file Rossby_front_2d_initialization.F90.
subroutine, public rossby_front_2d_initialization::rossby_front_initialize_thickness | ( | real, dimension(szi_(g),szj_(g),szk_(gv)), intent(out) | h, |
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(param_file_type), intent(in) | param_file, | ||
logical, intent(in), optional | just_read_params | ||
) |
Initialization of thicknesses in 2D Rossby front test.
[in] | g | Grid structure |
[in] | gv | Vertical grid structure |
[out] | h | The thickness that is being initialized, in m. |
[in] | param_file | A structure indicating the open file to parse for model parameter values. |
[in] | just_read_params | If present and true, this call will only read parameters without changing h. |
Definition at line 42 of file Rossby_front_2d_initialization.F90.
References hml(), mdl, mom_error_handler::mom_error(), mom_error_handler::mom_mesg(), regrid_consts::regridding_rho, and regrid_consts::regridding_sigma.
subroutine, public rossby_front_2d_initialization::rossby_front_initialize_velocity | ( | real, dimension(szib_(g),szj_(g),szk_(g)), intent(out) | u, |
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(out) | v, | ||
real, dimension(szi_(g),szj_(g), szk_(g)), intent(in) | h, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(param_file_type), intent(in) | param_file, | ||
logical, intent(in), optional | just_read_params | ||
) |
Initialization of u and v in the Rossby front test.
[in] | g | Grid structure |
[in] | gv | Vertical grid structure |
[out] | u | i-component of velocity [m/s] |
[out] | v | j-component of velocity [m/s] |
[in] | h | Thickness [H] |
[in] | param_file | A structure indicating the open file to parse for model parameter values. |
[in] | just_read_params | If present and true, this call will only read parameters without changing h. |
Definition at line 165 of file Rossby_front_2d_initialization.F90.
References dtdy(), hml(), and mdl.
Referenced by mom_state_initialization::mom_initialize_state().
|
private |
Pseudo coordinate across domain used by Hml() and dTdy() returns a coordinate from -PI/2 .. PI/2 squashed towards the center of the domain.
[in] | g | Grid structure |
[in] | lat | Latitude |
Definition at line 221 of file Rossby_front_2d_initialization.F90.
References frontfractionalwidth.
Referenced by dtdy(), and hml().
|
private |
Width of front as fraction of domain.
Definition at line 34 of file Rossby_front_2d_initialization.F90.
|
private |
Deepest ML as fractional depth of ocean.
Definition at line 36 of file Rossby_front_2d_initialization.F90.
|
private |
Shallowest ML as fractional depth of ocean.
Definition at line 35 of file Rossby_front_2d_initialization.F90.
|
private |
This module's name.
Definition at line 25 of file Rossby_front_2d_initialization.F90.
Referenced by rossby_front_initialize_temperature_salinity(), rossby_front_initialize_thickness(), and rossby_front_initialize_velocity().