MOM6
circle_obcs_initialization.F90
1 !> Configures the model for the "circle_obcs" experiment which tests
2 !! Open Boundary Conditions radiating an SSH anomaly.
4 
5 ! This file is part of MOM6. See LICENSE.md for the license.
6 
7 use mom_sponge, only : sponge_cs, set_up_sponge_field, initialize_sponge
8 use mom_error_handler, only : mom_mesg, mom_error, fatal, is_root_pe
10 use mom_get_input, only : directories
11 use mom_grid, only : ocean_grid_type
16 
17 implicit none ; private
18 
19 #include <MOM_memory.h>
20 
21 public circle_obcs_initialize_thickness
22 
23 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
24 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
25 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
26 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
27 
28 contains
29 
30 !> This subroutine initializes layer thicknesses for the circle_obcs experiment.
31 subroutine circle_obcs_initialize_thickness(h, G, GV, param_file, just_read_params)
32  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
33  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
34  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)), &
35  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
36  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
37  !! to parse for model parameter values.
38  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
39  !! only read parameters without changing h.
40 
41  real :: e0(szk_(gv)+1) ! The resting interface heights, in depth units [Z ~> m], usually
42  ! negative because it is positive upward.
43  real :: eta1D(szk_(gv)+1)! Interface height relative to the sea surface
44  ! positive upward, in depth units [Z ~> m].
45  real :: IC_amp ! The amplitude of the initial height displacement [H ~> m or kg m-2].
46  real :: diskrad, rad, xCenter, xRadius, lonC, latC, xOffset
47  logical :: just_read
48  ! This include declares and sets the variable "version".
49 # include "version_variable.h"
50  character(len=40) :: mdl = "circle_obcs_initialization" ! This module's name.
51  integer :: i, j, k, is, ie, js, je, nz
52 
53  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
54 
55  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
56 
57  if (.not.just_read) &
58  call mom_mesg(" circle_obcs_initialization.F90, circle_obcs_initialize_thickness: setting thickness", 5)
59 
60  if (.not.just_read) call log_version(param_file, mdl, version, "")
61  ! Parameters read by cartesian grid initialization
62  call get_param(param_file, mdl, "DISK_RADIUS", diskrad, &
63  "The radius of the initially elevated disk in the "//&
64  "circle_obcs test case.", units=g%x_axis_units, &
65  fail_if_missing=.not.just_read, do_not_log=just_read)
66  call get_param(param_file, mdl, "DISK_X_OFFSET", xoffset, &
67  "The x-offset of the initially elevated disk in the "//&
68  "circle_obcs test case.", units=g%x_axis_units, &
69  default = 0.0, do_not_log=just_read)
70  call get_param(param_file, mdl, "DISK_IC_AMPLITUDE", ic_amp, &
71  "Initial amplitude of interface height displacements "//&
72  "in the circle_obcs test case.", &
73  units='m', default=5.0, scale=gv%m_to_H, do_not_log=just_read)
74 
75  if (just_read) return ! All run-time parameters have been read, so return.
76 
77  do k=1,nz
78  e0(k) = -g%max_depth * real(k-1) / real(nz)
79  enddo
80 
81  ! Uniform thicknesses for base state
82  do j=js,je ; do i=is,ie !
83  eta1d(nz+1) = -g%bathyT(i,j)
84  do k=nz,1,-1
85  eta1d(k) = e0(k)
86  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z)) then
87  eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
88  h(i,j,k) = gv%Angstrom_H
89  else
90  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
91  endif
92  enddo
93  enddo ; enddo
94 
95  ! Perturb base state by circular anomaly in center
96  k=nz
97  latc = g%south_lat + 0.5*g%len_lat
98  lonc = g%west_lon + 0.5*g%len_lon + xoffset
99  do j=js,je ; do i=is,ie
100  rad = sqrt((g%geoLonT(i,j)-lonc)**2+(g%geoLatT(i,j)-latc)**2)/(diskrad)
101  ! if (rad <= 6.*diskrad) h(i,j,k) = h(i,j,k)+10.0*exp( -0.5*( rad**2 ) )
102  rad = min( rad, 1. ) ! Flatten outside radius of diskrad
103  rad = rad*(2.*asin(1.)) ! Map 0-1 to 0-pi
104  if (nz==1) then
105  ! The model is barotropic
106  h(i,j,k) = h(i,j,k) + ic_amp * 0.5*(1.+cos(rad)) ! cosine bell
107  else
108  ! The model is baroclinic
109  do k = 1, nz
110  h(i,j,k) = h(i,j,k) - 0.5*(1.+cos(rad)) * ic_amp * real( 2*k-nz )
111  enddo
112  endif
113  enddo ; enddo
114 
115 end subroutine circle_obcs_initialize_thickness
116 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
Reads the only Fortran name list needed to boot-strap the model.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
The MOM6 facility to parse input files for runtime parameters.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Type to carry basic tracer information.
Routines for error handling and I/O management.
This control structure holds memory and parameters for the MOM_sponge module.
Definition: MOM_sponge.F90:41
Implements sponge regions in isopycnal mode.
Definition: MOM_sponge.F90:2
Configures the model for the "circle_obcs" experiment which tests Open Boundary Conditions radiating ...
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Definition: MOM_EOS.F90:81
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Container for paths and parameter file names.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108