MOM6
rossby_front_2d_initialization Module Reference

Detailed Description

Initial conditions for the 2D Rossby front test.

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

\[ \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...
 

Function/Subroutine Documentation

◆ dtdy()

real function rossby_front_2d_initialization::dtdy ( type(ocean_grid_type), intent(in)  G,
real, intent(in)  dT,
real, intent(in)  lat 
)
private

Analytic prescription of mixed layer temperature gradient in 2d Rossby front test.

Parameters
[in]gGrid structure
[in]dtTop to bottom temperature difference
[in]latLatitude

Definition at line 247 of file Rossby_front_2d_initialization.F90.

References frontfractionalwidth, hmlmax, hmlmin, and ypseudo().

Referenced by rossby_front_initialize_velocity().

247  type(ocean_grid_type), intent(in) :: g !< Grid structure
248  real, intent(in) :: dt !< Top to bottom temperature difference
249  real, intent(in) :: lat !< Latitude
250  ! Local
251  real :: pi, dhml, dhdy
252  real :: km = 1.e3 ! AXIS_UNITS = 'k' (1000 m)
253 
254  pi = 4.0 * atan(1.0)
255  dhml = 0.5 * ( hmlmax - hmlmin ) * g%max_depth
256  dhdy = dhml * ( pi / ( frontfractionalwidth * g%len_lat * km ) ) * cos( ypseudo(g, lat) )
257  dtdy = -( dt / g%max_depth ) * dhdy
258 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ hml()

real function rossby_front_2d_initialization::hml ( type(ocean_grid_type), intent(in)  G,
real, intent(in)  lat 
)
private

Analytic prescription of mixed layer depth in 2d Rossby front test.

Parameters
[in]gGrid structure
[in]latLatitude

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

234  type(ocean_grid_type), intent(in) :: g !< Grid structure
235  real, intent(in) :: lat !< Latitude
236  ! Local
237  real :: dhml, hmlmean
238 
239  dhml = 0.5 * ( hmlmax - hmlmin ) * g%max_depth
240  hmlmean = 0.5 * ( hmlmin + hmlmax ) * g%max_depth
241  hml = hmlmean + dhml * sin( ypseudo(g, lat) )
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ rossby_front_initialize_temperature_salinity()

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.

Parameters
[in]gGrid structure
[out]tPotential temperature [deg C]
[out]sSalinity [ppt]
[in]hThickness
[in]param_fileParameter file handle
eqn_of_stateEquation of state structure
[in]just_read_paramsIf present and true, this call will only read parameters without changing h.

Definition at line 114 of file Rossby_front_2d_initialization.F90.

References hml(), and mdl.

114  type(ocean_grid_type), intent(in) :: g !< Grid structure
115  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: t !< Potential temperature [deg C]
116  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: s !< Salinity [ppt]
117  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Thickness
118  type(param_file_type), intent(in) :: param_file !< Parameter file handle
119  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
120  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
121  !! only read parameters without changing h.
122 
123  integer :: i, j, k, is, ie, js, je, nz
124  real :: t_ref, s_ref ! Reference salinity and temerature within surface layer
125  real :: t_range ! Range of salinities and temperatures over the vertical
126  real :: y, zc, zi, dtdz
127  logical :: just_read ! If true, just read parameters but set nothing.
128  character(len=40) :: verticalcoordinate
129  real :: pi ! 3.1415926... calculated as 4*atan(1)
130 
131  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
132 
133  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
134 
135  call get_param(param_file, mdl,"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
136  default=default_coordinate_mode, do_not_log=just_read)
137  call get_param(param_file, mdl,"S_REF",s_ref,'Reference salinity', units='1e-3', &
138  fail_if_missing=.not.just_read, do_not_log=just_read)
139  call get_param(param_file, mdl,"T_REF",t_ref,'Reference temperature',units='C',&
140  fail_if_missing=.not.just_read, do_not_log=just_read)
141  call get_param(param_file, mdl,"T_RANGE",t_range,'Initial temperature range',&
142  units='C', default=0.0, do_not_log=just_read)
143 
144  if (just_read) return ! All run-time parameters have been read, so return.
145 
146  t(:,:,:) = 0.0
147  s(:,:,:) = s_ref
148  dtdz = t_range / g%max_depth
149 
150  do j = g%jsc,g%jec ; do i = g%isc,g%iec
151  zi = 0.
152  do k = 1, nz
153  zi = zi - h(i,j,k) ! Bottom interface position
154  zc = zi - 0.5*h(i,j,k) ! Position of middle of cell
155  zc = min( zc, -hml(g, g%geoLatT(i,j)) ) ! Bound by depth of mixed layer
156  t(i,j,k) = t_ref + dtdz * zc ! Linear temperature profile
157  enddo
158  end do ; end do
159 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:

◆ rossby_front_initialize_thickness()

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.

Parameters
[in]gGrid structure
[in]gvVertical grid structure
[out]hThe thickness that is being initialized, in m.
[in]param_fileA structure indicating the open file to parse for model parameter values.
[in]just_read_paramsIf 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.

42  type(ocean_grid_type), intent(in) :: g !< Grid structure
43  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
44  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
45  intent(out) :: h !< The thickness that is being initialized, in m.
46  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
47  !! to parse for model parameter values.
48  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
49  !! only read parameters without changing h.
50 
51  integer :: i, j, k, is, ie, js, je, nz
52  real :: tz, dml, eta, stretch, h0
53  real :: min_thickness, t_range, drho_dt
54  logical :: just_read ! If true, just read parameters but set nothing.
55  character(len=40) :: verticalcoordinate
56 
57  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
58 
59  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
60 
61  if (.not.just_read) &
62  call mom_mesg("Rossby_front_2d_initialization.F90, Rossby_front_initialize_thickness: setting thickness")
63 
64  if (.not.just_read) call log_version(param_file, mdl, version, "")
65  ! Read parameters needed to set thickness
66  call get_param(param_file, mdl, "MIN_THICKNESS", min_thickness, &
67  'Minimum layer thickness',units='m',default=1.e-3, do_not_log=just_read)
68  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
69  default=default_coordinate_mode, do_not_log=just_read)
70  call get_param(param_file, mdl, "T_RANGE", t_range, 'Initial temperature range', &
71  units='C', default=0.0, do_not_log=just_read)
72  call get_param(param_file, mdl, "DRHO_DT", drho_dt, default=-0.2, do_not_log=.true.)
73 
74  if (just_read) return ! All run-time parameters have been read, so return.
75 
76  tz = t_range / g%max_depth
77 
78  select case ( coordinatemode(verticalcoordinate) )
79 
80  case (regridding_layer, regridding_rho)
81  do j = g%jsc,g%jec ; do i = g%isc,g%iec
82  dml = hml( g, g%geoLatT(i,j) )
83  eta = -( -drho_dt / gv%Rho0 ) * tz * 0.5 * ( dml * dml )
84  stretch = ( ( g%max_depth + eta ) / g%max_depth )
85  h0 = ( g%max_depth / real(nz) ) * stretch
86  do k = 1, nz
87  h(i,j,k) = h0
88  enddo
89  end do ; end do
90 
91  case (regridding_zstar, regridding_sigma)
92  do j = g%jsc,g%jec ; do i = g%isc,g%iec
93  dml = hml( g, g%geoLatT(i,j) )
94  eta = -( -drho_dt / gv%Rho0 ) * tz * 0.5 * ( dml * dml )
95  stretch = ( ( g%max_depth + eta ) / g%max_depth )
96  h0 = ( g%max_depth / real(nz) ) * stretch
97  do k = 1, nz
98  h(i,j,k) = h0
99  enddo
100  end do ; end do
101 
102  case default
103  call mom_error(fatal,"Rossby_front_initialize: "// &
104  "Unrecognized i.c. setup - set REGRIDDING_COORDINATE_MODE")
105 
106  end select
107 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:

