MOM6
seamount_initialization Module Reference

Detailed Description

The module configures the model for the idealized seamount test case.

Functions/Subroutines

subroutine, public seamount_initialize_topography (D, G, param_file, max_depth)
 Initialization of topography. More...
 
subroutine, public seamount_initialize_thickness (h, G, GV, param_file, just_read_params)
 Initialization of thicknesses. This subroutine initializes the layer thicknesses to be uniform. More...
 
subroutine, public seamount_initialize_temperature_salinity (T, S, h, G, GV, param_file, eqn_of_state, just_read_params)
 Initial values for temperature and salinity. More...
 

Variables

character(len=40) mdl = "seamount_initialization"
 

Function/Subroutine Documentation

◆ seamount_initialize_temperature_salinity()

subroutine, public seamount_initialization::seamount_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.

Parameters
[in]gOcean grid structure
[in]gvVertical grid structure
[out]tPotential temperature (degC)
[out]sSalinity (ppt)
[in]hLayer thickness (m or Pa)
[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 209 of file seamount_initialization.F90.

References mdl, mom_error_handler::mom_error(), regrid_consts::regridding_rho, and regrid_consts::regridding_sigma.

Referenced by mom_state_initialization::mom_initialize_state().

209  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
210  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
211  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: t !< Potential temperature (degC)
212  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: s !< Salinity (ppt)
213  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Layer thickness (m or Pa)
214  type(param_file_type), intent(in) :: param_file !< Parameter file structure
215  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
216  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
217  !! only read parameters without changing h.
218 
219  ! Local variables
220  integer :: i, j, k, is, ie, js, je, nz, k_light
221  real :: xi0, xi1, dxi, r, s_surf, t_surf, s_range, t_range
222  real :: t_ref, t_light, t_dense, s_ref, s_light, s_dense, a1, frac_dense, k_frac, res_rat
223  logical :: just_read ! If true, just read parameters but set nothing.
224  character(len=20) :: verticalcoordinate, density_profile
225 
226  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
227 
228  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
229 
230  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
231  default=default_coordinate_mode, do_not_log=just_read)
232  call get_param(param_file, mdl,"INITIAL_DENSITY_PROFILE", density_profile, &
233  'Initial profile shape. Valid values are "linear", "parabolic"\n'// &
234  'and "exponential".', default='linear', do_not_log=just_read)
235  call get_param(param_file, mdl,"INITIAL_SSS", s_surf, &
236  'Initial surface salinity', units='1e-3', default=34., do_not_log=just_read)
237  call get_param(param_file, mdl,"INITIAL_SST", t_surf, &
238  'Initial surface temperature', units='C', default=0., do_not_log=just_read)
239  call get_param(param_file, mdl,"INITIAL_S_RANGE", s_range, &
240  'Initial salinity range (bottom - surface)', units='1e-3', &
241  default=2., do_not_log=just_read)
242  call get_param(param_file, mdl,"INITIAL_T_RANGE", t_range, &
243  'Initial temperature range (bottom - surface)', units='C', &
244  default=0., do_not_log=just_read)
245 
246  select case ( coordinatemode(verticalcoordinate) )
247  case ( regridding_layer ) ! Initial thicknesses for layer isopycnal coordinates
248  ! These parameters are used in MOM_fixed_initialization.F90 when CONFIG_COORD="ts_range"
249  call get_param(param_file, mdl, "T_REF", t_ref, default=10.0, do_not_log=.true.)
250  call get_param(param_file, mdl, "TS_RANGE_T_LIGHT", t_light, default=t_ref, do_not_log=.true.)
251  call get_param(param_file, mdl, "TS_RANGE_T_DENSE", t_dense, default=t_ref, do_not_log=.true.)
252  call get_param(param_file, mdl, "S_REF", s_ref, default=35.0, do_not_log=.true.)
253  call get_param(param_file, mdl, "TS_RANGE_S_LIGHT", s_light, default = s_ref, do_not_log=.true.)
254  call get_param(param_file, mdl, "TS_RANGE_S_DENSE", s_dense, default = s_ref, do_not_log=.true.)
255  call get_param(param_file, mdl, "TS_RANGE_RESOLN_RATIO", res_rat, default=1.0, do_not_log=.true.)
256  if (just_read) return ! All run-time parameters have been read, so return.
257 
258  ! Emulate the T,S used in the "ts_range" coordinate configuration code
259  k_light = gv%nk_rho_varies + 1
260  do j=js,je ; do i=is,ie
261  t(i,j,k_light) = t_light ; s(i,j,k_light) = s_light
262  enddo ; enddo
263  a1 = 2.0 * res_rat / (1.0 + res_rat)
264  do k=k_light+1,nz
265  k_frac = real(k-k_light)/real(nz-k_light)
266  frac_dense = a1 * k_frac + (1.0 - a1) * k_frac**2
267  do j=js,je ; do i=is,ie
268  t(i,j,k) = frac_dense * (t_dense - t_light) + t_light
269  s(i,j,k) = frac_dense * (s_dense - s_light) + s_light
270  enddo ; enddo
271  enddo
272  case ( regridding_sigma, regridding_zstar, regridding_rho ) ! All other coordinate use FV initialization
273  if (just_read) return ! All run-time parameters have been read, so return.
274  do j=js,je ; do i=is,ie
275  xi0 = 0.0
276  do k = 1,nz
277  xi1 = xi0 + h(i,j,k) / g%max_depth
278  select case ( trim(density_profile) )
279  case ('linear')
280  !S(i,j,k) = S_surf + S_range * 0.5 * (xi0 + xi1)
281  s(i,j,k) = s_surf + ( 0.5 * s_range ) * (xi0 + xi1) ! Coded this way to reproduce old hard-coded answers
282  t(i,j,k) = t_surf + t_range * 0.5 * (xi0 + xi1)
283  case ('parabolic')
284  s(i,j,k) = s_surf + s_range * (2.0 / 3.0) * (xi1**3 - xi0**3) / (xi1 - xi0)
285  t(i,j,k) = t_surf + t_range * (2.0 / 3.0) * (xi1**3 - xi0**3) / (xi1 - xi0)
286  case ('exponential')
287  s(i,j,k) = s_surf + s_range * (exp(xi1/r)-exp(xi0/r)) / (xi1 - xi0)
288  t(i,j,k) = t_surf + t_range * (exp(xi1/r)-exp(xi0/r)) / (xi1 - xi0)
289  case default
290  call mom_error(fatal, 'Unknown value for "INITIAL_DENSITY_PROFILE"')
291  end select
292  xi0 = xi1
293  enddo
294  enddo ; enddo
295  end select
296 
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:

