Skip to content

Commit

Permalink
+Added a run-time parameter for circle_obcs
Browse files Browse the repository at this point in the history
  Added a new runtime parameter, DISK_IC_AMPLITUDE, that is used in setting up
the circle_obcs test case.  By default all answers are bitwise identical, but
the MOM_parameter_doc.short files change for the benchmark test case.
  • Loading branch information
Hallberg-NOAA committed Nov 6, 2018
1 parent fa1a762 commit 6ccb6b6
Showing 1 changed file with 9 additions and 5 deletions.
14 changes: 9 additions & 5 deletions src/user/circle_obcs_initialization.F90
Original file line number Diff line number Diff line change
Expand Up @@ -37,10 +37,11 @@ subroutine circle_obcs_initialize_thickness(h, G, GV, param_file, just_read_para
! negative because it is positive upward.
real :: eta1D(SZK_(GV)+1)! Interface height relative to the sea surface
! positive upward, in in depth units (Z).
real :: IC_amp ! The amplitude of the initial height displacement, in H.
real :: diskrad, rad, xCenter, xRadius, lonC, latC, xOffset
logical :: just_read
! This include declares and sets the variable "version".
#include "version_variable.h"
! This include declares and sets the variable "version".
# include "version_variable.h"
character(len=40) :: mdl = "circle_obcs_initialization" ! This module's name.
integer :: i, j, k, is, ie, js, je, nz

Expand All @@ -61,6 +62,10 @@ subroutine circle_obcs_initialize_thickness(h, G, GV, param_file, just_read_para
"The x-offset of the initially elevated disk in the \n"//&
"circle_obcs test case.", units=G%x_axis_units, &
default = 0.0, do_not_log=just_read)
call get_param(param_file, mdl, "DISK_IC_AMPLITUDE", IC_amp, &
"Initial amplitude of interface height displacements \n"//&
"in the circle_obcs test case.", &
units='m', default=5.0, scale=GV%m_to_H, do_not_log=just_read)

if (just_read) return ! All run-time parameters have been read, so return.

Expand Down Expand Up @@ -93,12 +98,11 @@ subroutine circle_obcs_initialize_thickness(h, G, GV, param_file, just_read_para
rad = rad*(2.*asin(1.)) ! Map 0-1 to 0-pi
if (nz==1) then
! The model is barotropic
h(i,j,k) = h(i,j,k) + GV%m_to_H * 1.0*0.5*(1.+cos(rad)) ! cosine bell
h(i,j,k) = h(i,j,k) + IC_amp * 0.5*(1.+cos(rad)) ! cosine bell
else
! The model is baroclinic
do k = 1, nz
h(i,j,k) = h(i,j,k) - GV%m_to_H * 0.5*(1.+cos(rad)) & ! cosine bell
* 5.0 * real( 2*k-nz )
h(i,j,k) = h(i,j,k) - 0.5*(1.+cos(rad)) * IC_amp * real( 2*k-nz )
enddo
endif
enddo ; enddo
Expand Down

0 comments on commit 6ccb6b6

Please sign in to comment.