Skip to content

Commit f835394

Browse files
committed
Calculate volo in scaled units
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.
1 parent d6f3fa0 commit f835394

File tree

1 file changed

+5
-5
lines changed

1 file changed

+5
-5
lines changed

src/diagnostics/MOM_diagnostics.F90

Lines changed: 5 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -204,7 +204,7 @@ subroutine calculate_diagnostic_fields(u, v, h, uh, vh, tv, ADp, CDp, p_surf, &
204204
! including [nondim] and [H ~> m or kg m-2].
205205
real :: uh_tmp(SZIB_(G),SZJ_(G),SZK_(GV)) ! A temporary zonal transport [H L2 T-1 ~> m3 s-1 or kg s-1]
206206
real :: vh_tmp(SZI_(G),SZJB_(G),SZK_(GV)) ! A temporary meridional transport [H L2 T-1 ~> m3 s-1 or kg s-1]
207-
real :: mass_cell(SZI_(G),SZJ_(G)) ! The vertically integrated mass in a grid cell [kg]
207+
real :: mass_cell(SZI_(G),SZJ_(G)) ! The vertically integrated mass in a grid cell [R Z L2 ~> kg]
208208
real :: rho_in_situ(SZI_(G)) ! In situ density [R ~> kg m-3]
209209
real :: cg1(SZI_(G),SZJ_(G)) ! First baroclinic gravity wave speed [L T-1 ~> m s-1]
210210
real :: Rd1(SZI_(G),SZJ_(G)) ! First baroclinic deformation radius [L ~> m]
@@ -226,7 +226,7 @@ subroutine calculate_diagnostic_fields(u, v, h, uh, vh, tv, ADp, CDp, p_surf, &
226226
real, dimension(SZK_(GV)) :: salt_layer_ave ! The average salinity in a layer [S ~> ppt]
227227
real :: thetaoga ! The volume mean potential temperature [C ~> degC]
228228
real :: soga ! The volume mean ocean salinity [S ~> ppt]
229-
real :: masso ! The total mass of the ocean [kg]
229+
real :: masso ! The total mass of the ocean [R Z L2 ~> kg]
230230
real :: tosga ! The area mean sea surface temperature [C ~> degC]
231231
real :: sosga ! The area mean sea surface salinity [S ~> ppt]
232232

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

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

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

0 commit comments

Comments
 (0)