MOM6
mom_spatial_means Module Reference

Detailed Description

Functions and routines to take area, volume, mass-weighted, layerwise, zonal or meridional means.

Functions/Subroutines

real function, public global_area_mean (var, G, scale)
 Return the global area mean of a variable. This uses reproducing sums. More...
 
real function, public global_area_integral (var, G, scale, area)
 Return the global area integral of a variable, by default using the masked area from the grid, but an alternate could be used instead. This uses reproducing sums. More...
 
real function, dimension(gv %ke), public global_layer_mean (var, h, G, GV, scale)
 Return the layerwise global thickness-weighted mean of a variable. This uses reproducing sums. More...
 
real function, public global_volume_mean (var, h, G, GV, scale)
 Find the global thickness-weighted mean of a variable. This uses reproducing sums. More...
 
real function, public global_mass_integral (h, G, GV, var, on_PE_only, scale)
 Find the global mass-weighted integral of a variable. This uses reproducing sums. More...
 
subroutine, public global_i_mean (array, i_mean, G, mask, scale, tmp_scale)
 Determine the global mean of a field along rows of constant i, returning it in a 1-d array using the local indexing. This uses reproducing sums. More...
 
subroutine, public global_j_mean (array, j_mean, G, mask, scale, tmp_scale)
 Determine the global mean of a field along rows of constant j, returning it in a 1-d array using the local indexing. This uses reproducing sums. More...
 
subroutine, public adjust_area_mean_to_zero (array, G, scaling, unit_scale)
 Adjust 2d array such that area mean is zero without moving the zero contour. More...
 

Function/Subroutine Documentation

◆ adjust_area_mean_to_zero()

subroutine, public mom_spatial_means::adjust_area_mean_to_zero ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(inout)  array,
type(ocean_grid_type), intent(in)  G,
real, intent(out), optional  scaling,
real, intent(in), optional  unit_scale 
)

Adjust 2d array such that area mean is zero without moving the zero contour.

Parameters
[in]gGrid structure
[in,out]array2D array to be adjusted
[out]scalingThe scaling factor used
[in]unit_scaleA rescaling factor for the variable

Definition at line 367 of file MOM_spatial_means.F90.

368  type(ocean_grid_type), intent(in) :: G !< Grid structure
369  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: array !< 2D array to be adjusted
370  real, optional, intent(out) :: scaling !< The scaling factor used
371  real, optional, intent(in) :: unit_scale !< A rescaling factor for the variable
372  ! Local variables
373  real, dimension(G%isc:G%iec, G%jsc:G%jec) :: posVals, negVals, areaXposVals, areaXnegVals
374  integer :: i,j
375  type(EFP_type), dimension(2) :: areaInt_EFP
376  real :: scalefac ! A scaling factor for the variable.
377  real :: I_scalefac ! The Adcroft reciprocal of scalefac
378  real :: areaIntPosVals, areaIntNegVals, posScale, negScale
379 
380  scalefac = 1.0 ; if (present(unit_scale)) scalefac = unit_scale
381  i_scalefac = 0.0 ; if (scalefac /= 0.0) i_scalefac = 1.0 / scalefac
382 
383  ! areaXposVals(:,:) = 0. ! This zeros out halo points.
384  ! areaXnegVals(:,:) = 0. ! This zeros out halo points.
385 
386  do j=g%jsc,g%jec ; do i=g%isc,g%iec
387  posvals(i,j) = max(0., scalefac*array(i,j))
388  areaxposvals(i,j) = g%US%L_to_m**2*g%areaT(i,j) * posvals(i,j)
389  negvals(i,j) = min(0., scalefac*array(i,j))
390  areaxnegvals(i,j) = g%US%L_to_m**2*g%areaT(i,j) * negvals(i,j)
391  enddo ; enddo
392 
393  ! Combining the sums like this avoids separate blocking global sums.
394  areaint_efp(1) = reproducing_sum_efp( areaxposvals, only_on_pe=.true. )
395  areaint_efp(2) = reproducing_sum_efp( areaxnegvals, only_on_pe=.true. )
396  call efp_sum_across_pes(areaint_efp, 2)
397  areaintposvals = efp_to_real( areaint_efp(1) )
398  areaintnegvals = efp_to_real( areaint_efp(2) )
399 
400  posscale = 0.0 ; negscale = 0.0
401  if ((areaintposvals>0.).and.(areaintnegvals<0.)) then ! Only adjust if possible
402  if (areaintposvals>-areaintnegvals) then ! Scale down positive values
403  posscale = - areaintnegvals / areaintposvals
404  do j=g%jsc,g%jec ; do i=g%isc,g%iec
405  array(i,j) = ((posscale * posvals(i,j)) + negvals(i,j)) * i_scalefac
406  enddo ; enddo
407  elseif (areaintposvals<-areaintnegvals) then ! Scale down negative values
408  negscale = - areaintposvals / areaintnegvals
409  do j=g%jsc,g%jec ; do i=g%isc,g%iec
410  array(i,j) = (posvals(i,j) + (negscale * negvals(i,j))) * i_scalefac
411  enddo ; enddo
412  endif
413  endif
414  if (present(scaling)) scaling = posscale - negscale
415 

