Skip to content

Commit d6f594e

Browse files
committed
Use reproducing_sum with unscale in 5 places
Use reproducing_sum with unscale to calculate the global mass diagnostic, globally integrated ocean surface area, offline tracer transport residuals and in calculating the OBC inflow area in tidal_bay_set_OBC_data(). This change allows for the elimination or replacement of 6 rescaling factors and one added instance with a consistent conversion factor and diagnostic units occur on the same line. All answers are bitwise identical.
1 parent 1c0d3d1 commit d6f594e

File tree

4 files changed

+21
-31
lines changed

4 files changed

+21
-31
lines changed

src/diagnostics/MOM_diagnostics.F90

+5-4
Original file line numberDiff line numberDiff line change
@@ -332,9 +332,9 @@ subroutine calculate_diagnostic_fields(u, v, h, uh, vh, tv, ADp, CDp, p_surf, &
332332
if (CS%id_masso > 0) then
333333
mass_cell(:,:) = 0.0
334334
do k=1,nz ; do j=js,je ; do i=is,ie
335-
mass_cell(i,j) = mass_cell(i,j) + (GV%H_to_kg_m2*h(i,j,k)) * US%L_to_m**2*G%areaT(i,j)
335+
mass_cell(i,j) = mass_cell(i,j) + (GV%H_to_RZ*h(i,j,k)) * G%areaT(i,j)
336336
enddo ; enddo ; enddo
337-
masso = reproducing_sum(mass_cell)
337+
masso = reproducing_sum(mass_cell, unscale=US%RZL2_to_kg)
338338
call post_data(CS%id_masso, masso, CS%diag)
339339
endif
340340

@@ -1644,8 +1644,9 @@ subroutine MOM_diagnostics_init(MIS, ADp, CDp, Time, G, GV, US, param_file, diag
16441644
Time, 'Mass per unit area of liquid ocean grid cell', 'kg m-2', conversion=GV%H_to_kg_m2, &
16451645
standard_name='sea_water_mass_per_unit_area', v_extensive=.true.)
16461646

1647-
CS%id_masso = register_scalar_field('ocean_model', 'masso', Time, &
1648-
diag, 'Mass of liquid ocean', 'kg', standard_name='sea_water_mass')
1647+
CS%id_masso = register_scalar_field('ocean_model', 'masso', Time, diag, &
1648+
'Mass of liquid ocean', units='kg', conversion=US%RZL2_to_kg, &
1649+
standard_name='sea_water_mass')
16491650

16501651
CS%id_thkcello = register_diag_field('ocean_model', 'thkcello', diag%axesTL, Time, &
16511652
long_name='Cell Thickness', standard_name='cell_thickness', &

src/initialization/MOM_shared_initialization.F90

+3-6
Original file line numberDiff line numberDiff line change
@@ -1314,17 +1314,14 @@ subroutine compute_global_grid_integrals(G, US)
13141314

13151315
! Local variables
13161316
real, dimension(G%isc:G%iec, G%jsc:G%jec) :: tmpForSumming ! Masked and unscaled cell areas [m2]
1317-
real :: area_scale ! A scaling factor for area into MKS units [m2 L-2 ~> 1]
1318-
integer :: i,j
1319-
1320-
area_scale = US%L_to_m**2
1317+
integer :: i, j
13211318

13221319
tmpForSumming(:,:) = 0.
13231320
G%areaT_global = 0.0 ; G%IareaT_global = 0.0
13241321
do j=G%jsc,G%jec ; do i=G%isc,G%iec
1325-
tmpForSumming(i,j) = area_scale*G%areaT(i,j) * G%mask2dT(i,j)
1322+
tmpForSumming(i,j) = G%areaT(i,j) * G%mask2dT(i,j)
13261323
enddo ; enddo
1327-
G%areaT_global = reproducing_sum(tmpForSumming)
1324+
G%areaT_global = US%L_to_m**2 * reproducing_sum(tmpForSumming, unscale=US%L_to_m**2)
13281325

13291326
if (G%areaT_global == 0.0) &
13301327
call MOM_error(FATAL, "compute_global_grid_integrals: zero ocean area (check topography?)")

src/tracer/MOM_offline_main.F90

+10-14
Original file line numberDiff line numberDiff line change
@@ -625,27 +625,24 @@ real function remaining_transport_sum(G, GV, US, uhtr, vhtr, h_new)
625625

626626
! Local variables
627627
real, dimension(SZI_(G),SZJ_(G)) :: trans_rem_col !< The vertical sum of the absolute value of
628-
!! transports through the faces of a column, in MKS units [kg].
628+
!! transports through the faces of a column [R Z L2 ~> kg].
629629
real :: trans_cell !< The sum of the absolute value of the remaining transports through the faces
630630
!! of a tracer cell [H L2 ~> m3 or kg]
631-
real :: HL2_to_kg_scale !< Unit conversion factor to cell mass [kg H-1 L-2 ~> kg m-3 or 1]
632631
integer :: i, j, k, is, ie, js, je, nz
633632

634633
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec ; nz = GV%ke
635634

636-
HL2_to_kg_scale = GV%H_to_kg_m2 * US%L_to_m**2
637-
638635
trans_rem_col(:,:) = 0.0
639636
do k=1,nz ; do j=js,je ; do i=is,ie
640637
trans_cell = (ABS(uhtr(I-1,j,k)) + ABS(uhtr(I,j,k))) + &
641638
(ABS(vhtr(i,J-1,k)) + ABS(vhtr(i,J,k)))
642639
if (trans_cell > max(1.0e-16*h_new(i,j,k), GV%H_subroundoff) * G%areaT(i,j)) &
643-
trans_rem_col(i,j) = trans_rem_col(i,j) + HL2_to_kg_scale * trans_cell
640+
trans_rem_col(i,j) = trans_rem_col(i,j) + GV%H_to_RZ * trans_cell
644641
enddo ; enddo ; enddo
645642

646643
! The factor of 0.5 here is to avoid double-counting because two cells share a face.
647-
remaining_transport_sum = 0.5 * GV%kg_m2_to_H*US%m_to_L**2 * &
648-
reproducing_sum(trans_rem_col, is+(1-G%isd), ie+(1-G%isd), js+(1-G%jsd), je+(1-G%jsd))
644+
remaining_transport_sum = 0.5 * GV%RZ_to_H * reproducing_sum(trans_rem_col, &
645+
is+(1-G%isd), ie+(1-G%isd), js+(1-G%jsd), je+(1-G%jsd), unscale=US%RZL2_to_kg)
649646

650647
end function remaining_transport_sum
651648

@@ -876,8 +873,8 @@ subroutine offline_advection_layer(fluxes, Time_start, time_interval, G, GV, US,
876873
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)) :: vhtr_sub ! Remaining meridional mass transports [H L2 ~> m3 or kg]
877874

878875
real, dimension(SZI_(G),SZJB_(G)) :: rem_col_flux ! The summed absolute value of the remaining
879-
! fluxes through the faces of a column or within a column, in mks units [kg]
880-
real :: sum_flux ! Globally summed absolute value of fluxes in mks units [kg], which is
876+
! mass fluxes through the faces of a column or within a column [R Z L2 ~> kg]
877+
real :: sum_flux ! Globally summed absolute value of fluxes [R Z L2 ~> kg], which is
881878
! used to keep track of how close to convergence we are.
882879

883880
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: &
@@ -890,7 +887,6 @@ subroutine offline_advection_layer(fluxes, Time_start, time_interval, G, GV, US,
890887
! Work arrays for temperature and salinity
891888
integer :: iter
892889
real :: dt_iter ! The timestep of each iteration [T ~> s]
893-
real :: HL2_to_kg_scale ! Unit conversion factors to cell mass [kg H-1 L-2 ~> kg m-3 or 1]
894890
character(len=160) :: mesg ! The text of an error message
895891
integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
896892
integer :: IsdB, IedB, JsdB, JedB
@@ -993,22 +989,22 @@ subroutine offline_advection_layer(fluxes, Time_start, time_interval, G, GV, US,
993989
call pass_vector(uhtr,vhtr,G%Domain)
994990

995991
! Calculate how close we are to converging by summing the remaining fluxes at each point
996-
HL2_to_kg_scale = US%L_to_m**2*GV%H_to_kg_m2
997992
rem_col_flux(:,:) = 0.0
998993
do k=1,nz ; do j=js,je ; do i=is,ie
999-
rem_col_flux(i,j) = rem_col_flux(i,j) + HL2_to_kg_scale * &
994+
rem_col_flux(i,j) = rem_col_flux(i,j) + GV%H_to_RZ * &
1000995
( (abs(eatr(i,j,k)) + abs(ebtr(i,j,k))) + &
1001996
((abs(uhtr(I-1,j,k)) + abs(uhtr(I,j,k))) + &
1002997
(abs(vhtr(i,J-1,k)) + abs(vhtr(i,J,k))) ) )
1003998
enddo ; enddo ; enddo
1004-
sum_flux = reproducing_sum(rem_col_flux, is+(1-G%isd), ie+(1-G%isd), js+(1-G%jsd), je+(1-G%jsd))
999+
sum_flux = reproducing_sum(rem_col_flux, is+(1-G%isd), ie+(1-G%isd), js+(1-G%jsd), je+(1-G%jsd), &
1000+
unscale=US%RZL2_to_kg)
10051001

10061002
if (sum_flux==0) then
10071003
write(mesg,*) 'offline_advection_layer: Converged after iteration', iter
10081004
call MOM_mesg(mesg)
10091005
exit
10101006
else
1011-
write(mesg,*) "offline_advection_layer: Iteration ", iter, " remaining total fluxes: ", sum_flux
1007+
write(mesg,*) "offline_advection_layer: Iteration ", iter, " remaining total fluxes: ", sum_flux*US%RZL2_to_kg
10121008
call MOM_mesg(mesg)
10131009
endif
10141010

src/user/tidal_bay_initialization.F90

+3-7
Original file line numberDiff line numberDiff line change
@@ -78,8 +78,7 @@ subroutine tidal_bay_set_OBC_data(OBC, CS, G, GV, US, h, Time)
7878
real :: my_flux ! The vlume flux through the face [L2 Z T-1 ~> m3 s-1]
7979
real :: total_area ! The total face area of the OBCs [L Z ~> m2]
8080
real :: PI ! The ratio of the circumference of a circle to its diameter [nondim]
81-
real :: flux_scale ! A scaling factor for the areas [m2 H-1 L-1 ~> nondim or m3 kg-1]
82-
real, allocatable :: my_area(:,:) ! The total OBC inflow area [m2]
81+
real, allocatable :: my_area(:,:) ! The total OBC inflow area [L Z ~> m2]
8382
integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, n
8483
integer :: IsdB, IedB, JsdB, JedB
8584
type(OBC_segment_type), pointer :: segment => NULL()
@@ -94,8 +93,6 @@ subroutine tidal_bay_set_OBC_data(OBC, CS, G, GV, US, h, Time)
9493

9594
allocate(my_area(1:1,js:je))
9695

97-
flux_scale = GV%H_to_m*US%L_to_m
98-
9996
time_sec = US%s_to_T*time_type_to_real(Time)
10097
cff_eta = CS%tide_ssh_amp * sin(2.0*PI*time_sec / CS%tide_period)
10198
my_area = 0.0
@@ -105,12 +102,11 @@ subroutine tidal_bay_set_OBC_data(OBC, CS, G, GV, US, h, Time)
105102
do j=segment%HI%jsc,segment%HI%jec ; do I=segment%HI%IscB,segment%HI%IecB
106103
if (OBC%segnum_u(I,j) /= OBC_NONE) then
107104
do k=1,nz
108-
! This area has to be in MKS units to work with reproducing_sum.
109-
my_area(1,j) = my_area(1,j) + h(I,j,k)*flux_scale*G%dyCu(I,j)
105+
my_area(1,j) = my_area(1,j) + h(I,j,k)*(GV%H_to_m*US%m_to_Z)*G%dyCu(I,j)
110106
enddo
111107
endif
112108
enddo ; enddo
113-
total_area = US%m_to_Z*US%m_to_L * reproducing_sum(my_area)
109+
total_area = reproducing_sum(my_area, unscale=US%Z_to_m*US%L_to_m)
114110
my_flux = - CS%tide_flow * SIN(2.0*PI*time_sec / CS%tide_period)
115111

116112
do n = 1, OBC%number_of_segments

0 commit comments

Comments
 (0)