MOM6
dumbbell_initialization Module Reference

Detailed Description

Configures the model for the idealized dumbbell test case.

Functions/Subroutines

subroutine, public dumbbell_initialize_topography (D, G, param_file, max_depth)
 Initialization of topography. More...
 
subroutine, public dumbbell_initialize_thickness (h, G, GV, US, param_file, just_read_params)
 Initializes the layer thicknesses to be uniform in the dumbbell test case. More...
 
subroutine, public dumbbell_initialize_temperature_salinity (T, S, h, G, GV, param_file, eqn_of_state, just_read_params)
 Initial values for temperature and salinity for the dumbbell test case. More...
 
subroutine, public dumbbell_initialize_sponges (G, GV, US, tv, param_file, use_ALE, CSp, ACSp)
 Initialize the restoring sponges for the dumbbell test case. More...
 

Variables

character(len=40) mdl = "dumbbell_initialization"
 This module's name.
 

Function/Subroutine Documentation

◆ dumbbell_initialize_sponges()

subroutine, public dumbbell_initialization::dumbbell_initialize_sponges ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(thermo_var_ptrs), intent(in)  tv,
type(param_file_type), intent(in)  param_file,
logical, intent(in)  use_ALE,
type(sponge_cs), pointer  CSp,
type(ale_sponge_cs), pointer  ACSp 
)

Initialize the restoring sponges for the dumbbell test case.

Parameters
[in]gHorizontal grid control structure
[in]gvVertical grid control structure
[in]usA dimensional unit scaling type
[in]tvThermodynamic variables
[in]param_fileParameter file structure
[in]use_aleALE flag
cspLayered sponge control structure pointer
acspALE sponge control structure pointer

Definition at line 288 of file dumbbell_initialization.F90.

288  type(ocean_grid_type), intent(in) :: g !< Horizontal grid control structure
289  type(verticalgrid_type), intent(in) :: gv !< Vertical grid control structure
290  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
291  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables
292  type(param_file_type), intent(in) :: param_file !< Parameter file structure
293  logical, intent(in) :: use_ale !< ALE flag
294  type(sponge_cs), pointer :: csp !< Layered sponge control structure pointer
295  type(ale_sponge_cs), pointer :: acsp !< ALE sponge control structure pointer
296 
297  real :: sponge_time_scale ! The damping time scale [T ~> s]
298 
299  real, dimension(SZI_(G),SZJ_(G)) :: idamp ! inverse damping timescale [T-1 ~> s-1]
300  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h, t, s ! sponge thicknesses, temp and salt
301  real, dimension(SZK_(GV)+1) :: e0, eta1d ! interface positions for ALE sponge
302 
303  integer :: i, j, k, nz
304  real :: x, zi, zmid, dist, min_thickness, dblen
305  real :: mld, s_ref, s_range, s_dense, t_ref, sill_height
306  logical :: dbrotate ! If true, rotate the domain.
307 
308  call get_param(param_file, mdl,"DUMBBELL_LEN",dblen, &
309  'Lateral Length scale for dumbbell ',&
310  units='k', default=600., do_not_log=.true.)
311  call get_param(param_file, mdl, "DUMBBELL_ROTATION", dbrotate, &
312  'Logical for rotation of dumbbell domain.',&
313  units='nondim', default=.false., do_not_log=.true.)
314 
315  if (g%x_axis_units == 'm') then
316  dblen=dblen*1.e3
317  endif
318 
319  nz = gv%ke
320 
321  call get_param(param_file, mdl, "DUMBBELL_SPONGE_TIME_SCALE", sponge_time_scale, &
322  "The time scale in the reservoir for restoring. If zero, the sponge is disabled.", &
323  units="s", default=0., scale=us%s_to_T)
324  call get_param(param_file, mdl, "DUMBBELL_SREF", s_ref, do_not_log=.true.)
325  call get_param(param_file, mdl, "DUMBBELL_S_RANGE", s_range, do_not_log=.true.)
326  call get_param(param_file, mdl,"MIN_THICKNESS", min_thickness, &
327  'Minimum thickness for layer',&
328  units='m', default=1.0e-3, do_not_log=.true., scale=us%m_to_Z)
329 
330  ! no active sponges
331  if (sponge_time_scale <= 0.) return
332 
333  ! everywhere is initially unsponged
334  idamp(:,:) = 0.0
335 
336  do j = g%jsc, g%jec
337  do i = g%isc,g%iec
338  if (g%mask2dT(i,j) > 0.) then
339  ! nondimensional x position
340  if (dbrotate) then
341  ! This is really y in the rotated case
342  x = ( g%geoLatT(i,j) ) / dblen
343  else
344  x = ( g%geoLonT(i,j) ) / dblen
345  endif
346  if (x > 0.25 .or. x < -0.25) then
347  ! scale restoring by depth into sponge
348  idamp(i,j) = 1. / sponge_time_scale
349  endif
350  endif
351  enddo
352  enddo
353 
354  if (use_ale) then
355  ! construct a uniform grid for the sponge
356  do j=g%jsc,g%jec ; do i=g%isc,g%iec
357  eta1d(nz+1) = -g%bathyT(i,j)
358  do k=nz,1,-1
359  eta1d(k) = -g%max_depth * real(k-1) / real(nz)
360  if (eta1d(k) < (eta1d(k+1) + min_thickness)) then
361  eta1d(k) = eta1d(k+1) + min_thickness
362  h(i,j,k) = gv%Z_to_H * min_thickness
363  else
364  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
365  endif
366  enddo
367  enddo ; enddo
368 
369  call initialize_ale_sponge(idamp, g, param_file, acsp, h, nz)
370 
371  ! construct temperature and salinity for the sponge
372  ! start with initial condition
373  s(:,:,:) = 0.0
374 
375  do j=g%jsc,g%jec ; do i=g%isc,g%iec
376  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
377  if (dbrotate) then
378  ! This is really y in the rotated case
379  x = ( g%geoLatT(i,j) ) / dblen
380  else
381  x = ( g%geoLonT(i,j) ) / dblen
382  endif
383  if (x>=0.25 ) then
384  do k=1,nz
385  s(i,j,k)=s_ref + 0.5*s_range
386  enddo
387  endif
388  if (x<=-0.25 ) then
389  do k=1,nz
390  s(i,j,k)=s_ref - 0.5*s_range
391  enddo
392  endif
393  enddo ; enddo
394  endif
395 
396  if (associated(tv%S)) call set_up_ale_sponge_field(s, g, tv%S, acsp)
397 