◆ global_area_integral()

real function, public mom_spatial_means::global_area_integral ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  var,
type(ocean_grid_type), intent(in)  G,
real, intent(in), optional  scale,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in), optional  area 
)

Return the global area integral of a variable, by default using the masked area from the grid, but an alternate could be used instead. This uses reproducing sums.

Parameters
[in]gThe ocean's grid structure
[in]varThe variable to integrate
[in]scaleA rescaling factor for the variable
[in]areaThe alternate area to use, including any required masking [L2 ~> m2].
Returns
The returned area integral, usually in the units of var times [m2].

Definition at line 52 of file MOM_spatial_means.F90.

53  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
54  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: var !< The variable to integrate
55  real, optional, intent(in) :: scale !< A rescaling factor for the variable
56  real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: area !< The alternate area to use, including
57  !! any required masking [L2 ~> m2].
58  real :: global_area_integral !< The returned area integral, usually in the units of var times [m2].
59 
60  ! Local variables
61  real, dimension(SZI_(G),SZJ_(G)) :: tmpForSumming
62  real :: scalefac ! An overall scaling factor for the areas and variable.
63  integer :: i, j, is, ie, js, je
64  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
65 
66  scalefac = g%US%L_to_m**2 ; if (present(scale)) scalefac = g%US%L_to_m**2*scale
67 
68  tmpforsumming(:,:) = 0.
69  if (present(area)) then
70  do j=js,je ; do i=is,ie
71  tmpforsumming(i,j) = var(i,j) * (scalefac * area(i,j))
72  enddo ; enddo
73  else
74  do j=js,je ; do i=is,ie
75  tmpforsumming(i,j) = var(i,j) * (scalefac * g%areaT(i,j) * g%mask2dT(i,j))
76  enddo ; enddo
77  endif
78  global_area_integral = reproducing_sum(tmpforsumming)
79 

◆ global_area_mean()

real function, public mom_spatial_means::global_area_mean ( real, dimension(szi_(g), szj_(g)), intent(in)  var,
type(ocean_grid_type), intent(in)  G,
real, intent(in), optional  scale 
)

Return the global area mean of a variable. This uses reproducing sums.

Parameters
[in]gThe ocean's grid structure
[in]varThe variable to average
[in]scaleA rescaling factor for the variable

Definition at line 28 of file MOM_spatial_means.F90.

29  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
30  real, dimension(SZI_(G), SZJ_(G)), intent(in) :: var !< The variable to average
31  real, optional, intent(in) :: scale !< A rescaling factor for the variable
32 
33  real, dimension(SZI_(G), SZJ_(G)) :: tmpForSumming
34  real :: global_area_mean
35 
36  real :: scalefac ! An overall scaling factor for the areas and variable.
37  integer :: i, j, is, ie, js, je
38  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
39 
40  scalefac = g%US%L_to_m**2 ; if (present(scale)) scalefac = g%US%L_to_m**2*scale
41 
42  tmpforsumming(:,:) = 0.
43  do j=js,je ; do i=is,ie
44  tmpforsumming(i,j) = var(i,j) * (scalefac * g%areaT(i,j) * g%mask2dT(i,j))
45  enddo ; enddo
46  global_area_mean = reproducing_sum(tmpforsumming) * g%IareaT_global
47 

