MOM6
DOME2d_initialization.F90
Go to the documentation of this file.
2 
5 use mom_error_handler, only : mom_mesg, mom_error, fatal
7 use mom_get_input, only : directories
8 use mom_grid, only : ocean_grid_type
9 use mom_io, only : close_file, fieldtype, file_exists
10 use mom_io, only : open_file, read_data, read_axis_data, single_file
11 use mom_io, only : write_field, slasher, vardesc
19 
20 implicit none ; private
21 
22 #include <MOM_memory.h>
23 
24 ! Public functions
29 
30 character(len=40) :: mdl = "DOME2D_initialization" !< This module's name.
31 
32 contains
33 
34 !> Initialize topography with a shelf and slope in a 2D domain
35 subroutine dome2d_initialize_topography ( D, G, param_file, max_depth )
36  ! Arguments
37  type(dyn_horgrid_type), intent(in) :: G !< The dynamic horizontal grid type
38  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
39  intent(out) :: D !< Ocean bottom depth in m
40  type(param_file_type), intent(in) :: param_file !< Parameter file structure
41  real, intent(in) :: max_depth !< Maximum depth of model in m
42  ! Local variables
43  integer :: i, j
44  real :: x, bay_depth, l1, l2
45  real :: dome2d_width_bay, dome2d_width_bottom, dome2d_depth_bay
46 ! This include declares and sets the variable "version".
47 #include "version_variable.h"
48 
49  call log_version(param_file, mdl, version, "")
50  call get_param(param_file, mdl, "DOME2D_SHELF_WIDTH", dome2d_width_bay, &
51  'Width of shelf, as fraction of domain, in 2d DOME configuration.', &
52  units='nondim',default=0.1)
53  call get_param(param_file, mdl, "DOME2D_BASIN_WIDTH", dome2d_width_bottom, &
54  'Width of deep ocean basin, as fraction of domain, in 2d DOME configuration.', &
55  units='nondim',default=0.3)
56  call get_param(param_file, mdl, "DOME2D_SHELF_DEPTH", dome2d_depth_bay, &
57  'Depth of shelf, as fraction of basin depth, in 2d DOME configuration.', &
58  units='nondim',default=0.2)
59 
60  ! location where downslope starts
61  l1 = dome2d_width_bay
62 
63  ! location where downslope reaches maximum depth
64  l2 = 1.0 - dome2d_width_bottom
65 
66  bay_depth = dome2d_depth_bay
67 
68  do i=g%isc,g%iec
69  do j=g%jsc,g%jec
70 
71  ! Compute normalized zonal coordinate
72  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon;
73 
74  if ( x .le. l1 ) then
75  d(i,j) = bay_depth * max_depth
76  else if (( x .gt. l1 ) .and. ( x .lt. l2 )) then
77  d(i,j) = bay_depth * max_depth + (1.0-bay_depth) * max_depth * &
78  ( x - l1 ) / (l2 - l1)
79  else
80  d(i,j) = max_depth
81  end if
82 
83  enddo
84  enddo
85 end subroutine dome2d_initialize_topography
86 
87 !> Initialize thicknesses according to coordinate mode
88 subroutine dome2d_initialize_thickness ( h, G, GV, param_file, just_read_params )
89  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
90  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
91  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
92  intent(out) :: h !< The thickness that is being initialized, in m.
93  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
94  !! to parse for model parameter values.
95  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
96  !! only read parameters without changing h.
97 
98  ! Local variables
99  real :: e0(szk_(gv)) ! The resting interface heights, in m, usually !
100  ! negative because it is positive upward. !
101  real :: eta1D(szk_(gv)+1)! Interface height relative to the sea surface !
102  ! positive upward, in m. !
103  integer :: i, j, k, is, ie, js, je, nz
104  real :: x
105  real :: delta_h
106  real :: min_thickness
107  real :: dome2d_width_bay, dome2d_width_bottom, dome2d_depth_bay
108  logical :: just_read ! If true, just read parameters but set nothing.
109  character(len=40) :: verticalCoordinate
110 
111  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
112 
113  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
114 
115  if (.not.just_read) &
116  call mom_mesg("MOM_initialization.F90, DOME2d_initialize_thickness: setting thickness")
117 
118  call get_param(param_file, mdl,"MIN_THICKNESS",min_thickness, &
119  default=1.e-3, do_not_log=.true.)
120  call get_param(param_file, mdl,"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
121  default=default_coordinate_mode, do_not_log=.true.)
122  call get_param(param_file, mdl, "DOME2D_SHELF_WIDTH", dome2d_width_bay, &
123  default=0.1, do_not_log=.true.)
124  call get_param(param_file, mdl, "DOME2D_BASIN_WIDTH", dome2d_width_bottom, &
125  default=0.3, do_not_log=.true.)
126  call get_param(param_file, mdl, "DOME2D_SHELF_DEPTH", dome2d_depth_bay, &
127  default=0.2, do_not_log=.true.)
128 
129  if (just_read) return ! All run-time parameters have been read, so return.
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
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 )
145 
146  do j=js,je ; do i=is,ie
147  eta1d(nz+1) = -1.0*g%bathyT(i,j)
148  do k=nz,1,-1
149  eta1d(k) = e0(k)
150  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_z)) then
151  eta1d(k) = eta1d(k+1) + gv%Angstrom_z
152  h(i,j,k) = gv%Angstrom_z
153  else
154  h(i,j,k) = eta1d(k) - eta1d(k+1)
155  endif
156  enddo
157 
158  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon;
159  if ( x .le. dome2d_width_bay ) then
160  h(i,j,1:nz-1) = gv%Angstrom;
161  h(i,j,nz) = dome2d_depth_bay * g%max_depth - (nz-1) * gv%Angstrom;
162  end if
163 
164  end do ; end do
165 
166  ! case ( IC_RHO_C )
167  !
168  ! do j=js,je ; do i=is,ie
169  ! eta1D(nz+1) = -1.0*G%bathyT(i,j)
170  ! do k=nz,1,-1
171  ! eta1D(k) = e0(k)
172  ! if (eta1D(k) < (eta1D(k+1) + min_thickness)) then
173  ! eta1D(k) = eta1D(k+1) + min_thickness
174  ! h(i,j,k) = min_thickness
175  ! else
176  ! h(i,j,k) = eta1D(k) - eta1D(k+1)
177  ! endif
178  ! enddo
179  !
180  ! x = G%geoLonT(i,j) / G%len_lon;
181  ! if ( x .le. dome2d_width_bay ) then
182  ! h(i,j,1:nz-1) = min_thickness;
183  ! h(i,j,nz) = dome2d_depth_bay * G%max_depth - (nz-1) * min_thickness;
184  ! end if
185  !
186  ! enddo ; enddo
187 
188  case ( regridding_zstar )
189 
190  do j=js,je ; do i=is,ie
191  eta1d(nz+1) = -1.0*g%bathyT(i,j)
192  do k=nz,1,-1
193  eta1d(k) = e0(k)
194  if (eta1d(k) < (eta1d(k+1) + min_thickness)) then
195  eta1d(k) = eta1d(k+1) + min_thickness
196  h(i,j,k) = min_thickness
197  else
198  h(i,j,k) = eta1d(k) - eta1d(k+1)
199  endif
200  enddo
201  enddo ; enddo
202 
203  case ( regridding_sigma )
204  do j=js,je ; do i=is,ie
205  delta_h = g%bathyT(i,j) / nz;
206  h(i,j,:) = delta_h;
207  end do ; end do
208 
209  case default
210  call mom_error(fatal,"dome2d_initialize: "// &
211  "Unrecognized i.c. setup - set REGRIDDING_COORDINATE_MODE")
212 
213  end select
214 
215 end subroutine dome2d_initialize_thickness
216 
217 
218 !> Initialize temperature and salinity in the 2d DOME configuration
219 subroutine dome2d_initialize_temperature_salinity ( T, S, h, G, param_file, &
220  eqn_of_state, just_read_params)
221  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
222  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: T !< Potential temperature (degC)
223  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: S !< Salinity (ppt)
224  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Layer thickness (m or Pa)
225  type(param_file_type), intent(in) :: param_file !< Parameter file structure
226  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
227  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
228  !! only read parameters without changing h.
229 
230  ! Local variables
231  integer :: i, j, k, is, ie, js, je, nz
232  real :: x;
233  integer :: index_bay_z;
234  real :: delta_S, delta_T;
235  real :: S_ref, T_ref; ! Reference salinity and temperature within surface layer
236  real :: S_range, T_range; ! Range of salinities and temperatures over the vertical
237  real :: xi0, xi1;
238  logical :: just_read ! If true, just read parameters but set nothing.
239  character(len=40) :: verticalCoordinate
240  real :: dome2d_width_bay, dome2d_width_bottom, dome2d_depth_bay
241 
242  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
243 
244  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
245 
246  call get_param(param_file, mdl,"REGRIDDING_COORDINATE_MODE", verticalcoordinate, &
247  default=default_coordinate_mode, do_not_log=.true.)
248  call get_param(param_file, mdl, "DOME2D_SHELF_WIDTH", dome2d_width_bay, &
249  default=0.1, do_not_log=.true.)
250  call get_param(param_file, mdl, "DOME2D_BASIN_WIDTH", dome2d_width_bottom, &
251  default=0.3, do_not_log=.true.)
252  call get_param(param_file, mdl, "DOME2D_SHELF_DEPTH", dome2d_depth_bay, &
253  default=0.2, do_not_log=.true.)
254  call get_param(param_file, mdl,"S_REF",s_ref,'Reference salinity',units='1e-3', &
255  fail_if_missing=.not.just_read, do_not_log=just_read)
256  call get_param(param_file, mdl,"T_REF",t_ref,'Refernce temperature',units='C', &
257  fail_if_missing=.not.just_read, do_not_log=just_read)
258  call get_param(param_file, mdl,"S_RANGE",s_range,'Initial salinity range', &
259  units='1e-3', default=2.0, do_not_log=just_read)
260  call get_param(param_file, mdl,"T_RANGE",t_range,'Initial temperature range', &
261  units='1e-3', default=0.0, do_not_log=just_read)
262 
263  if (just_read) return ! All run-time parameters have been read, so return.
264 
265  t(:,:,:) = 0.0
266  s(:,:,:) = 0.0
267 
268  ! Linear salinity profile
269 
270  select case ( coordinatemode(verticalcoordinate) )
271 
272  case ( regridding_zstar, regridding_sigma )
273 
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  s(i,j,k) = 34.0 + 0.5 * s_range * (xi0 + xi1);
279  xi0 = xi1;
280  enddo
281  enddo ; enddo
282 
283  case ( regridding_rho )
284 
285  do j=js,je ; do i=is,ie
286  xi0 = 0.0;
287  do k = 1,nz
288  xi1 = xi0 + h(i,j,k) / g%max_depth;
289  s(i,j,k) = 34.0 + 0.5 * s_range * (xi0 + xi1);
290  xi0 = xi1;
291  enddo
292  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon;
293  if ( x .le. dome2d_width_bay ) then
294  s(i,j,nz) = 34.0 + s_range;
295  endif
296  enddo ; enddo
297 
298  case ( regridding_layer )
299 
300  delta_s = s_range / ( g%ke - 1.0 );
301  s(:,:,1) = s_ref;
302  do k = 2,g%ke
303  s(:,:,k) = s(:,:,k-1) + delta_s;
304  enddo
305 
306  case default
307  call mom_error(fatal,"dome2d_initialize: "// &
308  "Unrecognized i.c. setup - set REGRIDDING_COORDINATE_MODE")
309 
310  end select
311 
312  ! Modify salinity and temperature when z coordinates are used
313  if ( coordinatemode(verticalcoordinate) .eq. regridding_zstar ) then
314  index_bay_z = nint( dome2d_depth_bay * g%ke );
315  do j = g%jsc,g%jec ; do i = g%isc,g%iec
316  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon;
317  if ( x .le. dome2d_width_bay ) then
318  s(i,j,1:index_bay_z) = s_ref + s_range; ! Use for z coordinates
319  t(i,j,1:index_bay_z) = 1.0; ! Use for z coordinates
320  endif
321  enddo ; enddo ! i and j loops
322  endif ! Z initial conditions
323 
324  ! Modify salinity and temperature when sigma coordinates are used
325  if ( coordinatemode(verticalcoordinate) .eq. regridding_sigma ) then
326  do i = g%isc,g%iec ; do j = g%jsc,g%jec
327  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon;
328  if ( x .le. dome2d_width_bay ) then
329  s(i,j,1:g%ke) = s_ref + s_range; ! Use for sigma coordinates
330  t(i,j,1:g%ke) = 1.0; ! Use for sigma coordinates
331  endif
332  enddo ; enddo
333  endif
334 
335  ! Modify temperature when rho coordinates are used
336  t(g%isc:g%iec,g%jsc:g%jec,1:g%ke) = 0.0;
337  if (( coordinatemode(verticalcoordinate) .eq. regridding_rho ) .or. ( coordinatemode(verticalcoordinate) .eq. regridding_layer )) then
338  do i = g%isc,g%iec ; do j = g%jsc,g%jec
339  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon;
340  if ( x .le. dome2d_width_bay ) then
341  t(i,j,g%ke) = 1.0;
342  end if
343  end do ; end do
344  end if
345 
347 
348 !> Set up sponges in 2d DOME configuration
349 subroutine dome2d_initialize_sponges(G, GV, tv, param_file, use_ALE, CSp, ACSp)
350  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
351  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
352  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
353  type(param_file_type), intent(in) :: param_file !< Parameter file structure
354  logical, intent(in) :: use_ALE !< If true, indicates model is in ALE mode
355  type(sponge_cs), pointer :: CSp !< Layer-mode sponge structure
356  type(ale_sponge_cs), pointer :: ACSp !< ALE-mode sponge structure
357  ! Local variables
358  real :: T(szi_(g),szj_(g),szk_(g)) ! A temporary array for temp
359  real :: S(szi_(g),szj_(g),szk_(g)) ! A temporary array for salt
360  real :: RHO(szi_(g),szj_(g),szk_(g)) ! A temporary array for RHO
361  real :: h(szi_(g),szj_(g),szk_(g)) ! A temporary array for thickness
362  real :: eta(szi_(g),szj_(g),szk_(g)+1) ! A temporary array for thickness
363  real :: Idamp(szi_(g),szj_(g)) ! The inverse damping rate, in s-1.
364  real :: S_ref, T_ref ! Reference salinity and temerature within surface layer
365  real :: S_range, T_range ! Range of salinities and temperatures over the vertical
366  real :: e0(szk_(g)+1) ! The resting interface heights, in m, usually !
367  ! negative because it is positive upward. !
368  real :: eta1D(szk_(g)+1) ! Interface height relative to the sea surface !
369  ! positive upward, in m.
370  real :: dome2d_width_bay, dome2d_width_bottom, dome2d_depth_bay
371  real :: dome2d_west_sponge_time_scale, dome2d_east_sponge_time_scale
372  real :: dome2d_west_sponge_width, dome2d_east_sponge_width
373  real :: dummy1, x, z
374  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
375 
376  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
377  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
378 
379  call get_param(param_file, mdl, "DOME2D_WEST_SPONGE_TIME_SCALE", dome2d_west_sponge_time_scale, &
380  'The time-scale on the west edge of the domain for restoring T/S\n' //&
381  'in the sponge. If zero, the western sponge is disabled', &
382  units='s', default=0.)
383  call get_param(param_file, mdl, "DOME2D_EAST_SPONGE_TIME_SCALE", dome2d_east_sponge_time_scale, &
384  'The time-scale on the east edge of the domain for restoring T/S\n' //&
385  'in the sponge. If zero, the eastern sponge is disabled', &
386  units='s', default=0.)
387  call get_param(param_file, mdl, "DOME2D_WEST_SPONGE_WIDTH", dome2d_west_sponge_width, &
388  'The fraction of the domain in which the western sponge for restoring T/S\n' //&
389  'is active.', &
390  units='nondim', default=0.1)
391  call get_param(param_file, mdl, "DOME2D_EAST_SPONGE_WIDTH", dome2d_east_sponge_width, &
392  'The fraction of the domain in which the eastern sponge for restoring T/S\n' //&
393  'is active.', &
394  units='nondim', default=0.1)
395 
396  ! Return if sponges are not in use
397  if (dome2d_west_sponge_time_scale <= 0. .and. dome2d_east_sponge_time_scale <= 0.) return
398 
399  if (associated(csp)) call mom_error(fatal, &
400  "DOME2d_initialize_sponges called with an associated control structure.")
401  if (associated(acsp)) call mom_error(fatal, &
402  "DOME2d_initialize_sponges called with an associated ALE-sponge control structure.")
403 
404  call get_param(param_file, mdl, "DOME2D_SHELF_WIDTH", dome2d_width_bay, &
405  default=0.1, do_not_log=.true.)
406  call get_param(param_file, mdl, "DOME2D_BASIN_WIDTH", dome2d_width_bottom, &
407  default=0.3, do_not_log=.true.)
408  call get_param(param_file, mdl, "DOME2D_SHELF_DEPTH", dome2d_depth_bay, &
409  default=0.2, do_not_log=.true.)
410  call get_param(param_file, mdl,"S_REF",s_ref)
411  call get_param(param_file, mdl,"T_REF",t_ref)
412  call get_param(param_file, mdl,"S_RANGE",s_range,default=2.0)
413  call get_param(param_file, mdl,"T_RANGE",t_range,default=0.0)
414 
415 
416  ! Set the inverse damping rate as a function of position
417  idamp(:,:) = 0.0
418  do j=js,je ; do i=is,ie
419  if (g%mask2dT(i,j) > 0.) then ! Only set damping rate for wet points
420  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon ! Non-dimensional position within domain (0,1)
421  if ( dome2d_west_sponge_time_scale > 0. .and. x < dome2d_west_sponge_width ) then
422  ! Within half the shelf width from the left edge
423  dummy1 = 1. - x / dome2d_west_sponge_width
424  idamp(i,j) = 1./dome2d_west_sponge_time_scale * max(0., min(1., dummy1))
425  elseif ( dome2d_east_sponge_time_scale > 0. .and. x > ( 1. - dome2d_east_sponge_width ) ) then
426  ! Within a quarter of the basin width from the right
427  dummy1 = 1. - ( 1. - x ) / dome2d_east_sponge_width
428  idamp(i,j) = 1./dome2d_east_sponge_time_scale * max(0., min(1., dummy1))
429  else
430  idamp(i,j) = 0.
431  endif
432  else
433  idamp(i,j) = 0.
434  endif
435  enddo ; enddo
436 
437 
438  if (use_ale) then
439 
440  ! Construct a grid (somewhat arbitrarily) to describe the sponge T/S on
441  do k=1,nz
442  e0(k) = -g%max_depth * ( real(k-1) / real(nz) )
443  enddo
444  e0(nz+1) = -g%max_depth
445  do j=js,je ; do i=is,ie
446  eta1d(nz+1) = -1.0*g%bathyT(i,j)
447  do k=nz,1,-1
448  eta1d(k) = e0(k)
449  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_z)) then
450  eta1d(k) = eta1d(k+1) + gv%Angstrom_z
451  h(i,j,k) = gv%Angstrom_z
452  else
453  h(i,j,k) = eta1d(k) - eta1d(k+1)
454  endif
455  enddo
456  enddo; enddo
457  ! Store the grid on which the T/S sponge data will reside
458  call initialize_ale_sponge(idamp, h, nz, g, param_file, acsp)
459 
460  ! Construct temperature and salinity on the arbitrary grid
461  t(:,:,:) = 0.0 ; s(:,:,:) = 0.0
462  do j=js,je ; do i=is,ie
463  z = -g%bathyT(i,j)
464  do k = nz,1,-1
465  z = z + 0.5 * h(i,j,k) ! Position of the center of layer k
466  s(i,j,k) = 34.0 - 1.0 * (z/g%max_depth)
467  if ( ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon < dome2d_west_sponge_width ) s(i,j,k) = s_ref + s_range
468  z = z + 0.5 * h(i,j,k) ! Position of the interface k
469  enddo
470  enddo ; enddo
471 
472  if ( associated(tv%T) ) then
473  call set_up_ale_sponge_field(t,g,tv%T,acsp)
474  endif
475  if ( associated(tv%S) ) then
476  call set_up_ale_sponge_field(s,g,tv%S,acsp)
477  endif
478 
479  else
480 
481  ! Construct thicknesses to restore to
482  do j=js,je ; do i=is,ie
483  eta1d(nz+1) = -1.0*g%bathyT(i,j)
484  do k=nz,1,-1
485  eta1d(k) = -g%max_depth * real(k-1) / real(nz)
486  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_z)) then
487  eta1d(k) = eta1d(k+1) + gv%Angstrom_z
488  h(i,j,k) = gv%Angstrom_z
489  else
490  h(i,j,k) = eta1d(k) - eta1d(k+1)
491  endif
492  enddo
493 
494  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon;
495  if ( x .le. dome2d_width_bay ) then
496  h(i,j,1:nz-1) = gv%Angstrom;
497  h(i,j,nz) = dome2d_depth_bay * g%max_depth - (nz-1) * gv%Angstrom;
498  end if
499 
500  eta(i,j,nz+1) = -g%bathyT(i,j)
501  do k=nz,1,-1
502  eta(i,j,k) = eta(i,j,k+1) + h(i,j,k)
503  enddo
504  enddo ; enddo
505  call initialize_sponge(idamp, eta, g, param_file, csp)
506 
507  endif
508 
509 end subroutine dome2d_initialize_sponges
510 
511 end module dome2d_initialization
integer, parameter regridding_layer
Layer mode.
subroutine, public read_axis_data(filename, axis_name, var)
Definition: MOM_io.F90:424
integer function coordinatemode(string)
Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer...
integer, parameter regridding_sigma
Sigma coordinates.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:45
Provides the ocean grid type.
Definition: MOM_grid.F90:2
character(len=40) mdl
This module&#39;s name.
character(len= *), parameter default_coordinate_mode
Default coordinate mode.
This module contains I/O framework code.
Definition: MOM_io.F90:2
subroutine, public set_up_ale_sponge_field(sp_val, G, f_ptr, CS)
This subroutine stores the reference profile at h points for the variable.
This module contains the routines used to apply sponge layers when using the ALE mode. Applying sponges requires the following: (1) initialize_ALE_sponge (2) set_up_ALE_sponge_field (tracers) and set_up_ALE_sponge_vel_field (vel) (3) apply_ALE_sponge (4) init_ALE_sponge_diags (not being used for now) (5) ALE_sponge_end (not being used for now)
SPONGE control structure.
subroutine, public dome2d_initialize_thickness(h, G, GV, param_file, just_read_params)
Initialize thicknesses according to coordinate mode.
subroutine, public calculate_density_derivs(T, S, pressure, drho_dT, drho_dS, start, npts, EOS)
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
Definition: MOM_EOS.F90:214
subroutine, public initialize_sponge(Iresttime, int_height, G, param_file, CS, Iresttime_i_mean, int_height_i_mean)
Definition: MOM_sponge.F90:142
subroutine, public set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
Definition: MOM_sponge.F90:271
subroutine, public dome2d_initialize_sponges(G, GV, tv, param_file, use_ALE, CSp, ACSp)
Set up sponges in 2d DOME configuration.
subroutine, public initialize_ale_sponge(Iresttime, data_h, nz_data, G, param_file, CS)
This subroutine determines the number of points which are within.
subroutine, public mom_mesg(message, verb, all_print)
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:51
integer, parameter regridding_zstar
z* coordinates
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
Contains constants for interpreting input parameters that control regridding.
subroutine, public dome2d_initialize_temperature_salinity(T, S, h, G, param_file, eqn_of_state, just_read_params)
Initialize temperature and salinity in the 2d DOME configuration.
subroutine, public dome2d_initialize_topography(D, G, param_file, max_depth)
Initialize topography with a shelf and slope in a 2D domain.
subroutine, public mom_error(level, message, all_print)
integer, parameter regridding_rho
Target interface densities.
A control structure for the equation of state.
Definition: MOM_EOS.F90:55