◆ dumbbell_initialize_temperature_salinity()

subroutine, public dumbbell_initialization::dumbbell_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(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 
)

Initial values for temperature and salinity for the dumbbell test case.

Parameters
[in]gOcean grid structure
[in]gvVertical grid structure
[out]tPotential temperature [degC]
[out]sSalinity [ppt]
[in]hLayer thickness [H ~> m or kg m-2]
[in]param_fileParameter file structure
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 212 of file dumbbell_initialization.F90.

212  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
213  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
214  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: t !< Potential temperature [degC]
215  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: s !< Salinity [ppt]
216  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
217  type(param_file_type), intent(in) :: param_file !< Parameter file structure
218  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
219  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
220  !! only read parameters without changing h.
221 
222  ! Local variables
223  integer :: i, j, k, is, ie, js, je, nz, k_light
224  real :: xi0, xi1, dxi, r, s_surf, t_surf, s_range, t_range
225  real :: x, y, dblen
226  real :: t_ref, t_light, t_dense, s_ref, s_light, s_dense, a1, frac_dense, k_frac, res_rat
227  logical :: just_read ! If true, just read parameters but set nothing.
228  logical :: dbrotate ! If true, rotate the domain.
229  character(len=20) :: verticalcoordinate, density_profile
230 
231  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
232 
233  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
234 
235  t_surf = 20.0
236 
237  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
238  default=default_coordinate_mode, do_not_log=just_read)
239  call get_param(param_file, mdl,"INITIAL_DENSITY_PROFILE", density_profile, &
240  'Initial profile shape. Valid values are "linear", "parabolic" '// &
241  'and "exponential".', default='linear', do_not_log=just_read)
242  call get_param(param_file, mdl,"DUMBBELL_SREF", s_surf, &
243  'DUMBBELL REFERENCE SALINITY', units='1e-3', default=34., do_not_log=just_read)
244  call get_param(param_file, mdl,"DUMBBELL_S_RANGE", s_range, &
245  'DUMBBELL salinity range (right-left)', units='1e-3', &
246  default=2., do_not_log=just_read)
247  call get_param(param_file, mdl,"DUMBBELL_LEN",dblen, &
248  'Lateral Length scale for dumbbell ',&
249  units='k', default=600., do_not_log=just_read)
250  call get_param(param_file, mdl, "DUMBBELL_ROTATION", dbrotate, &
251  'Logical for rotation of dumbbell domain.',&
252  units='nondim', default=.false., do_not_log=just_read)
253 
254  if (g%x_axis_units == 'm') then
255  dblen=dblen*1.e3
256  endif
257 
258  do j=g%jsc,g%jec
259  do i=g%isc,g%iec
260  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
261  if (dbrotate) then
262  ! This is really y in the rotated case
263  x = ( g%geoLatT(i,j) ) / dblen
264  else
265  x = ( g%geoLonT(i,j) ) / dblen
266  endif
267  do k=1,nz
268  t(i,j,k)=t_surf
269  enddo
270  if (x>=0. ) then
271  do k=1,nz
272  s(i,j,k)=s_surf + 0.5*s_range
273  enddo
274  endif
275  if (x<0. ) then
276  do k=1,nz
277  s(i,j,k)=s_surf - 0.5*s_range
278  enddo
279  endif
280 
281  enddo
282  enddo
283 