◆ global_i_mean()

subroutine, public mom_spatial_means::global_i_mean ( real, dimension(szi_(g),szj_(g)), intent(in)  array,
real, dimension(szj_(g)), intent(out)  i_mean,
type(ocean_grid_type), intent(inout)  G,
real, dimension(szi_(g),szj_(g)), intent(in), optional  mask,
real, intent(in), optional  scale,
real, intent(in), optional  tmp_scale 
)

Determine the global mean of a field along rows of constant i, returning it in a 1-d array using the local indexing. This uses reproducing sums.

Parameters
[in,out]gThe ocean's grid structure
[in]arrayThe variable being averaged
[out]i_meanGlobal mean of array along its i-axis
[in]maskAn array used for weighting the i-mean
[in]scaleA rescaling factor for the output variable
[in]tmp_scaleA rescaling factor for the internal calculations that is removed from the output

Definition at line 196 of file MOM_spatial_means.F90.

197  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
198  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: array !< The variable being averaged
199  real, dimension(SZJ_(G)), intent(out) :: i_mean !< Global mean of array along its i-axis
200  real, dimension(SZI_(G),SZJ_(G)), &
201  optional, intent(in) :: mask !< An array used for weighting the i-mean
202  real, optional, intent(in) :: scale !< A rescaling factor for the output variable
203  real, optional, intent(in) :: tmp_scale !< A rescaling factor for the internal
204  !! calculations that is removed from the output
205 
206  ! Local variables
207  type(EFP_type), allocatable, dimension(:) :: asum, mask_sum
208  real :: scalefac ! A scaling factor for the variable.
209  real :: unscale ! A factor for undoing any internal rescaling before output.
210  real :: mask_sum_r
211  integer :: is, ie, js, je, idg_off, jdg_off
212  integer :: i, j
213 
214  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
215  idg_off = g%idg_offset ; jdg_off = g%jdg_offset
216 
217  scalefac = 1.0 ; if (present(scale)) scalefac = scale
218  unscale = 1.0
219  if (present(tmp_scale)) then ; if (tmp_scale /= 0.0) then
220  scalefac = scalefac * tmp_scale ; unscale = 1.0 / tmp_scale
221  endif ; endif
222  call reset_efp_overflow_error()
223 
224  allocate(asum(g%jsg:g%jeg))
225  if (present(mask)) then
226  allocate(mask_sum(g%jsg:g%jeg))
227 
228  do j=g%jsg,g%jeg
229  asum(j) = real_to_efp(0.0) ; mask_sum(j) = real_to_efp(0.0)
230  enddo
231 
232  do i=is,ie ; do j=js,je
233  asum(j+jdg_off) = asum(j+jdg_off) + real_to_efp(scalefac*array(i,j)*mask(i,j))
234  mask_sum(j+jdg_off) = mask_sum(j+jdg_off) + real_to_efp(mask(i,j))
235  enddo ; enddo
236 
237  if (query_efp_overflow_error()) call mom_error(fatal, &
238  "global_i_mean overflow error occurred before sums across PEs.")
239 
240  call efp_sum_across_pes(asum(g%jsg:g%jeg), g%jeg-g%jsg+1)
241  call efp_sum_across_pes(mask_sum(g%jsg:g%jeg), g%jeg-g%jsg+1)
242 
243  if (query_efp_overflow_error()) call mom_error(fatal, &
244  "global_i_mean overflow error occurred during sums across PEs.")
245 
246  do j=js,je
247  mask_sum_r = efp_to_real(mask_sum(j+jdg_off))
248  if (mask_sum_r == 0.0 ) then ; i_mean(j) = 0.0 ; else
249  i_mean(j) = efp_to_real(asum(j+jdg_off)) / mask_sum_r
250  endif
251  enddo
252 
253  deallocate(mask_sum)
254  else
255  do j=g%jsg,g%jeg ; asum(j) = real_to_efp(0.0) ; enddo
256 
257  do i=is,ie ; do j=js,je
258  asum(j+jdg_off) = asum(j+jdg_off) + real_to_efp(scalefac*array(i,j))
259  enddo ; enddo
260 
261  if (query_efp_overflow_error()) call mom_error(fatal, &
262  "global_i_mean overflow error occurred before sum across PEs.")
263 
264  call efp_sum_across_pes(asum(g%jsg:g%jeg), g%jeg-g%jsg+1)
265 
266  if (query_efp_overflow_error()) call mom_error(fatal, &
267  "global_i_mean overflow error occurred during sum across PEs.")
268 
269  do j=js,je
270  i_mean(j) = efp_to_real(asum(j+jdg_off)) / real(g%ieg-g%isg+1)
271  enddo
272  endif
273 
274  if (unscale /= 1.0) then ; do j=js,je ; i_mean(j) = unscale*i_mean(j) ; enddo ; endif
275 
276  deallocate(asum)
277 

