Skip to content

Commit

Permalink
*Altered initial conditions for circle_obcs
Browse files Browse the repository at this point in the history
- The circle_obcs experiment ostensibly tests the Flather open
  boundary conditions except it was modified to use 10 layers
  and also test the Orlanski boundary conditions too. The initialization
  code appears to have been written for a 2-layer system and
  essentially only forces a very high-mode. This change produces a
  mixed mode initial condition (thickness anomaly linear with z) but
  with most energy in the first baroclinic mode.
- Changes circle_obcs experiment only.
  • Loading branch information
adcroft committed Sep 13, 2016
1 parent 2577c8a commit 36275ae
Showing 1 changed file with 10 additions and 2 deletions.
12 changes: 10 additions & 2 deletions src/user/circle_obcs_initialization.F90
Original file line number Diff line number Diff line change
Expand Up @@ -93,8 +93,16 @@ subroutine circle_obcs_initialize_thickness(h, G, GV, param_file)
! if (rad <= 6.*diskrad) h(i,j,k) = h(i,j,k)+10.0*exp( -0.5*( rad**2 ) )
rad = min( rad, 1. ) ! Flatten outside radius of diskrad
rad = rad*(2.*asin(1.)) ! Map 0-1 to 0-pi
h(i,j,k) = h(i,j,k) + 10.0*0.5*(1.+cos(rad)) ! cosine bell
h(i,j,k-1) = h(i,j,k-1) - 10.0*0.5*(1.+cos(rad)) ! cosine bell
if (Nz==1) then
! The model is barotropic
h(i,j,k) = h(i,j,k) + 1.0*0.5*(1.+cos(rad)) ! cosine bell
else
! The model is baroclinic
do k = 1, Nz
h(i,j,k) = h(i,j,k) - 0.5*(1.+cos(rad)) & ! cosine bell
* 5.0 * real( 2*k-Nz )
enddo
endif
enddo ; enddo

end subroutine circle_obcs_initialize_thickness
Expand Down

0 comments on commit 36275ae

Please sign in to comment.