◆ dumbbell_initialize_thickness()

subroutine, public dumbbell_initialization::dumbbell_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(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
logical, intent(in), optional  just_read_params 
)

Initializes the layer thicknesses to be uniform in the dumbbell test case.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[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 94 of file dumbbell_initialization.F90.

94  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
95  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
96  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
97  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
98  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
99  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
100  !! to parse for model parameter values.
101  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
102  !! only read parameters without changing h.
103 
104  real :: e0(szk_(g)+1) ! The resting interface heights [Z ~> m], usually
105  ! negative because it is positive upward.
106  real :: eta1d(szk_(g)+1)! Interface height relative to the sea surface
107  ! positive upward [Z ~> m].
108  real :: min_thickness ! The minimum layer thicknesses [Z ~> m].
109  real :: s_surf, s_range, s_ref, s_light, s_dense ! Various salinities [ppt].
110  real :: eta_ic_quanta ! The granularity of quantization of intial interface heights [Z-1 ~> m-1].
111  ! This include declares and sets the variable "version".
112 # include "version_variable.h"
113  character(len=20) :: verticalcoordinate
114  logical :: just_read ! If true, just read parameters but set nothing.
115  integer :: i, j, k, is, ie, js, je, nz
116 
117  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
118 
119  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
120 
121  if (.not.just_read) &
122  call mom_mesg("MOM_initialization.F90, initialize_thickness_uniform: setting thickness")
123 
124  if (.not.just_read) call log_version(param_file, mdl, version, "")
125  call get_param(param_file, mdl,"MIN_THICKNESS", min_thickness, &
126  'Minimum thickness for layer',&
127  units='m', default=1.0e-3, do_not_log=just_read, scale=us%m_to_Z)
128  call get_param(param_file, mdl,"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
129  default=default_coordinate_mode, do_not_log=just_read)
130 
131  ! WARNING: this routine specifies the interface heights so that the last layer
132  ! is vanished, even at maximum depth. In order to have a uniform
133  ! layer distribution, use this line of code within the loop:
134  ! e0(k) = -G%max_depth * real(k-1) / real(nz)
135  ! To obtain a thickness distribution where the last layer is
136  ! vanished and the other thicknesses uniformly distributed, use:
137  ! e0(k) = -G%max_depth * real(k-1) / real(nz-1)
138  !do k=1,nz+1
139  ! e0(k) = -G%max_depth * real(k-1) / real(nz)
140  !enddo
141 
142  select case ( coordinatemode(verticalcoordinate) )
143 
144  case ( regridding_layer, regridding_rho ) ! Initial thicknesses for isopycnal coordinates
145  call get_param(param_file, mdl,"INITIAL_SSS", s_surf, default=34., do_not_log=.true.)
146  call get_param(param_file, mdl,"INITIAL_S_RANGE", s_range, default=2., do_not_log=.true.)
147  call get_param(param_file, mdl, "S_REF", s_ref, default=35.0, do_not_log=.true.)
148  call get_param(param_file, mdl, "TS_RANGE_S_LIGHT", s_light, default = s_ref, do_not_log=.true.)
149  call get_param(param_file, mdl, "TS_RANGE_S_DENSE", s_dense, default = s_ref, do_not_log=.true.)
150  call get_param(param_file, mdl, "INTERFACE_IC_QUANTA", eta_ic_quanta, &
151  "The granularity of initial interface height values "//&
152  "per meter, to avoid sensivity to order-of-arithmetic changes.", &
153  default=2048.0, units="m-1", scale=us%Z_to_m, do_not_log=just_read)
154  if (just_read) return ! All run-time parameters have been read, so return.
155 
156  do k=1,nz+1
157  ! Salinity of layer k is S_light + (k-1)/(nz-1) * (S_dense - S_light)
158  ! Salinity of interface K is S_light + (K-3/2)/(nz-1) * (S_dense - S_light)
159  ! Salinity at depth z should be S(z) = S_surf - S_range * z/max_depth
160  ! Equating: S_surf - S_range * z/max_depth = S_light + (K-3/2)/(nz-1) * (S_dense - S_light)
161  ! Equating: - S_range * z/max_depth = S_light - S_surf + (K-3/2)/(nz-1) * (S_dense - S_light)
162  ! Equating: z/max_depth = - ( S_light - S_surf + (K-3/2)/(nz-1) * (S_dense - S_light) ) / S_range
163  e0(k) = - g%max_depth * ( ( s_light - s_surf ) + ( s_dense - s_light ) * &
164  ( (real(k)-1.5) / real(nz-1) ) ) / s_range
165  ! Force round numbers ... the above expression has irrational factors ...
166  if (eta_ic_quanta > 0.0) &
167  e0(k) = nint(eta_ic_quanta*e0(k)) / eta_ic_quanta
168  e0(k) = min(real(1-k)*gv%angstrom_z, e0(k)) ! Bound by surface
169  e0(k) = max(-g%max_depth, e0(k)) ! Bound by bottom
170  enddo
171  do j=js,je ; do i=is,ie
172  eta1d(nz+1) = -g%bathyT(i,j)
173  do k=nz,1,-1
174  eta1d(k) = e0(k)
175  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z)) then
176  eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
177  h(i,j,k) = gv%Angstrom_H
178  else
179  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
180  endif
181  enddo
182  enddo ; enddo
183 
184  case ( regridding_zstar ) ! Initial thicknesses for z coordinates
185  if (just_read) return ! All run-time parameters have been read, so return.
186  do j=js,je ; do i=is,ie
187  eta1d(nz+1) = -g%bathyT(i,j)
188  do k=nz,1,-1
189  eta1d(k) = -g%max_depth * real(k-1) / real(nz)
190  if (eta1d(k) < (eta1d(k+1) + min_thickness)) then
191  eta1d(k) = eta1d(k+1) + min_thickness
192  h(i,j,k) = gv%Z_to_H * min_thickness
193  else
194  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
195  endif
196  enddo
197  enddo ; enddo
198 
199  case ( regridding_sigma ) ! Initial thicknesses for sigma coordinates
200  if (just_read) return ! All run-time parameters have been read, so return.
201  do j=js,je ; do i=is,ie
202  h(i,j,:) = gv%Z_to_H * g%bathyT(i,j) / dfloat(nz)
203  enddo ; enddo
204 
205 end select
206 