◆ global_j_mean()

subroutine, public mom_spatial_means::global_j_mean ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  array,
real, dimension( g %isd: g %ied), intent(out)  j_mean,
type(ocean_grid_type), intent(inout)  G,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in), optional  mask,
real, intent(in), optional  scale,
real, intent(in), optional  tmp_scale 
)

Determine the global mean of a field along rows of constant j, returning it in a 1-d array using the local indexing. This uses reproducing sums.

Parameters
[in,out]gThe ocean's grid structure
[in]arrayThe variable being averaged
[out]j_meanGlobal mean of array along its j-axis
[in]maskAn array used for weighting the j-mean
[in]scaleA rescaling factor for the output variable
[in]tmp_scaleA rescaling factor for the internal calculations that is removed from the output

Definition at line 282 of file MOM_spatial_means.F90.

283  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
284  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: array !< The variable being averaged
285  real, dimension(SZI_(G)), intent(out) :: j_mean !< Global mean of array along its j-axis
286  real, dimension(SZI_(G),SZJ_(G)), &
287  optional, intent(in) :: mask !< An array used for weighting the j-mean
288  real, optional, intent(in) :: scale !< A rescaling factor for the output variable
289  real, optional, intent(in) :: tmp_scale !< A rescaling factor for the internal
290  !! calculations that is removed from the output
291 
292  ! Local variables
293  type(EFP_type), allocatable, dimension(:) :: asum, mask_sum
294  real :: mask_sum_r
295  real :: scalefac ! A scaling factor for the variable.
296  real :: unscale ! A factor for undoing any internal rescaling before output.
297  integer :: is, ie, js, je, idg_off, jdg_off
298  integer :: i, j
299 
300  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
301  idg_off = g%idg_offset ; jdg_off = g%jdg_offset
302 
303  scalefac = 1.0 ; if (present(scale)) scalefac = scale
304  unscale = 1.0
305  if (present(tmp_scale)) then ; if (tmp_scale /= 0.0) then
306  scalefac = scalefac * tmp_scale ; unscale = 1.0 / tmp_scale
307  endif ; endif
308  call reset_efp_overflow_error()
309 
310  allocate(asum(g%isg:g%ieg))
311  if (present(mask)) then
312  allocate (mask_sum(g%isg:g%ieg))
313 
314  do i=g%isg,g%ieg
315  asum(i) = real_to_efp(0.0) ; mask_sum(i) = real_to_efp(0.0)
316  enddo
317 
318  do i=is,ie ; do j=js,je
319  asum(i+idg_off) = asum(i+idg_off) + real_to_efp(scalefac*array(i,j)*mask(i,j))
320  mask_sum(i+idg_off) = mask_sum(i+idg_off) + real_to_efp(mask(i,j))
321  enddo ; enddo
322 
323  if (query_efp_overflow_error()) call mom_error(fatal, &
324  "global_j_mean overflow error occurred before sums across PEs.")
325 
326  call efp_sum_across_pes(asum(g%isg:g%ieg), g%ieg-g%isg+1)
327  call efp_sum_across_pes(mask_sum(g%isg:g%ieg), g%ieg-g%isg+1)
328 
329  if (query_efp_overflow_error()) call mom_error(fatal, &
330  "global_j_mean overflow error occurred during sums across PEs.")
331 
332  do i=is,ie
333  mask_sum_r = efp_to_real(mask_sum(i+idg_off))
334  if (mask_sum_r == 0.0 ) then ; j_mean(i) = 0.0 ; else
335  j_mean(i) = efp_to_real(asum(i+idg_off)) / mask_sum_r
336  endif
337  enddo
338 
339  deallocate(mask_sum)
340  else
341  do i=g%isg,g%ieg ; asum(i) = real_to_efp(0.0) ; enddo
342 
343  do i=is,ie ; do j=js,je
344  asum(i+idg_off) = asum(i+idg_off) + real_to_efp(scalefac*array(i,j))
345  enddo ; enddo
346 
347  if (query_efp_overflow_error()) call mom_error(fatal, &
348  "global_j_mean overflow error occurred before sum across PEs.")
349 
350  call efp_sum_across_pes(asum(g%isg:g%ieg), g%ieg-g%isg+1)
351 
352  if (query_efp_overflow_error()) call mom_error(fatal, &
353  "global_j_mean overflow error occurred during sum across PEs.")
354 
355  do i=is,ie
356  j_mean(i) = efp_to_real(asum(i+idg_off)) / real(g%jeg-g%jsg+1)
357  enddo
358  endif
359 
360  if (unscale /= 1.0) then ; do i=is,ie ; j_mean(i) = unscale*j_mean(i) ; enddo ; endif
361 
362  deallocate(asum)
363 

