MOM6
sloshing_initialization Module Reference

Detailed Description

Initialization for the "sloshing" internal waves configuration.

The module configures the model for the non-rotating sloshing test case.

Functions/Subroutines

subroutine, public sloshing_initialize_topography (D, G, param_file, max_depth)
 Initialization of topography. More...
 
subroutine, public sloshing_initialize_thickness (h, G, GV, US, param_file, just_read_params)
 Initialization of thicknesses This routine is called when THICKNESS_CONFIG is set to 'sloshing'. More...
 
subroutine, public sloshing_initialize_temperature_salinity (T, S, h, G, GV, param_file, eqn_of_state, just_read_params)
 Initialization of temperature and salinity. More...
 

Function/Subroutine Documentation

◆ sloshing_initialize_temperature_salinity()

subroutine, public sloshing_initialization::sloshing_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.

This subroutine initializes linear profiles for T and S according to reference surface layer salinity and temperature and a specified range. Note that the linear distribution is set up with respect to the layer number, not the physical position).

Parameters
[in]gOcean grid structure.
[in]gvThe ocean's vertical grid structure.
[out]tPotential temperature [degC].
[out]sSalinity [ppt].
[in]hLayer thickness [H ~> m or kg m-2].
[in]param_fileA structure indicating the open file to parse for model parameter values.
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 181 of file sloshing_initialization.F90.

181  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure.
182  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
183  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: T !< Potential temperature [degC].
184  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: S !< Salinity [ppt].
185  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
186  type(param_file_type), intent(in) :: param_file !< A structure indicating the
187  !! open file to parse for model
188  !! parameter values.
189  type(EOS_type), pointer :: eqn_of_state !< Equation of state structure.
190  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
191  !! only read parameters without changing h.
192 
193  integer :: i, j, k, is, ie, js, je, nz
194  real :: delta_S, delta_T
195  real :: S_ref, T_ref; ! Reference salinity and temerature within
196  ! surface layer
197  real :: S_range, T_range; ! Range of salinities and temperatures over the
198  ! vertical
199  integer :: kdelta
200  real :: deltah
201  real :: xi0, xi1
202  logical :: just_read ! If true, just read parameters but set nothing.
203  character(len=40) :: mdl = "initialize_temp_salt_linear" ! This subroutine's
204  ! name.
205 
206  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
207 
208  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
209 
210  call get_param(param_file, mdl,"S_REF",s_ref,'Reference value for salinity', &
211  units='1e-3', fail_if_missing=.not.just_read, do_not_log=just_read)
212  call get_param(param_file, mdl,"T_REF",t_ref,'Refernce value for temperature', &
213  units='C', fail_if_missing=.not.just_read, do_not_log=just_read)
214 
215  ! The default is to assume an increase by 2 for the salinity and a uniform
216  ! temperature
217  call get_param(param_file, mdl,"S_RANGE",s_range,'Initial salinity range.', &
218  units='1e-3', default=2.0, do_not_log=just_read)
219  call get_param(param_file, mdl,"T_RANGE",t_range,'Initial temperature range', &
220  units='C', default=0.0, do_not_log=just_read)
221 
222  if (just_read) return ! All run-time parameters have been read, so return.
223 
224  ! Prescribe salinity
225  !delta_S = S_range / ( G%ke - 1.0 )
226 
227  !S(:,:,1) = S_ref
228  !do k = 2,G%ke
229  ! S(:,:,k) = S(:,:,k-1) + delta_S
230  !enddo
231 
232  deltah = g%max_depth / nz
233  do j=js,je ; do i=is,ie
234  xi0 = 0.0
235  do k = 1,nz
236  xi1 = xi0 + deltah / g%max_depth ! = xi0 + 1.0 / real(nz)
237  s(i,j,k) = 34.0 + 0.5 * s_range * (xi0 + xi1)
238  xi0 = xi1
239  enddo
240  enddo ; enddo
241 
242  ! Prescribe temperature
243  delta_t = t_range / ( g%ke - 1.0 )
244 
245  t(:,:,1) = t_ref
246  do k = 2,g%ke
247  t(:,:,k) = t(:,:,k-1) + delta_t
248  enddo
249  kdelta = 2
250  t(:,:,g%ke/2 - (kdelta-1):g%ke/2 + kdelta) = 1.0
251 

Referenced by mom_state_initialization::mom_initialize_state().

Here is the caller graph for this function:

◆ sloshing_initialize_thickness()

subroutine, public sloshing_initialization::sloshing_initialize_thickness ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), 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 
)

Initialization of thicknesses This routine is called when THICKNESS_CONFIG is set to 'sloshing'.

This routine initializes layer positions to set off a sloshing motion in the zonal direction in a rectangular basin. All layers have initially the same thickness but all interfaces (except bottom and sea surface) are displaced according to a half-period cosine, with maximum value on the left and minimum value on the right. This sets off a regular sloshing motion.

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 57 of file sloshing_initialization.F90.