◆ dumbbell_initialize_topography()

subroutine, public dumbbell_initialization::dumbbell_initialize_topography ( real, dimension(g%isd:g%ied,g%jsd:g%jed), intent(out)  D,
type(dyn_horgrid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
real, intent(in)  max_depth 
)

Initialization of topography.

Parameters
[in]gThe dynamic horizontal grid type
[out]dOcean bottom depth in the units of depth_max
[in]param_fileParameter file structure
[in]max_depthMaximum ocean depth in arbitrary units

Definition at line 43 of file dumbbell_initialization.F90.

43  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
44  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
45  intent(out) :: d !< Ocean bottom depth in the units of depth_max
46  type(param_file_type), intent(in) :: param_file !< Parameter file structure
47  real, intent(in) :: max_depth !< Maximum ocean depth in arbitrary units
48 
49  ! Local variables
50  integer :: i, j
51  real :: x, y, delta, dblen, dbfrac
52  logical :: dbrotate
53 
54  call get_param(param_file, mdl,"DUMBBELL_LEN",dblen, &
55  'Lateral Length scale for dumbbell.',&
56  units='k', default=600., do_not_log=.false.)
57  call get_param(param_file, mdl,"DUMBBELL_FRACTION",dbfrac, &
58  'Meridional fraction for narrow part of dumbbell.',&
59  units='nondim', default=0.5, do_not_log=.false.)
60  call get_param(param_file, mdl, "DUMBBELL_ROTATION", dbrotate, &
61  'Logical for rotation of dumbbell domain.',&
62  units='nondim', default=.false., do_not_log=.false.)
63 
64  if (g%x_axis_units == 'm') then
65  dblen=dblen*1.e3
66  endif
67 
68  if (dbrotate) then
69  do j=g%jsc,g%jec ; do i=g%isc,g%iec
70  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
71  x = ( g%geoLonT(i,j) ) / g%len_lon
72  y = ( g%geoLatT(i,j) ) / dblen
73  d(i,j) = g%max_depth
74  if ((y>=-0.25 .and. y<=0.25) .and. (x <= -0.5*dbfrac .or. x >= 0.5*dbfrac)) then
75  d(i,j) = 0.0
76  endif
77  enddo ; enddo
78  else
79  do j=g%jsc,g%jec ; do i=g%isc,g%iec
80  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
81  x = ( g%geoLonT(i,j) ) / dblen
82  y = ( g%geoLatT(i,j) ) / g%len_lat
83  d(i,j) = g%max_depth
84  if ((x>=-0.25 .and. x<=0.25) .and. (y <= -0.5*dbfrac .or. y >= 0.5*dbfrac)) then
85  d(i,j) = 0.0
86  endif
87  enddo ; enddo
88  endif
89