◆ global_layer_mean()

real function, dimension( gv %ke), public mom_spatial_means::global_layer_mean ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  var,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, intent(in), optional  scale 
)

Return the layerwise global thickness-weighted mean of a variable. This uses reproducing sums.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]varThe variable to average
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]scaleA rescaling factor for the variable

Definition at line 83 of file MOM_spatial_means.F90.

84  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
85  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
86  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: var !< The variable to average
87  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
88  real, optional, intent(in) :: scale !< A rescaling factor for the variable
89  real, dimension(SZK_(GV)) :: global_layer_mean
90 
91  real, dimension(G%isc:G%iec, G%jsc:G%jec, SZK_(GV)) :: tmpForSumming, weight
92  type(EFP_type), dimension(2*SZK_(GV)) :: laysums
93  real, dimension(SZK_(GV)) :: global_temp_scalar, global_weight_scalar
94  real :: scalefac ! A scaling factor for the variable.
95  integer :: i, j, k, is, ie, js, je, nz
96  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
97 
98  scalefac = 1.0 ; if (present(scale)) scalefac = scale
99  tmpforsumming(:,:,:) = 0. ; weight(:,:,:) = 0.
100 
101  do k=1,nz ; do j=js,je ; do i=is,ie
102  weight(i,j,k) = (gv%H_to_m * h(i,j,k)) * (g%US%L_to_m**2*g%areaT(i,j) * g%mask2dT(i,j))
103  tmpforsumming(i,j,k) = scalefac * var(i,j,k) * weight(i,j,k)
104  enddo ; enddo ; enddo
105 
106  global_temp_scalar = reproducing_sum(tmpforsumming, efp_lay_sums=laysums(1:nz), only_on_pe=.true.)
107  global_weight_scalar = reproducing_sum(weight, efp_lay_sums=laysums(nz+1:2*nz), only_on_pe=.true.)
108  call efp_sum_across_pes(laysums, 2*nz)
109 
110  do k=1,nz
111  global_layer_mean(k) = efp_to_real(laysums(k)) / efp_to_real(laysums(nz+k))
112  enddo
113 

◆ global_mass_integral()

real function, public mom_spatial_means::global_mass_integral ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in), optional  var,
logical, intent(in), optional  on_PE_only,
real, intent(in), optional  scale 
)

Find the global mass-weighted integral of a variable. This uses reproducing sums.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]varThe variable being integrated
[in]on_pe_onlyIf present and true, the sum is only done on the local PE, and it is not order invariant.
[in]scaleA rescaling factor for the variable
Returns
The mass-weighted integral of var (or 1) in kg times the units of var

Definition at line 148 of file MOM_spatial_means.F90.

