MOM6
sloshing_initialization.F90
Go to the documentation of this file.
2 !***********************************************************************
3 !* GNU General Public License *
4 !* This file is a part of MOM. *
5 !* *
6 !* MOM is free software; you can redistribute it and/or modify it and *
7 !* are expected to follow the terms of the GNU General Public License *
8 !* as published by the Free Software Foundation; either version 2 of *
9 !* the License, or (at your option) any later version. *
10 !* *
11 !* MOM is distributed in the hope that it will be useful, but WITHOUT *
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
13 !* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
14 !* License for more details. *
15 !* *
16 !* For the full text of the GNU General Public License, *
17 !* write to: Free Software Foundation, Inc., *
18 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
19 !* or see: http://www.gnu.org/licenses/gpl.html *
20 !***********************************************************************
21 
22 use mom_domains, only : sum_across_pes
26 use mom_get_input, only : directories
27 use mom_grid, only : ocean_grid_type
28 use mom_io, only : close_file, fieldtype, file_exists
29 use mom_io, only : open_file, read_data, read_axis_data, single_file
30 use mom_io, only : write_field, slasher
36 
37 implicit none ; private
38 
39 #include <MOM_memory.h>
40 
41 ! -----------------------------------------------------------------------------
42 ! The following routines are visible to the outside world
43 ! -----------------------------------------------------------------------------
47 
48 character(len=40) :: mdl = "sloshing_initialization" ! This module's name.
49 
50 ! -----------------------------------------------------------------------------
51 ! This module contains the following routines
52 ! -----------------------------------------------------------------------------
53 contains
54 
55 !> Initialization of topography.
56 subroutine sloshing_initialize_topography ( D, G, param_file, max_depth )
57  type(dyn_horgrid_type), intent(in) :: G !< The dynamic horizontal grid type
58  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
59  intent(out) :: D !< Ocean bottom depth in m
60  type(param_file_type), intent(in) :: param_file !< Parameter file structure
61  real, intent(in) :: max_depth !< Maximum depth of model in m
62 
63  ! Local variables
64  integer :: i, j
65 
66  do i=g%isc,g%iec
67  do j=g%jsc,g%jec
68 
69  d(i,j) = max_depth
70 
71  enddo
72  enddo
73 
74 end subroutine sloshing_initialize_topography
75 
76 
77 !> Initialization of thicknesses
78 !! This routine is called when THICKNESS_CONFIG is set to 'sloshing'
79 !!
80 !! This routine initializes layer positions to set off a sloshing motion in
81 !! the zonal direction in a rectangular basin. All layers have initially the
82 !! same thickness but all interfaces (except bottom and sea surface) are
83 !! displaced according to a half-period cosine, with maximum value on the
84 !! left and minimum value on the right. This sets off a regular sloshing motion.
85 subroutine sloshing_initialize_thickness ( h, G, GV, param_file, just_read_params)
86  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
87  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
88  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
89  intent(out) :: h !< The thickness that is being initialized, in m.
90  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
91  !! to parse for model parameter values.
92  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
93  !! only read parameters without changing h.
94 
95  real :: displ(szk_(g)+1)
96  real :: z_unif(szk_(g)+1)
97  real :: z_inter(szk_(g)+1)
98  real :: x
99  real :: a0
100  real :: deltah
101  real :: total_height
102  real :: weight_z
103  real :: x1, y1, x2, y2
104  real :: t
105  logical :: just_read ! If true, just read parameters but set nothing.
106  integer :: n
107 
108  integer :: i, j, k, is, ie, js, je, nx, nz
109 
110  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
111 
112  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
113 
114  if (just_read) return ! This subroutine has no run-time parameters.
115 
116  deltah = g%max_depth / nz
117 
118  ! Define thicknesses
119  do j=g%jsc,g%jec ; do i=g%isc,g%iec
120 
121  ! Define uniform interfaces
122  do k = 0,nz
123  z_unif(k+1) = -real(k)/real(nz)
124  end do
125 
126  ! 1. Define stratification
127  n = 3
128  do k = 1,nz+1
129 
130  ! Thin pycnocline in the middle
131  !z_inter(k) = (2.0**(n-1)) * (z_unif(k) + 0.5)**n - 0.5
132 
133  ! Thin pycnocline in the middle (piecewise linear profile)
134  x1 = 0.30; y1 = 0.48; x2 = 0.70; y2 = 0.52
135 
136  x = -z_unif(k)
137 
138  if ( x .le. x1 ) then
139  t = y1*x/x1;
140  else if ( (x .gt. x1 ) .and. ( x .lt. x2 )) then
141  t = y1 + (y2-y1) * (x-x1) / (x2-x1)
142  else
143  t = y2 + (1.0-y2) * (x-x2) / (1.0-x2)
144  end if
145 
146  t = - z_unif(k)
147 
148  z_inter(k) = -t * g%max_depth
149 
150  end do
151 
152  ! 2. Define displacement
153  a0 = 75.0; ! Displacement amplitude (meters)
154  do k = 1,nz+1
155 
156  weight_z = - 4.0 * ( z_unif(k) + 0.5 )**2 + 1
157 
158  x = g%geoLonT(i,j) / g%len_lon
159  displ(k) = a0 * cos(acos(-1.0)*x) + weight_z;
160 
161  if ( k .EQ. 1 ) then
162  displ(k) = 0.0
163  end if
164 
165  if ( k .EQ. nz+1 ) then
166  displ(k) = 0.0
167  end if
168 
169  z_inter(k) = z_inter(k) + displ(k)
170 
171  end do
172 
173  ! 3. The last interface must coincide with the seabed
174  z_inter(nz+1) = -g%bathyT(i,j)
175 
176  ! Modify interface heights to make sure all thicknesses
177  ! are strictly positive
178  do k = nz,1,-1
179 
180  if ( z_inter(k) .LT. (z_inter(k+1) + gv%Angstrom) ) then
181  z_inter(k) = z_inter(k+1) + gv%Angstrom
182  end if
183 
184  end do
185 
186  ! 4. Define layers
187  total_height = 0.0
188  do k = 1,nz
189  h(i,j,k) = z_inter(k) - z_inter(k+1)
190 
191  total_height = total_height + h(i,j,k)
192  end do
193 
194  enddo ; enddo
195 
196 end subroutine sloshing_initialize_thickness
197 
198 
199 !------------------------------------------------------------------------------
200 !> Initialization of temperature and salinity
201 !!
202 !! This subroutine initializes linear profiles for T and S according to
203 !! reference surface layer salinity and temperature and a specified range.
204 !! Note that the linear distribution is set up with respect to the layer
205 !! number, not the physical position).
206 subroutine sloshing_initialize_temperature_salinity ( T, S, h, G, param_file, &
207  eqn_of_state, just_read_params)
208  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure.
209  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: T !< Potential temperature (degC).
210  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: S !< Salinity (ppt).
211  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Layer thickness (m or Pa).
212  type(param_file_type), intent(in) :: param_file !< A structure indicating the
213  !! open file to parse for model
214  !! parameter values.
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  integer :: i, j, k, is, ie, js, je, nz
220  real :: delta_S, delta_T
221  real :: S_ref, T_ref; ! Reference salinity and temerature within
222  ! surface layer
223  real :: S_range, T_range; ! Range of salinities and temperatures over the
224  ! vertical
225  integer :: kdelta
226  real :: deltah
227  real :: xi0, xi1
228  logical :: just_read ! If true, just read parameters but set nothing.
229  character(len=40) :: mdl = "initialize_temp_salt_linear" ! This subroutine's
230  ! name.
231 
232  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
233 
234  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
235 
236  call get_param(param_file, mdl,"S_REF",s_ref,'Reference value for salinity', &
237  units='1e-3', fail_if_missing=.not.just_read, do_not_log=just_read)
238  call get_param(param_file, mdl,"T_REF",t_ref,'Refernce value for temperature', &
239  units='C', fail_if_missing=.not.just_read, do_not_log=just_read)
240 
241  ! The default is to assume an increase by 2 for the salinity and a uniform
242  ! temperature
243  call get_param(param_file, mdl,"S_RANGE",s_range,'Initial salinity range.', &
244  units='1e-3', default=2.0, do_not_log=just_read)
245  call get_param(param_file, mdl,"T_RANGE",t_range,'Initial temperature range', &
246  units='C', default=0.0, do_not_log=just_read)
247 
248  if (just_read) return ! All run-time parameters have been read, so return.
249 
250  ! Prescribe salinity
251  !delta_S = S_range / ( G%ke - 1.0 )
252 
253  !S(:,:,1) = S_ref
254  !do k = 2,G%ke
255  ! S(:,:,k) = S(:,:,k-1) + delta_S
256  !end do
257 
258  deltah = g%max_depth / nz;
259  do j=js,je ; do i=is,ie
260  xi0 = 0.0
261  do k = 1,nz
262  xi1 = xi0 + deltah / g%max_depth
263  s(i,j,k) = 34.0 + 0.5 * s_range * (xi0 + xi1)
264  xi0 = xi1
265  enddo
266  enddo ; enddo
267 
268  ! Prescribe temperature
269  delta_t = t_range / ( g%ke - 1.0 )
270 
271  t(:,:,1) = t_ref
272  do k = 2,g%ke
273  t(:,:,k) = t(:,:,k-1) + delta_t
274  end do
275  kdelta = 2
276  t(:,:,g%ke/2 - (kdelta-1):g%ke/2 + kdelta) = 1.0
277 
279 
280 !> \namespace sloshing_initialization
281 !!
282 !! The module configures the model for the non-rotating sloshing
283 !! test case.
284 end module sloshing_initialization
subroutine, public read_axis_data(filename, axis_name, var)
Definition: MOM_io.F90:424
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
This module contains I/O framework code.
Definition: MOM_io.F90:2
subroutine, public sloshing_initialize_thickness(h, G, GV, param_file, just_read_params)
Initialization of thicknesses This routine is called when THICKNESS_CONFIG is set to &#39;sloshing&#39;...
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
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 sloshing_initialize_temperature_salinity(T, S, h, G, param_file, eqn_of_state, just_read_params)
Initialization of temperature and salinity.
subroutine, public set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
Definition: MOM_sponge.F90:271
Type to carry basic tracer information.
logical function, public is_root_pe()
The module configures the model for the non-rotating sloshing test case.
subroutine, public mom_mesg(message, verb, all_print)
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
subroutine, public mom_error(level, message, all_print)
subroutine, public sloshing_initialize_topography(D, G, param_file, max_depth)
Initialization of topography.
A control structure for the equation of state.
Definition: MOM_EOS.F90:55