57  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
58  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
59  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
60  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
61  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
62  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
63  !! to parse for model parameter values.
64  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
65  !! only read parameters without changing h.
66 
67  real :: displ(SZK_(G)+1) ! The interface displacement in depth units.
68  real :: z_unif(SZK_(G)+1) ! Fractional uniform interface heights [nondim].
69  real :: z_inter(SZK_(G)+1) ! Interface heights, in depth units.
70  real :: a0 ! The displacement amplitude in depth units.
71  real :: weight_z ! A (misused?) depth-space weighting, in inconsistent units.
72  real :: x1, y1, x2, y2 ! Dimensonless parameters.
73  real :: x, t ! Dimensionless depth coordinates?
74  logical :: use_IC_bug ! If true, set the initial conditions retaining an old bug.
75  logical :: just_read ! If true, just read parameters but set nothing.
76  ! This include declares and sets the variable "version".
77 # include "version_variable.h"
78  character(len=40) :: mdl = "sloshing_initialization" !< This module's name.
79 
80  integer :: i, j, k, is, ie, js, je, nx, nz
81 
82  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
83 
84  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
85  if (.not.just_read) call log_version(param_file, mdl, version, "")
86  call get_param(param_file, mdl, "SLOSHING_IC_AMPLITUDE", a0, &
87  "Initial amplitude of sloshing internal interface height "//&
88  "displacements it the sloshing test case.", &
89  units='m', default=75.0, scale=us%m_to_Z, do_not_log=just_read)
90  call get_param(param_file, mdl, "SLOSHING_IC_BUG", use_ic_bug, &
91  "If true, use code with a bug to set the sloshing initial conditions.", &
92  default=.true., do_not_log=just_read)
93 
94  if (just_read) return ! All run-time parameters have been read, so return.
95 
96  ! Define thicknesses
97  do j=g%jsc,g%jec ; do i=g%isc,g%iec
98 
99  ! Define uniform interfaces
100  do k = 0,nz
101  z_unif(k+1) = -real(k)/real(nz)
102  enddo
103 
104  ! 1. Define stratification
105  do k = 1,nz+1
106 
107  ! Thin pycnocline in the middle
108  !z_inter(k) = (2.0**(n-1)) * (z_unif(k) + 0.5)**n - 0.5
109 
110  ! Thin pycnocline in the middle (piecewise linear profile)
111  x1 = 0.30; y1 = 0.48; x2 = 0.70; y2 = 0.52
112 
113  x = -z_unif(k)
114 
115  if ( x <= x1 ) then
116  t = y1*x/x1
117  elseif ( (x > x1 ) .and. ( x < x2 )) then
118  t = y1 + (y2-y1) * (x-x1) / (x2-x1)
119  else
120  t = y2 + (1.0-y2) * (x-x2) / (1.0-x2)
121  endif
122 
123  t = - z_unif(k)
124 
125  z_inter(k) = -t * g%max_depth
126 
127  enddo
128 
129  ! 2. Define displacement
130  ! a0 is set via get_param; by default a0 is a 75m Displacement amplitude in depth units.
131  do k = 1,nz+1
132 
133  weight_z = - 4.0 * ( z_unif(k) + 0.5 )**2 + 1.0
134 
135  x = g%geoLonT(i,j) / g%len_lon
136  if (use_ic_bug) then
137  displ(k) = a0 * cos(acos(-1.0)*x) + weight_z * us%m_to_Z
138  else
139  displ(k) = a0 * cos(acos(-1.0)*x) * weight_z
140  endif
141 
142  if ( k == 1 ) then
143  displ(k) = 0.0
144  endif
145 
146  if ( k == nz+1 ) then
147  displ(k) = 0.0
148  endif
149 
150  z_inter(k) = z_inter(k) + displ(k)
151 
152  enddo
153 
154  ! 3. The last interface must coincide with the seabed
155  z_inter(nz+1) = -g%bathyT(i,j)
156  ! Modify interface heights to make sure all thicknesses are strictly positive
157  do k = nz,1,-1
158  if ( z_inter(k) < (z_inter(k+1) + gv%Angstrom_Z) ) then
159  z_inter(k) = z_inter(k+1) + gv%Angstrom_Z
160  endif
161  enddo
162 
163  ! 4. Define layers
164  do k = 1,nz
165  h(i,j,k) = gv%Z_to_H * (z_inter(k) - z_inter(k+1))
166  enddo
167 
168  enddo ; enddo
169 

◆ sloshing_initialize_topography()

subroutine, public sloshing_initialization::sloshing_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 32 of file sloshing_initialization.F90.

32  type(dyn_horgrid_type), intent(in) :: G !< The dynamic horizontal grid type
33  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
34  intent(out) :: D !< Ocean bottom depth in the units of depth_max
35  type(param_file_type), intent(in) :: param_file !< Parameter file structure
36  real, intent(in) :: max_depth !< Maximum ocean depth in arbitrary units
37 
38  ! Local variables
39  integer :: i, j
40 
41  do i=g%isc,g%iec ; do j=g%jsc,g%jec
42  d(i,j) = max_depth
43  enddo ; enddo
44 

Referenced by mom_fixed_initialization::mom_initialize_topography().

Here is the caller graph for this function: