MOM6
adjustment_initialization.F90
1 !> Configures the model for the geostrophic adjustment test case.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_mesg, mom_error, fatal, is_root_pe
8 use mom_get_input, only : directories
9 use mom_grid, only : ocean_grid_type
14 use regrid_consts, only : coordinatemode, default_coordinate_mode
15 use regrid_consts, only : regridding_layer, regridding_zstar
16 use regrid_consts, only : regridding_rho, regridding_sigma
17 
18 implicit none ; private
19 
20 character(len=40) :: mdl = "adjustment_initialization" !< This module's name.
21 
22 #include <MOM_memory.h>
23 
24 public adjustment_initialize_thickness
25 public adjustment_initialize_temperature_salinity
26 
27 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
28 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
29 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
30 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
31 
32 contains
33 
34 !> Initializes the layer thicknesses in the adjustment test case
35 subroutine adjustment_initialize_thickness ( h, G, GV, US, param_file, just_read_params)
36  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
37  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
38  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
39  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
40  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
41  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
42  !! to parse for model parameter values.
43  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
44  !! only read parameters without changing h.
45  ! Local variables
46  real :: e0(szk_(g)+1) ! The resting interface heights, in depth units [Z ~> m], usually
47  ! negative because it is positive upward.
48  real :: eta1d(szk_(g)+1)! Interface height relative to the sea surface
49  ! positive upward, in depth units [Z ~> m].
50  real :: x, y, yy, delta_s_strat, dsdz, delta_s, s_ref
51  real :: min_thickness, adjustment_width, adjustment_delta, adjustment_deltas
52  real :: front_wave_amp, front_wave_length, front_wave_asym
53  real :: target_values(szk_(g)+1)
54  logical :: just_read ! If true, just read parameters but set nothing.
55  character(len=20) :: verticalcoordinate
56 ! This include declares and sets the variable "version".
57 #include "version_variable.h"
58  integer :: i, j, k, is, ie, js, je, nz
59 
60  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
61 
62  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
63 
64  if (.not.just_read) &
65  call mom_mesg("initialize_thickness_uniform: setting thickness")
66 
67  ! Parameters used by main model initialization
68  if (.not.just_read) call log_version(param_file, mdl, version, "")
69  call get_param(param_file, mdl,"S_REF",s_ref,fail_if_missing=.true.,do_not_log=.true.)
70  call get_param(param_file, mdl,"MIN_THICKNESS",min_thickness,'Minimum layer thickness', &
71  units='m', default=1.0e-3, do_not_log=just_read, scale=us%m_to_Z)
72 
73  ! Parameters specific to this experiment configuration
74  call get_param(param_file, mdl,"REGRIDDING_COORDINATE_MODE",verticalcoordinate, &
75  default=default_coordinate_mode, do_not_log=just_read)
76  call get_param(param_file, mdl,"ADJUSTMENT_WIDTH",adjustment_width, &
77  "Width of frontal zone", &
78  units="same as x,y", fail_if_missing=.not.just_read, do_not_log=just_read)
79  call get_param(param_file, mdl,"DELTA_S_STRAT",delta_s_strat, &
80  "Top-to-bottom salinity difference of stratification", &
81  units="1e-3", fail_if_missing=.not.just_read, do_not_log=just_read)
82  call get_param(param_file, mdl,"ADJUSTMENT_DELTAS",adjustment_deltas, &
83  "Salinity difference across front", &
84  units="1e-3", fail_if_missing=.not.just_read, do_not_log=just_read)
85  call get_param(param_file, mdl,"FRONT_WAVE_AMP",front_wave_amp, &
86  "Amplitude of trans-frontal wave perturbation", &
87  units="same as x,y",default=0., do_not_log=just_read)
88  call get_param(param_file, mdl,"FRONT_WAVE_LENGTH",front_wave_length, &
89  "Wave-length of trans-frontal wave perturbation", &
90  units="same as x,y",default=0., do_not_log=just_read)
91  call get_param(param_file, mdl,"FRONT_WAVE_ASYM",front_wave_asym, &
92  "Amplitude of frontal asymmetric perturbation", &
93  default=0., do_not_log=just_read)
94 
95  if (just_read) return ! All run-time parameters have been read, so return.
96 
97  ! WARNING: this routine specifies the interface heights so that the last layer
98  ! is vanished, even at maximum depth. In order to have a uniform
99  ! layer distribution, use this line of code within the loop:
100  ! e0(k) = -G%max_depth * real(k-1) / real(nz)
101  ! To obtain a thickness distribution where the last layer is
102  ! vanished and the other thicknesses uniformly distributed, use:
103  ! e0(k) = -G%max_depth * real(k-1) / real(nz-1)
104 
105  dsdz = -delta_s_strat / g%max_depth
106 
107  select case ( coordinatemode(verticalcoordinate) )
108 
109  case ( regridding_layer, regridding_rho )
110  if (delta_s_strat /= 0.) then
111  ! This was previously coded ambiguously.
112  adjustment_delta = (adjustment_deltas / delta_s_strat) * g%max_depth
113  do k=1,nz+1
114  e0(k) = adjustment_delta - (g%max_depth + 2*adjustment_delta) * (real(k-1) / real(nz))
115  enddo
116  else
117  adjustment_delta = 2.*g%max_depth
118  do k=1,nz+1
119  e0(k) = -g%max_depth * (real(k-1) / real(nz))
120  enddo
121  endif
122  target_values(1) = gv%Rlay(1) + 0.5*(gv%Rlay(1)-gv%Rlay(2))
123  target_values(nz+1) = gv%Rlay(nz) + 0.5*(gv%Rlay(nz)-gv%Rlay(nz-1))
124  do k = 2,nz
125  target_values(k) = target_values(k-1) + ( gv%Rlay(nz) - gv%Rlay(1) ) / (nz-1)
126  enddo
127  target_values(:) = target_values(:) - 1000.
128  do j=js,je ; do i=is,ie
129  if (front_wave_length /= 0.) then
130  y = ( 0.125 + g%geoLatT(i,j) / front_wave_length ) * ( 4. * acos(0.) )
131  yy = 2. * ( g%geoLatT(i,j) - 0.5 * g%len_lat ) / adjustment_width
132  yy = min(1.0, yy); yy = max(-1.0, yy)
133  yy = yy * 2. * acos( 0. )
134  y = front_wave_amp*sin(y) + front_wave_asym*sin(yy)
135  else
136  y = 0.
137  endif
138  x = ( ( g%geoLonT(i,j) - 0.5 * g%len_lon ) + y ) / adjustment_width
139  x = min(1.0, x); x = max(-1.0, x)
140  x = x * acos( 0. )
141  delta_s = adjustment_deltas * 0.5 * (1. - sin( x ) )
142  do k=2,nz
143  if (dsdz /= 0.) then
144  eta1d(k) = ( target_values(k) - ( s_ref + delta_s ) ) / dsdz
145  else
146  eta1d(k) = e0(k) - (0.5*adjustment_delta) * sin( x )
147  endif
148  eta1d(k) = max( eta1d(k), -g%max_depth )
149  eta1d(k) = min( eta1d(k), 0. )
150  enddo
151  eta1d(1) = 0.; eta1d(nz+1) = -g%max_depth
152  do k=nz,1,-1
153  if (eta1d(k) > 0.) then
154  eta1d(k) = max( eta1d(k+1) + min_thickness, 0. )
155  h(i,j,k) = gv%Z_to_H * max( eta1d(k) - eta1d(k+1), min_thickness )
156  elseif (eta1d(k) <= (eta1d(k+1) + min_thickness)) then
157  eta1d(k) = eta1d(k+1) + min_thickness
158  h(i,j,k) = gv%Z_to_H * min_thickness
159  else
160  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
161  endif
162  enddo
163  enddo ; enddo
164 
165  case ( regridding_zstar, regridding_sigma )
166  do k=1,nz+1
167  eta1d(k) = -g%max_depth * (real(k-1) / real(nz))
168  eta1d(k) = max(min(eta1d(k), 0.), -g%max_depth)
169  enddo
170  do j=js,je ; do i=is,ie
171  do k=nz,1,-1
172  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
173  enddo
174  enddo ; enddo
175 
176  case default
177  call mom_error(fatal,"adjustment_initialize_thickness: "// &
178  "Unrecognized i.c. setup - set ADJUSTMENT_IC")
179 
180  end select
181 
182 end subroutine adjustment_initialize_thickness
183 
184 !> Initialization of temperature and salinity in the adjustment test case
185 subroutine adjustment_initialize_temperature_salinity(T, S, h, G, GV, param_file, &
186  eqn_of_state, just_read_params)
187  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
188  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
189  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: t !< The temperature that is being initialized.
190  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: s !< The salinity that is being initialized.
191  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< The model thicknesses [H ~> m or kg m-2].
192  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to
193  !! parse for model parameter values.
194  type(eos_type), pointer :: eqn_of_state !< Equation of state.
195  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
196  !! only read parameters without changing T & S.
197 
198  integer :: i, j, k, is, ie, js, je, nz
199  real :: x, y, yy
200  integer :: index_bay_z
201  real :: s_ref, t_ref ! Reference salinity and temerature within
202  ! surface layer
203  real :: s_range, t_range ! Range of salinities and temperatures over the
204  ! vertical
205  real :: xi0, xi1, dsdz, delta_s, delta_s_strat
206  real :: adjustment_width, adjustment_deltas
207  real :: front_wave_amp, front_wave_length, front_wave_asym
208  real :: eta1d(szk_(g)+1)
209  logical :: just_read ! If true, just read parameters but set nothing.
210  character(len=20) :: verticalcoordinate
211 
212  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
213 
214  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
215 
216  ! Parameters used by main model initialization
217  call get_param(param_file, mdl,"S_REF",s_ref,'Reference salinity', units='1e-3', &
218  fail_if_missing=.not.just_read, do_not_log=just_read)
219  call get_param(param_file, mdl,"T_REF",t_ref,'Reference temperature', units='C', &
220  fail_if_missing=.not.just_read, do_not_log=just_read)
221  call get_param(param_file, mdl,"S_RANGE",s_range,'Initial salinity range', units='1e-3', &
222  default=2.0, do_not_log=just_read)
223  call get_param(param_file, mdl,"T_RANGE",t_range,'Initial temperature range', units='C', &
224  default=0.0, do_not_log=just_read)
225  ! Parameters specific to this experiment configuration BUT logged in previous s/r
226  call get_param(param_file, mdl,"REGRIDDING_COORDINATE_MODE",verticalcoordinate, &
227  default=default_coordinate_mode, do_not_log=just_read)
228  call get_param(param_file, mdl,"ADJUSTMENT_WIDTH", adjustment_width, &
229  fail_if_missing=.not.just_read, do_not_log=.true.)
230  call get_param(param_file, mdl,"ADJUSTMENT_DELTAS", adjustment_deltas, &
231  fail_if_missing=.not.just_read, do_not_log=.true.)
232  call get_param(param_file, mdl,"DELTA_S_STRAT", delta_s_strat, &
233  fail_if_missing=.not.just_read, do_not_log=.true.)
234  call get_param(param_file, mdl,"FRONT_WAVE_AMP", front_wave_amp, default=0., &
235  do_not_log=.true.)
236  call get_param(param_file, mdl,"FRONT_WAVE_LENGTH",front_wave_length, &
237  default=0., do_not_log=.true.)
238  call get_param(param_file, mdl,"FRONT_WAVE_ASYM", front_wave_asym, default=0., &
239  do_not_log=.true.)
240 
241  if (just_read) return ! All run-time parameters have been read, so return.
242 
243  t(:,:,:) = 0.0
244  s(:,:,:) = 0.0
245 
246  ! Linear salinity profile
247  select case ( coordinatemode(verticalcoordinate) )
248 
249  case ( regridding_zstar, regridding_sigma )
250  dsdz = -delta_s_strat / g%max_depth
251  do j=js,je ; do i=is,ie
252  eta1d(nz+1) = -g%bathyT(i,j)
253  do k=nz,1,-1
254  eta1d(k) = eta1d(k+1) + h(i,j,k)*gv%H_to_Z
255  enddo
256  if (front_wave_length /= 0.) then
257  y = ( 0.125 + g%geoLatT(i,j) / front_wave_length ) * ( 4. * acos(0.) )
258  yy = 2. * ( g%geoLatT(i,j) - 0.5 * g%len_lat ) / front_wave_length
259  yy = min(1.0, yy); yy = max(-1.0, yy)
260  yy = yy * 2. * acos( 0. )
261  y = front_wave_amp*sin(y) + front_wave_asym*sin(yy)
262  else
263  y = 0.
264  endif
265  x = ( ( g%geoLonT(i,j) - 0.5 * g%len_lon ) + y ) / adjustment_width
266  x = min(1.0, x); x = max(-1.0, x)
267  x = x * acos( 0. )
268  delta_s = adjustment_deltas * 0.5 * (1. - sin( x ) )
269  do k=1,nz
270  s(i,j,k) = s_ref + delta_s + 0.5 * ( eta1d(k)+eta1d(k+1) ) * dsdz
271  x = abs(s(i,j,k) - 0.5*real(nz-1)/real(nz)*s_range)/s_range*real(2*nz)
272  x = 1. - min(1., x)
273  t(i,j,k) = x
274  enddo
275  ! x = GV%H_to_Z*sum(T(i,j,:)*h(i,j,:))
276  ! T(i,j,:) = (T(i,j,:) / x) * (G%max_depth*1.5/real(nz))
277  enddo ; enddo
278 
279  case ( regridding_layer, regridding_rho )
280  do k = 1,nz
281  s(:,:,k) = s_ref + s_range * ( (real(k)-0.5) / real( nz ) )
282  ! x = abs(S(1,1,k) - 0.5*real(nz-1)/real(nz)*S_range)/S_range*real(2*nz)
283  ! x = 1.-min(1., x)
284  ! T(:,:,k) = x
285  enddo
286 
287  case default
288  call mom_error(fatal,"adjustment_initialize_temperature_salinity: "// &
289  "Unrecognized i.c. setup - set ADJUSTMENT_IC")
290 
291  end select
292 
293 end subroutine adjustment_initialize_temperature_salinity
294 
295 end module adjustment_initialization
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:82
mom_get_input::directories
Container for paths and parameter file names.
Definition: MOM_get_input.F90:20
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_get_input
Reads the only Fortran name list needed to boot-strap the model.
Definition: MOM_get_input.F90:6
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
regrid_consts
Contains constants for interpreting input parameters that control regridding.
Definition: regrid_consts.F90:2
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:86
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_eos::calculate_density_derivs
Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition: MOM_EOS.F90:70
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
adjustment_initialization
Configures the model for the geostrophic adjustment test case.
Definition: adjustment_initialization.F90:2
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:60