149  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
150  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
151  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
152  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
153  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
154  optional, intent(in) :: var !< The variable being integrated
155  logical, optional, intent(in) :: on_PE_only !< If present and true, the sum is only
156  !! done on the local PE, and it is _not_ order invariant.
157  real, optional, intent(in) :: scale !< A rescaling factor for the variable
158  real :: global_mass_integral !< The mass-weighted integral of var (or 1) in
159  !! kg times the units of var
160 
161  real, dimension(SZI_(G), SZJ_(G)) :: tmpForSumming
162  real :: scalefac ! An overall scaling factor for the areas and variable.
163  logical :: global_sum
164  integer :: i, j, k, is, ie, js, je, nz
165  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
166 
167  scalefac = g%US%L_to_m**2 ; if (present(scale)) scalefac = g%US%L_to_m**2*scale
168  tmpforsumming(:,:) = 0.0
169 
170  if (present(var)) then
171  do k=1,nz ; do j=js,je ; do i=is,ie
172  tmpforsumming(i,j) = tmpforsumming(i,j) + var(i,j,k) * &
173  ((gv%H_to_kg_m2 * h(i,j,k)) * (scalefac*g%areaT(i,j) * g%mask2dT(i,j)))
174  enddo ; enddo ; enddo
175  else
176  do k=1,nz ; do j=js,je ; do i=is,ie
177  tmpforsumming(i,j) = tmpforsumming(i,j) + &
178  ((gv%H_to_kg_m2 * h(i,j,k)) * (scalefac*g%areaT(i,j) * g%mask2dT(i,j)))
179  enddo ; enddo ; enddo
180  endif
181  global_sum = .true. ; if (present(on_pe_only)) global_sum = .not.on_pe_only
182  if (global_sum) then
183  global_mass_integral = reproducing_sum(tmpforsumming)
184  else
185  global_mass_integral = 0.0
186  do j=js,je ; do i=is,ie
187  global_mass_integral = global_mass_integral + tmpforsumming(i,j)
188  enddo ; enddo
189  endif
190 

◆ global_volume_mean()

real function, public mom_spatial_means::global_volume_mean ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  var,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, intent(in), optional  scale 
)

Find the global thickness-weighted mean of a variable. This uses reproducing sums.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]varThe variable being averaged
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]scaleA rescaling factor for the variable
Returns
The thickness-weighted average of var

Definition at line 117 of file MOM_spatial_means.F90.

118  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
119  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
120  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
121  intent(in) :: var !< The variable being averaged
122  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
123  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
124  real, optional, intent(in) :: scale !< A rescaling factor for the variable
125  real :: global_volume_mean !< The thickness-weighted average of var
126 
127  real :: scalefac ! A scaling factor for the variable.
128  real :: weight_here
129  real, dimension(SZI_(G), SZJ_(G)) :: tmpForSumming, sum_weight
130  integer :: i, j, k, is, ie, js, je, nz
131  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
132 
133  scalefac = 1.0 ; if (present(scale)) scalefac = scale
134  tmpforsumming(:,:) = 0. ; sum_weight(:,:) = 0.
135 
136  do k=1,nz ; do j=js,je ; do i=is,ie
137  weight_here = (gv%H_to_m * h(i,j,k)) * (g%US%L_to_m**2*g%areaT(i,j) * g%mask2dT(i,j))
138  tmpforsumming(i,j) = tmpforsumming(i,j) + scalefac * var(i,j,k) * weight_here
139  sum_weight(i,j) = sum_weight(i,j) + weight_here
140  enddo ; enddo ; enddo
141  global_volume_mean = (reproducing_sum(tmpforsumming)) / &
142  (reproducing_sum(sum_weight))
143 
mom_coms::reproducing_sum_efp
Find an accurate and order-invariant sum of a distributed 2d field, returning the result in the form ...
Definition: MOM_coms.F90:59
mom_coms::reproducing_sum
Find an accurate and order-invariant sum of a distributed 2d or 3d field.
Definition: MOM_coms.F90:53
mom_coms::efp_sum_across_pes
Sum a value or 1-d array of values across processors, returning the sums in place.
Definition: MOM_coms.F90:64