MOM6
circle_obcs_initialization Module Reference

Detailed Description

The module configures the model for the "circle_obcs" experiment. circle_obcs = Test of Open Boundary Conditions for an SSH anomaly.

Functions/Subroutines

subroutine, public circle_obcs_initialize_thickness (h, G, GV, param_file, just_read_params)
 This subroutine initializes layer thicknesses for the circle_obcs experiment. More...
 

Function/Subroutine Documentation

◆ circle_obcs_initialize_thickness()

subroutine, public circle_obcs_initialization::circle_obcs_initialize_thickness ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(out)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
logical, intent(in), optional  just_read_params 
)

This subroutine initializes layer thicknesses for the circle_obcs experiment.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical 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 42 of file circle_obcs_initialization.F90.

References mom_error_handler::mom_mesg().

Referenced by mom_state_initialization::mom_initialize_state().

42  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
43  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
44  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)), &
45  intent(out) :: h !< The thickness that is being initialized, in m.
46  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
47  !! to parse for model parameter values.
48  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
49  !! only read parameters without changing h.
50 
51  real :: e0(szk_(gv)+1) ! The resting interface heights, in m, usually !
52  ! negative because it is positive upward. !
53  real :: eta1d(szk_(gv)+1)! Interface height relative to the sea surface !
54  ! positive upward, in m. !
55  real :: diskrad, rad, xcenter, xradius, lonc, latc
56  logical :: just_read
57 ! This include declares and sets the variable "version".
58 #include "version_variable.h"
59  character(len=40) :: mdl = "circle_obcs_initialization" ! This module's name.
60  integer :: i, j, k, is, ie, js, je, nz
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(" circle_obcs_initialization.F90, circle_obcs_initialize_thickness: setting thickness", 5)
68 
69  if (.not.just_read) call log_version(param_file, mdl, version, "")
70  ! Parameters read by cartesian grid initialization
71  call get_param(param_file, mdl, "DISK_RADIUS", diskrad, &
72  "The radius of the initially elevated disk in the \n"//&
73  "circle_obcs test case.", units=g%x_axis_units, &
74  fail_if_missing=.not.just_read, do_not_log=just_read)
75 
76  if (just_read) return ! All run-time parameters have been read, so return.
77 
78  do k=1,nz
79  e0(k) = -g%max_depth * real(k-1) / real(nz)
80  enddo
81 
82  ! Uniform thicknesses for base state
83  do j=js,je ; do i=is,ie !
84  eta1d(nz+1) = -1.0*g%bathyT(i,j)
85  do k=nz,1,-1
86  eta1d(k) = e0(k)
87  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_z)) then
88  eta1d(k) = eta1d(k+1) + gv%Angstrom_z
89  h(i,j,k) = gv%Angstrom_z
90  else
91  h(i,j,k) = eta1d(k) - eta1d(k+1)
92  endif
93  enddo
94  enddo ; enddo
95 
96  ! Perturb base state by circular anomaly in center
97  k=nz
98  latc = g%south_lat + 0.5*g%len_lat
99  lonc = g%west_lon + 0.5*g%len_lon
100  do j=js,je ; do i=is,ie
101  rad = sqrt((g%geoLonT(i,j)-lonc)**2+(g%geoLatT(i,j)-latc)**2)/(diskrad)
102  ! if (rad <= 6.*diskrad) h(i,j,k) = h(i,j,k)+10.0*exp( -0.5*( rad**2 ) )
103  rad = min( rad, 1. ) ! Flatten outside radius of diskrad
104  rad = rad*(2.*asin(1.)) ! Map 0-1 to 0-pi
105  if (nz==1) then
106  ! The model is barotropic
107  h(i,j,k) = h(i,j,k) + 1.0*0.5*(1.+cos(rad)) ! cosine bell
108  else
109  ! The model is baroclinic
110  do k = 1, nz
111  h(i,j,k) = h(i,j,k) - 0.5*(1.+cos(rad)) & ! cosine bell
112  * 5.0 * real( 2*k-nz )
113  enddo
114  endif
115  enddo ; enddo
116 
Here is the call graph for this function:
Here is the caller graph for this function: