MOM6
external_gwave_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 
24 use mom_get_input, only : directories
25 use mom_grid, only : ocean_grid_type
28 implicit none ; private
29 
30 #include <MOM_memory.h>
31 
33 
34 contains
35 
36 ! -----------------------------------------------------------------------------
37 !> This subroutine initializes layer thicknesses for the external_gwave experiment.
38 subroutine external_gwave_initialize_thickness(h, G, param_file, just_read_params)
39  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
40  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
41  intent(out) :: h !< The thickness that is being initialized, in m.
42  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
43  !! to parse for model parameter values.
44  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
45  !! only read parameters without changing h.
46 
47  real :: e0(szk_(g)) ! The resting interface heights, in m, usually !
48  ! negative because it is positive upward. !
49  real :: e_pert(szk_(g)) ! Interface height perturbations, positive !
50  ! upward, in m. !
51  real :: eta1D(szk_(g)+1)! Interface height relative to the sea surface !
52  ! positive upward, in m. !
53  real :: ssh_anomaly_height ! Vertical height of ssh anomaly
54  real :: ssh_anomaly_width ! Lateral width of anomaly
55  logical :: just_read ! If true, just read parameters but set nothing.
56  character(len=40) :: mdl = "external_gwave_initialize_thickness" ! This subroutine's name.
57 ! This include declares and sets the variable "version".
58 #include "version_variable.h"
59  integer :: i, j, k, is, ie, js, je, nz
60  real :: PI, Xnondim
61 
62  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
63 
64  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
65 
66  if (.not.just_read) &
67  call mom_mesg(" external_gwave_initialization.F90, external_gwave_initialize_thickness: setting thickness", 5)
68 
69  if (.not.just_read) call log_version(param_file, mdl, version, "")
70  call get_param(param_file, mdl, "SSH_ANOMALY_HEIGHT", ssh_anomaly_height, &
71  "The vertical displacement of the SSH anomaly. ", units="m", &
72  fail_if_missing=.not.just_read, do_not_log=just_read)
73  call get_param(param_file, mdl, "SSH_ANOMALY_WIDTH", ssh_anomaly_width, &
74  "The lateral width of the SSH anomaly. ", units="coordinate", &
75  fail_if_missing=.not.just_read, do_not_log=just_read)
76 
77  if (just_read) return ! All run-time parameters have been read, so return.
78 
79  pi = 4.0*atan(1.0)
80  do j=g%jsc,g%jec ; do i=g%isc,g%iec
81  xnondim = (g%geoLonT(i,j)-g%west_lon-0.5*g%len_lon) / ssh_anomaly_width
82  xnondim = min(1., abs(xnondim))
83  eta1d(1) = ssh_anomaly_height * 0.5 * ( 1. + cos(pi*xnondim) ) ! Cosine bell
84  do k=2,nz
85  eta1d(k) = -g%max_depth & ! Stretch interior interfaces with SSH
86  + (eta1d(1)+g%max_depth) * ( real(nz+1-k)/real(nz) ) ! Stratification
87  enddo
88  eta1d(nz+1) = -g%max_depth ! Force bottom interface to bottom
89  do k=1,nz
90  h(i,j,k) = eta1d(k) - eta1d(k+1)
91  enddo
92  enddo ; enddo
93 
95 ! -----------------------------------------------------------------------------
96 
97 !> \namespace external_gwave_initialization
98 !!
99 !! The module configures the model for the "external_gwave" experiment.
100 !! external_gwave = External Gravity Wave
The module configures the model for the "external_gwave" experiment. external_gwave = External Gravit...
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Provides the ocean grid type.
Definition: MOM_grid.F90:2
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Type to carry basic tracer information.
logical function, public is_root_pe()
subroutine, public mom_mesg(message, verb, all_print)
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
subroutine, public external_gwave_initialize_thickness(h, G, param_file, just_read_params)
This subroutine initializes layer thicknesses for the external_gwave experiment.
subroutine, public mom_error(level, message, all_print)