◆ rossby_front_initialize_velocity()

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.

Parameters
[in]gGrid structure
[in]gvVertical grid structure
[out]ui-component of velocity [m/s]
[out]vj-component of velocity [m/s]
[in]hThickness [H]
[in]param_fileA structure indicating the open file to parse for model parameter values.
[in]just_read_paramsIf 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().

165  type(ocean_grid_type), intent(in) :: g !< Grid structure
166  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
167  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out) :: u !< i-component of velocity [m/s]
168  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: v !< j-component of velocity [m/s]
169  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Thickness [H]
170  type(param_file_type), intent(in) :: param_file !< A structure indicating the
171  !! open file to parse for model
172  !! parameter values.
173  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
174  !! only read parameters without changing h.
175 
176  real :: y ! Non-dimensional coordinate across channel, 0..pi
177  real :: t_range ! Range of salinities and temperatures over the vertical
178  real :: dudt ! Factor to convert dT/dy into dU/dz, g*alpha/f
179  real :: drho_dt, zi, zc, zm, f, ty, dml, hatu
180  integer :: i, j, k, is, ie, js, je, nz
181  logical :: just_read ! If true, just read parameters but set nothing.
182  character(len=40) :: verticalcoordinate
183 
184  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
185 
186  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
187 
188  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
189  default=default_coordinate_mode, do_not_log=just_read)
190  call get_param(param_file, mdl, "T_RANGE", t_range, 'Initial temperature range', &
191  units='C', default=0.0, do_not_log=just_read)
192  call get_param(param_file, mdl, "DRHO_DT", drho_dt, default=-0.2, do_not_log=.true.)
193 
194  if (just_read) return ! All run-time parameters have been read, so return.
195 
196  v(:,:,:) = 0.0
197  u(:,:,:) = 0.0
198 
199  do j = g%jsc,g%jec ; do i = g%isc-1,g%iec+1
200  f = 0.5*( g%CoriolisBu(i,j) + g%CoriolisBu(i,j-1) )
201  dudt = 0.0 ; if (abs(f) > 0.0) &
202  dudt = ( gv%g_Earth * drho_dt ) / ( f * gv%Rho0 )
203  dml = hml( g, g%geoLatT(i,j) )
204  ty = dtdy( g, t_range, g%geoLatT(i,j) )
205  zi = 0.
206  do k = 1, nz
207  hatu = 0.5*(h(i,j,k)+h(i+1,j,k))
208  zi = zi - hatu ! Bottom interface position
209  zc = zi - 0.5*hatu ! Position of middle of cell
210  zm = max( zc + dml, 0. ) ! Height above bottom of mixed layer
211  u(i,j,k) = dudt * ty * zm ! Thermal wind starting at base of ML
212  enddo
213  enddo ; enddo
214 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ypseudo()