◆ seamount_initialize_thickness()

subroutine, public seamount_initialization::seamount_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. This subroutine initializes the layer thicknesses to be uniform.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical 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 101 of file seamount_initialization.F90.

References mdl, mom_error_handler::mom_mesg(), regrid_consts::regridding_rho, and regrid_consts::regridding_sigma.

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

◆ seamount_initialize_topography()

subroutine, public seamount_initialization::seamount_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 m
[in]param_fileParameter file structure
[in]max_depthMaximum depth of model in m

Definition at line 60 of file seamount_initialization.F90.

References mdl.

Referenced by mom_fixed_initialization::mom_initialize_topography().

60  ! Arguments
61  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
62  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
63  intent(out) :: d !< Ocean bottom depth in m
64  type(param_file_type), intent(in) :: param_file !< Parameter file structure
65  real, intent(in) :: max_depth !< Maximum depth of model in m
66 
67  ! Local variables
68  integer :: i, j
69  real :: x, y, delta, lx, rlx, ly, rly
70 
71  call get_param(param_file, mdl,"SEAMOUNT_DELTA",delta, &
72  "Non-dimensional height of seamount.", &
73  units="non-dim", default=0.5)
74  call get_param(param_file, mdl,"SEAMOUNT_X_LENGTH_SCALE",lx, &
75  "Length scale of seamount in x-direction.\n"//&
76  "Set to zero make topography uniform in the x-direction.", &
77  units="Same as x,y", default=20.)
78  call get_param(param_file, mdl,"SEAMOUNT_Y_LENGTH_SCALE",ly, &
79  "Length scale of seamount in y-direction.\n"//&
80  "Set to zero make topography uniform in the y-direction.", &
81  units="Same as x,y", default=0.)
82 
83  lx = lx / g%len_lon
84  ly = ly / g%len_lat
85  rlx = 0. ; if (lx>0.) rlx = 1. / lx
86  rly = 0. ; if (ly>0.) rly = 1. / ly
87  do i=g%isc,g%iec
88  do j=g%jsc,g%jec
89  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
90  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon - 0.5
91  y = ( g%geoLatT(i,j) - g%south_lat ) / g%len_lat - 0.5
92  d(i,j) = g%max_depth * ( 1.0 - delta * exp(-(rlx*x)**2 -(rly*y)**2) )
93  enddo
94  enddo
95 
Here is the caller graph for this function:

Variable Documentation

◆ mdl

character(len=40) seamount_initialization::mdl = "seamount_initialization"
private

Definition at line 44 of file seamount_initialization.F90.

Referenced by seamount_initialize_temperature_salinity(), seamount_initialize_thickness(), and seamount_initialize_topography().

44 character(len=40) :: mdl = "seamount_initialization" ! This module's name.