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, GV, 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, US, 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, in the same units as Gmax_depth. 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.
 
real, parameter frontfractionalwidth = 0.5
 Width of front as fraction of domain [nondim].
 
real, parameter hmlmin = 0.25
 Shallowest ML as fractional depth of ocean [nondim].
 
real, parameter hmlmax = 0.75
 Deepest ML as fractional depth of ocean [nondim].
 

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 254 of file Rossby_front_2d_initialization.F90.

254  type(ocean_grid_type), intent(in) :: G !< Grid structure
255  real, intent(in) :: dT !< Top to bottom temperature difference
256  real, intent(in) :: lat !< Latitude
257  ! Local
258  real :: PI, dHML, dHdy
259  real :: km = 1.e3 ! AXIS_UNITS = 'k' (1000 m)
260 
261  pi = 4.0 * atan(1.0)
262  dhml = 0.5 * ( hmlmax - hmlmin ) * g%max_depth
263  dhdy = dhml * ( pi / ( frontfractionalwidth * g%len_lat * km ) ) * cos( ypseudo(g, lat) )
264  dtdy = -( dt / g%max_depth ) * dhdy
265 

◆ 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, in the same units as Gmax_depth.

Parameters
[in]gGrid structure
[in]latLatitude

Definition at line 241 of file Rossby_front_2d_initialization.F90.

241  type(ocean_grid_type), intent(in) :: G !< Grid structure
242  real, intent(in) :: lat !< Latitude
243  ! Local
244  real :: dHML, HMLmean
245 
246  dhml = 0.5 * ( hmlmax - hmlmin ) * g%max_depth
247  hmlmean = 0.5 * ( hmlmin + hmlmax ) * g%max_depth
248  hml = hmlmean + dhml * sin( ypseudo(g, lat) )

◆ rossby_front_initialize_temperature_salinity()

subroutine, public rossby_front_2d_initialization::rossby_front_initialize_temperature_salinity ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(out)  T,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(out)  S,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
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
[in]gvThe ocean's vertical grid structure.
[out]tPotential temperature [degC]
[out]sSalinity [ppt]
[in]hThickness [H ~> m or kg m-2]
[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 T & S.

Definition at line 112 of file Rossby_front_2d_initialization.F90.

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

◆ 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 [H ~> m or kg m-2]
[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 40 of file Rossby_front_2d_initialization.F90.

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

◆ rossby_front_initialize_velocity()

subroutine, public rossby_front_2d_initialization::rossby_front_initialize_velocity ( real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(out)  u,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(out)  v,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
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
[in]usA dimensional unit scaling type
[out]ui-component of velocity [m s-1]
[out]vj-component of velocity [m s-1]
[in]hThickness [H ~> m or kg m-2]
[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 setting u & v.

Definition at line 164 of file Rossby_front_2d_initialization.F90.

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

◆ 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 227 of file Rossby_front_2d_initialization.F90.

227  type(ocean_grid_type), intent(in) :: G !< Grid structure
228  real, intent(in) :: lat !< Latitude
229  ! Local
230  real :: y, PI
231 
232  pi = 4.0 * atan(1.0)
233  ypseudo = ( ( lat - g%south_lat ) / g%len_lat ) - 0.5 ! -1/2 .. 1/.2
234  ypseudo = pi * max(-0.5, min(0.5, ypseudo / frontfractionalwidth))