real function rossby_front_2d_initialization::ypseudo ( type(ocean_grid_type), intent(in)  G,
real, intent(in)  lat 
)
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.

Parameters
[in]gGrid structure
[in]latLatitude

Definition at line 221 of file Rossby_front_2d_initialization.F90.

References frontfractionalwidth.

Referenced by dtdy(), and hml().

221  type(ocean_grid_type), intent(in) :: g !< Grid structure
222  real, intent(in) :: lat !< Latitude
223  ! Local
224  real :: y, pi
225 
226  pi = 4.0 * atan(1.0)
227  ypseudo = ( ( lat - g%south_lat ) / g%len_lat ) - 0.5 ! -1/2 .. 1/.2
228  ypseudo = pi * max(-0.5, min(0.5, ypseudo / frontfractionalwidth))
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

Variable Documentation

◆ frontfractionalwidth

real, parameter rossby_front_2d_initialization::frontfractionalwidth = 0.5
private

Width of front as fraction of domain.

Definition at line 34 of file Rossby_front_2d_initialization.F90.

Referenced by dtdy(), and ypseudo().

34 real, parameter :: frontfractionalwidth = 0.5 !< Width of front as fraction of domain

◆ hmlmax

real, parameter rossby_front_2d_initialization::hmlmax = 0.75
private

Deepest ML as fractional depth of ocean.

Definition at line 36 of file Rossby_front_2d_initialization.F90.

Referenced by dtdy(), and hml().

36 real, parameter :: hmlmax = 0.75 !< Deepest ML as fractional depth of ocean

◆ hmlmin

real, parameter rossby_front_2d_initialization::hmlmin = 0.25
private

Shallowest ML as fractional depth of ocean.

Definition at line 35 of file Rossby_front_2d_initialization.F90.

Referenced by dtdy(), and hml().

35 real, parameter :: hmlmin = 0.25 !< Shallowest ML as fractional depth of ocean

◆ mdl

character(len=40) rossby_front_2d_initialization::mdl = "Rossby_front_2d_initialization"
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().

25 character(len=40) :: mdl = "Rossby_front_2d_initialization" !< This module's name.