Skip to content

Commit

Permalink
Calculate volo in scaled units
Browse files Browse the repository at this point in the history
  Calculate the volo diagnostic in scaled units using the tmp_scale argument to
global_area_integral and undo this scaling with a conversion argument on its
register_scalar_field call.  Also corrected the units in the comments describing
the mass_celland masso variables in calculate_diagnostic_fields.  All answers
are bitwise identical.
  • Loading branch information
Hallberg-NOAA committed Jan 31, 2025
1 parent d6f3fa0 commit f835394
Showing 1 changed file with 5 additions and 5 deletions.
10 changes: 5 additions & 5 deletions src/diagnostics/MOM_diagnostics.F90
Original file line number Diff line number Diff line change
Expand Up @@ -204,7 +204,7 @@ subroutine calculate_diagnostic_fields(u, v, h, uh, vh, tv, ADp, CDp, p_surf, &
! including [nondim] and [H ~> m or kg m-2].
real :: uh_tmp(SZIB_(G),SZJ_(G),SZK_(GV)) ! A temporary zonal transport [H L2 T-1 ~> m3 s-1 or kg s-1]
real :: vh_tmp(SZI_(G),SZJB_(G),SZK_(GV)) ! A temporary meridional transport [H L2 T-1 ~> m3 s-1 or kg s-1]
real :: mass_cell(SZI_(G),SZJ_(G)) ! The vertically integrated mass in a grid cell [kg]
real :: mass_cell(SZI_(G),SZJ_(G)) ! The vertically integrated mass in a grid cell [R Z L2 ~> kg]
real :: rho_in_situ(SZI_(G)) ! In situ density [R ~> kg m-3]
real :: cg1(SZI_(G),SZJ_(G)) ! First baroclinic gravity wave speed [L T-1 ~> m s-1]
real :: Rd1(SZI_(G),SZJ_(G)) ! First baroclinic deformation radius [L ~> m]
Expand All @@ -226,7 +226,7 @@ subroutine calculate_diagnostic_fields(u, v, h, uh, vh, tv, ADp, CDp, p_surf, &
real, dimension(SZK_(GV)) :: salt_layer_ave ! The average salinity in a layer [S ~> ppt]
real :: thetaoga ! The volume mean potential temperature [C ~> degC]
real :: soga ! The volume mean ocean salinity [S ~> ppt]
real :: masso ! The total mass of the ocean [kg]
real :: masso ! The total mass of the ocean [R Z L2 ~> kg]
real :: tosga ! The area mean sea surface temperature [C ~> degC]
real :: sosga ! The area mean sea surface salinity [S ~> ppt]

Expand Down Expand Up @@ -1341,7 +1341,7 @@ subroutine post_surface_thermo_diags(IDs, G, GV, US, diag, dt_int, sfc_state, tv
zos ! dynamic sea lev (zero area mean) from inverse-barometer adjusted ssh [Z ~> m]
real :: I_time_int ! The inverse of the time interval [T-1 ~> s-1].
real :: zos_area_mean ! Global area mean sea surface height [Z ~> m]
real :: volo ! Total volume of the ocean [m3]
real :: volo ! Total volume of the ocean [Z L2 ~> m3]
real :: ssh_ga ! Global ocean area weighted mean sea seaface height [Z ~> m]
integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state
integer :: i, j, is, ie, js, je
Expand Down Expand Up @@ -1375,7 +1375,7 @@ subroutine post_surface_thermo_diags(IDs, G, GV, US, diag, dt_int, sfc_state, tv
do j=js,je ; do i=is,ie
work_2d(i,j) = G%mask2dT(i,j) * (ssh(i,j) + G%bathyT(i,j))
enddo ; enddo
volo = global_area_integral(work_2d, G, unscale=US%Z_to_m)
volo = global_area_integral(work_2d, G, tmp_scale=US%Z_to_m)
call post_data(IDs%id_volo, volo, diag)
endif

Expand Down Expand Up @@ -1914,7 +1914,7 @@ subroutine register_surface_diags(Time, G, US, IDs, diag, tv)

! Vertically integrated, budget, and surface state diagnostics
IDs%id_volo = register_scalar_field('ocean_model', 'volo', Time, diag, &
long_name='Total volume of liquid ocean', units='m3', &
long_name='Total volume of liquid ocean', units='m3', conversion=US%Z_to_m*US%L_to_m**2, &
standard_name='sea_water_volume')
IDs%id_zos = register_diag_field('ocean_model', 'zos', diag%axesT1, Time, &
standard_name = 'sea_surface_height_above_geoid', &
Expand Down

0 comments on commit f835394

Please sign in to comment.