MOM6
external_gwave_initialization Module Reference

Detailed Description

The module configures the model for the "external_gwave" experiment. external_gwave = External Gravity Wave.

Functions/Subroutines

subroutine, public external_gwave_initialize_thickness (h, G, param_file, just_read_params)
 This subroutine initializes layer thicknesses for the external_gwave experiment. More...
 

Function/Subroutine Documentation

◆ external_gwave_initialize_thickness()

subroutine, public external_gwave_initialization::external_gwave_initialize_thickness ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(out)  h,
type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
logical, intent(in), optional  just_read_params 
)

This subroutine initializes layer thicknesses for the external_gwave experiment.

Parameters
[in]gThe ocean's 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 39 of file external_gwave_initialization.F90.

References mom_error_handler::mom_mesg().

Referenced by mom_state_initialization::mom_initialize_state().

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 
Here is the call graph for this function:
Here is the caller graph for this function: