Skip to content

Commit

Permalink
(*)Fixes gprime(1) when no layer coordinates are set
Browse files Browse the repository at this point in the history
- The "none" option for COORD_CONFIG literally did nothing which
  meant that Rlay and gprime were unset (they actually appeared to
  be zero).  This mode was added with the intent that the model
  should not need these coordinate variables in ALE mode. However,
  gprime(k=1) is used and should be set to GFS, which this commit
  now does.
- I also added "ALE" as an option for COORD_CONFIG which is
  equivalent to "none" since "ALE" is bit more intuitive.
  - We might even use this as a master switch or check that
    COORD_CONFIG is consistent with USE_REGRIDDING (for later)?
- This commit "appears" to changes answers but I don't think it
  does. What changes is the calculation of energy for some
  experiments in which Rlay was not set before this commit.
  A non-zero gprime(1) produces different APEs.
  Experiments affected are:
    - CVmix_SCM_tests/cooling_only/EPBL
    - CVmix_SCM_tests/cooling_only/KPP
    - CVmix_SCM_tests/mech_only/EPBL
    - CVmix_SCM_tests/mech_only/KPP
    - CVmix_SCM_tests/skin_warming_wind/EPBL
    - CVmix_SCM_tests/skin_warming_wind/KPP
    - CVmix_SCM_tests/wind_only/EPBL
    - CVmix_SCM_tests/wind_only/KPP
    - SCM_idealized_hurricane
    - single_column/EPBL
    - single_column/KPP
    - unit_tests
  All the above have COORD_CONFIG="none".
  • Loading branch information
adcroft committed Jul 28, 2017
1 parent fccb7ae commit be46fee
Showing 1 changed file with 35 additions and 2 deletions.
37 changes: 35 additions & 2 deletions src/initialization/MOM_coord_initialization.F90
Original file line number Diff line number Diff line change
Expand Up @@ -58,6 +58,7 @@ subroutine MOM_initialize_coord(GV, PF, write_geom, output_dir, tv, max_depth)
! Set-up the layer densities, GV%Rlay, and reduced gravities, GV%g_prime.
call get_param(PF, mod, "COORD_CONFIG", config, &
"This specifies how layers are to be defined: \n"//&
" \t ALE or none - used to avoid defining layers in ALE mode \n"//&
" \t file - read coordinate information from the file \n"//&
" \t\t specified by (COORD_FILE).\n"//&
" \t BFB - Custom coords for buoyancy-forced basin case \n"//&
Expand Down Expand Up @@ -93,7 +94,8 @@ subroutine MOM_initialize_coord(GV, PF, write_geom, output_dir, tv, max_depth)
call user_set_coord(GV%Rlay, GV%g_prime, GV, PF, eos)
case ("BFB")
call BFB_set_coord(GV%Rlay, GV%g_prime, GV, PF, eos)
case ("none")
case ("none", "ALE")
call set_coord_to_none(GV%Rlay, GV%g_prime, GV, PF)
case default ; call MOM_error(FATAL,"MOM_initialize_coord: "// &
"Unrecognized coordinate setup"//trim(config))
end select
Expand Down Expand Up @@ -493,7 +495,38 @@ subroutine set_coord_linear(Rlay, g_prime, GV, param_file)
call callTree_leave(trim(mod)//'()')
end subroutine set_coord_linear

! -----------------------------------------------------------------------------
!> Sets Rlay to Rho0 and g_prime to zero except for the free surface.
!! This is for use only in ALE mode where Rlay should not be used and g_prime(1) alone
!! might be used.
subroutine set_coord_to_none(Rlay, g_prime, GV, param_file)
real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
!! (potential density).
real, dimension(:), intent(out) :: g_prime !< A structure indicating the open file to
!! parse for model parameter values.
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
real :: g_fs ! Reduced gravity across the free surface, in m s-2.
character(len=40) :: mdl = "set_coord_to_none" ! This subroutine's name.
integer :: k, nz
nz = GV%ke

call callTree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")

call get_param(param_file, mdl, "GFS" , g_fs, &
"The reduced gravity at the free surface.", units="m s-2", &
default=GV%g_Earth)

g_prime(1) = g_fs
do k=2,nz ; g_prime(k) = 0. ; enddo
Rlay(1) = GV%Rho0
do k=2,nz ; Rlay(k) = Rlay(k-1) + g_prime(k)*(GV%Rho0/GV%g_Earth) ; enddo

call callTree_leave(trim(mdl)//'()')

end subroutine set_coord_to_none

!> This subroutine writes out a file containing any available data related
!! to the vertical grid used by the MOM ocean model.
subroutine write_vertgrid_file(GV, param_file, directory)
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
Expand Down

0 comments on commit be46fee

Please sign in to comment.