MOM6
Rossby_front_2d_initialization.F90
1 !> Initial conditions for the 2D Rossby front test
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_mesg, mom_error, fatal, is_root_pe
8 use mom_get_input, only : directories
9 use mom_grid, only : ocean_grid_type
14 use regrid_consts, only : coordinatemode, default_coordinate_mode
15 use regrid_consts, only : regridding_layer, regridding_zstar
16 use regrid_consts, only : regridding_rho, regridding_sigma
17 
18 implicit none ; private
19 
20 #include <MOM_memory.h>
21 
22 ! Private (module-wise) parameters
23 character(len=40) :: mdl = "Rossby_front_2d_initialization" !< This module's name.
24 ! This include declares and sets the variable "version".
25 #include "version_variable.h"
26 
27 public rossby_front_initialize_thickness
28 public rossby_front_initialize_temperature_salinity
29 public rossby_front_initialize_velocity
30 
31 ! Parameters defining the initial conditions of this test case
32 real, parameter :: frontfractionalwidth = 0.5 !< Width of front as fraction of domain [nondim]
33 real, parameter :: hmlmin = 0.25 !< Shallowest ML as fractional depth of ocean [nondim]
34 real, parameter :: hmlmax = 0.75 !< Deepest ML as fractional depth of ocean [nondim]
35 
36 contains
37 
38 !> Initialization of thicknesses in 2D Rossby front test
39 subroutine rossby_front_initialize_thickness(h, G, GV, US, param_file, just_read_params)
40  type(ocean_grid_type), intent(in) :: G !< Grid structure
41  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
42  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
43  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
44  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2]
45  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
46  !! to parse for model parameter values.
47  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
48  !! only read parameters without changing h.
49 
50  integer :: i, j, k, is, ie, js, je, nz
51  real :: Tz, Dml, eta, stretch, h0
52  real :: min_thickness, T_range
53  real :: dRho_dT ! The partial derivative of density with temperature [R degC-1 ~> kg m-3 degC-1]
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, scale=us%kg_m3_to_R, 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 * gv%Z_to_H
88  enddo
89  enddo ; enddo
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 * gv%Z_to_H
99  enddo
100  enddo ; enddo
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 
108 end subroutine rossby_front_initialize_thickness
109 
110 
111 !> Initialization of temperature and salinity in the Rossby front test
112 subroutine rossby_front_initialize_temperature_salinity(T, S, h, G, GV, &
113  param_file, eqn_of_state, just_read_params)
114  type(ocean_grid_type), intent(in) :: G !< Grid structure
115  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
116  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: T !< Potential temperature [degC]
117  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: S !< Salinity [ppt]
118  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Thickness [H ~> m or kg m-2]
119  type(param_file_type), intent(in) :: param_file !< Parameter file handle
120  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
121  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
122  !! only read parameters without changing T & S.
123 
124  integer :: i, j, k, is, ie, js, je, nz
125  real :: T_ref, S_ref ! Reference salinity and temerature within surface layer
126  real :: T_range ! Range of salinities and temperatures over the vertical
127  real :: y, zc, zi, dTdz
128  logical :: just_read ! If true, just read parameters but set nothing.
129  character(len=40) :: verticalCoordinate
130  real :: PI ! 3.1415926... calculated as 4*atan(1)
131 
132  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
133 
134  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
135 
136  call get_param(param_file, mdl,"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
137  default=default_coordinate_mode, do_not_log=just_read)
138  call get_param(param_file, mdl, "S_REF", s_ref, 'Reference salinity', &
139  default=35.0, units='1e-3', do_not_log=just_read)
140  call get_param(param_file, mdl,"T_REF",t_ref,'Reference temperature',units='C',&
141  fail_if_missing=.not.just_read, do_not_log=just_read)
142  call get_param(param_file, mdl,"T_RANGE",t_range,'Initial temperature range',&
143  units='C', default=0.0, do_not_log=just_read)
144 
145  if (just_read) return ! All run-time parameters have been read, so return.
146 
147  t(:,:,:) = 0.0
148  s(:,:,:) = s_ref
149  dtdz = t_range / g%max_depth
150 
151  do j = g%jsc,g%jec ; do i = g%isc,g%iec
152  zi = 0.
153  do k = 1, nz
154  zi = zi - h(i,j,k) ! Bottom interface position
155  zc = gv%H_to_Z * (zi - 0.5*h(i,j,k)) ! Position of middle of cell
156  zc = min( zc, -hml(g, g%geoLatT(i,j)) ) ! Bound by depth of mixed layer
157  t(i,j,k) = t_ref + dtdz * zc ! Linear temperature profile
158  enddo
159  enddo ; enddo
160 
161 end subroutine rossby_front_initialize_temperature_salinity
162 
163 
164 !> Initialization of u and v in the Rossby front test
165 subroutine rossby_front_initialize_velocity(u, v, h, G, GV, US, param_file, just_read_params)
166  type(ocean_grid_type), intent(in) :: G !< Grid structure
167  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
168  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
169  intent(out) :: u !< i-component of velocity [L T-1 ~> m s-1]
170  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
171  intent(out) :: v !< j-component of velocity [L T-1 ~> m s-1]
172  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), &
173  intent(in) :: h !< Thickness [H ~> m or kg m-2]
174  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
175  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
176  !! to parse for model parameter values.
177  logical, optional, intent(in) :: just_read_params !< If present and true, this call
178  !! will only read parameters without setting u & v.
179 
180  real :: y ! Non-dimensional coordinate across channel, 0..pi
181  real :: T_range ! Range of salinities and temperatures over the vertical
182  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]
183  real :: dRho_dT ! The partial derivative of density with temperature [R degC-1 ~> kg m-3 degC-1]
184  real :: Dml, zi, zc, zm ! Depths [Z ~> m].
185  real :: f ! The local Coriolis parameter [T-1 ~> s-1]
186  real :: Ty ! The meridional temperature gradient [degC L-1 ~> degC m-1]
187  real :: hAtU ! Interpolated layer thickness [Z ~> m].
188  integer :: i, j, k, is, ie, js, je, nz
189  logical :: just_read ! If true, just read parameters but set nothing.
190  character(len=40) :: verticalCoordinate
191 
192  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
193 
194  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
195 
196  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
197  default=default_coordinate_mode, do_not_log=just_read)
198  call get_param(param_file, mdl, "T_RANGE", t_range, 'Initial temperature range', &
199  units='C', default=0.0, do_not_log=just_read)
200  call get_param(param_file, mdl, "DRHO_DT", drho_dt, default=-0.2, scale=us%kg_m3_to_R, do_not_log=.true.)
201 
202  if (just_read) return ! All run-time parameters have been read, so return.
203 
204  v(:,:,:) = 0.0
205  u(:,:,:) = 0.0
206 
207  do j = g%jsc,g%jec ; do i = g%isc-1,g%iec+1
208  f = 0.5* (g%CoriolisBu(i,j) + g%CoriolisBu(i,j-1) )
209  dudt = 0.0 ; if (abs(f) > 0.0) &
210  dudt = ( gv%g_Earth*drho_dt ) / ( f * gv%Rho0 )
211  dml = hml( g, g%geoLatT(i,j) )
212  ty = us%L_to_m*dtdy( g, t_range, g%geoLatT(i,j) )
213  zi = 0.
214  do k = 1, nz
215  hatu = 0.5*(h(i,j,k)+h(i+1,j,k)) * gv%H_to_Z
216  zi = zi - hatu ! Bottom interface position
217  zc = zi - 0.5*hatu ! Position of middle of cell
218  zm = max( zc + dml, 0. ) ! Height above bottom of mixed layer
219  u(i,j,k) = dudt * ty * zm ! Thermal wind starting at base of ML
220  enddo
221  enddo ; enddo
222 
223 end subroutine rossby_front_initialize_velocity
224 
225 !> Pseudo coordinate across domain used by Hml() and dTdy()
226 !! returns a coordinate from -PI/2 .. PI/2 squashed towards the
227 !! center of the domain.
228 real function ypseudo( G, lat )
229  type(ocean_grid_type), intent(in) :: G !< Grid structure
230  real, intent(in) :: lat !< Latitude
231  ! Local
232  real :: y, PI
233 
234  pi = 4.0 * atan(1.0)
235  ypseudo = ( ( lat - g%south_lat ) / g%len_lat ) - 0.5 ! -1/2 .. 1/.2
236  ypseudo = pi * max(-0.5, min(0.5, ypseudo / frontfractionalwidth))
237 end function ypseudo
238 
239 
240 !> Analytic prescription of mixed layer depth in 2d Rossby front test,
241 !! in the same units as G%max_depth
242 real function hml( G, lat )
243  type(ocean_grid_type), intent(in) :: G !< Grid structure
244  real, intent(in) :: lat !< Latitude
245  ! Local
246  real :: dHML, HMLmean
247 
248  dhml = 0.5 * ( hmlmax - hmlmin ) * g%max_depth
249  hmlmean = 0.5 * ( hmlmin + hmlmax ) * g%max_depth
250  hml = hmlmean + dhml * sin( ypseudo(g, lat) )
251 end function hml
252 
253 
254 !> Analytic prescription of mixed layer temperature gradient in 2d Rossby front test
255 real function dtdy( G, dT, lat )
256  type(ocean_grid_type), intent(in) :: G !< Grid structure
257  real, intent(in) :: dT !< Top to bottom temperature difference
258  real, intent(in) :: lat !< Latitude
259  ! Local
260  real :: PI, dHML, dHdy
261  real :: km = 1.e3 ! AXIS_UNITS = 'k' (1000 m)
262 
263  pi = 4.0 * atan(1.0)
264  dhml = 0.5 * ( hmlmax - hmlmin ) * g%max_depth
265  dhdy = dhml * ( pi / ( frontfractionalwidth * g%len_lat * km ) ) * cos( ypseudo(g, lat) )
266  dtdy = -( dt / g%max_depth ) * dhdy
267 
268 end function dtdy
269 
270 
271 !> \namespace rossby_front_2d_initialization
272 !!
273 !! \section section_Rossby_front_2d Description of the 2d Rossby front initial conditions
274 !!
275 !! Consistent with a linear equation of state, the system has a constant stratification
276 !! below the mixed layer, stratified in temperature only. Isotherms are flat below the
277 !! mixed layer and vertical within. Salinity is constant. The mixed layer has a half sine
278 !! form so that there are no mixed layer or temperature gradients at the side walls.
279 !!
280 !! Below the mixed layer the potential temperature, \f$\theta(z)\f$, is given by
281 !! \f[ \theta(z) = \theta_0 - \Delta \theta \left( z + h_{ML} \right) \f]
282 !! where \f$ \theta_0 \f$ and \f$ \Delta \theta \f$ are external model parameters.
283 !!
284 !! The depth of the mixed layer, \f$H_{ML}\f$ is
285 !! \f[ h_{ML}(y) = h_{min} + \left( h_{max} - h_{min} \right) \cos{\pi y/L} \f].
286 !! The temperature in mixed layer is given by the reference temperature at \f$z=h_{ML}\f$
287 !! so that
288 !! \f{eqnarray} \theta(y,z) =
289 !! \theta_0 - \Delta \theta \left( z + h_{ML} \right) & \forall & z < h_{ML}(y) T.B.D.
290 !! \f}
291 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
Reads the only Fortran name list needed to boot-strap the model.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
The MOM6 facility to parse input files for runtime parameters.
Initial conditions for the 2D Rossby front test.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Definition: MOM_EOS.F90:81
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Container for paths and parameter file names.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Contains constants for interpreting input parameters that control regridding.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108