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)
 Return the global area mean of a variable. This uses reproducing sums. More...
 
real function, public global_area_integral (var, G)
 Return the global area integral of a variable. This uses reproducing sums. More...
 
real function, dimension(gv %ke), public global_layer_mean (var, h, G, GV)
 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)
 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)
 Find the global mass-weighted integral of a variable. This uses reproducing sums. More...
 
subroutine, public global_i_mean (array, i_mean, G, mask)
 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)
 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)
 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 
)

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

Definition at line 315 of file MOM_spatial_means.F90.

315  type(ocean_grid_type), intent(in) :: G !< Grid structure
316  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: array !< 2D array to be adjusted
317  real, optional, intent(out) :: scaling !< The scaling factor used
318  ! Local variables
319  real, dimension(SZI_(G), SZJ_(G)) :: posVals, negVals, areaXposVals, areaXnegVals
320  integer :: i,j
321  real :: areaIntPosVals, areaIntNegVals, posScale, negScale
322 
323  areaxposvals(:,:) = 0.
324  areaxnegvals(:,:) = 0.
325 
326  do j=g%jsc,g%jec ; do i=g%isc,g%iec
327  posvals(i,j) = max(0., array(i,j))
328  areaxposvals(i,j) = g%areaT(i,j) * posvals(i,j)
329  negvals(i,j) = min(0., array(i,j))
330  areaxnegvals(i,j) = g%areaT(i,j) * negvals(i,j)
331  enddo ; enddo
332 
333  areaintposvals = reproducing_sum( areaxposvals )
334  areaintnegvals = reproducing_sum( areaxnegvals )
335 
336  posscale = 0.0 ; negscale = 0.0
337  if ((areaintposvals>0.).and.(areaintnegvals<0.)) then ! Only adjust if possible
338  if (areaintposvals>-areaintnegvals) then ! Scale down positive values
339  posscale = - areaintnegvals / areaintposvals
340  do j=g%jsc,g%jec ; do i=g%isc,g%iec
341  array(i,j) = (posscale * posvals(i,j)) + negvals(i,j)
342  enddo ; enddo
343  elseif (areaintposvals<-areaintnegvals) then ! Scale down negative values
344  negscale = - areaintposvals / areaintnegvals
345  do j=g%jsc,g%jec ; do i=g%isc,g%iec
346  array(i,j) = posvals(i,j) + (negscale * negvals(i,j))
347  enddo ; enddo
348  endif
349  endif
350  if (present(scaling)) scaling = posscale - negscale
351 

◆ 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 
)

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

Parameters
[in]gThe ocean's grid structure
[in]varThe variable to integrate

Definition at line 47 of file MOM_spatial_means.F90.

47  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
48  real, dimension(SZI_(G), SZJ_(G)), intent(in) :: var !< The variable to integrate
49  real, dimension(SZI_(G), SZJ_(G)) :: tmpForSumming
50  real :: global_area_integral
51 
52  integer :: i, j, is, ie, js, je
53  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
54 
55  tmpforsumming(:,:) = 0.
56  do j=js,je ; do i=is, ie
57  tmpforsumming(i,j) = ( var(i,j) * (g%areaT(i,j) * g%mask2dT(i,j)) )
58  enddo ; enddo
59  global_area_integral = reproducing_sum( tmpforsumming )
60 

◆ 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 
)

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

Parameters
[in]gThe ocean's grid structure
[in]varThe variable to average

Definition at line 29 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, dimension(SZI_(G), SZJ_(G)) :: tmpForSumming
32  real :: global_area_mean
33 
34  integer :: i, j, is, ie, js, je
35  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
36 
37  tmpforsumming(:,:) = 0.
38  do j=js,je ; do i=is, ie
39  tmpforsumming(i,j) = ( var(i,j) * (g%areaT(i,j) * g%mask2dT(i,j)) )
40  enddo ; enddo
41  global_area_mean = reproducing_sum( tmpforsumming ) * g%IareaT_global
42 

◆ 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 
)

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

Definition at line 168 of file MOM_spatial_means.F90.

168  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
169  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: array !< The variable being averaged
170  real, dimension(SZJ_(G)), intent(out) :: i_mean !< Global mean of array along its i-axis
171  real, dimension(SZI_(G),SZJ_(G)), &
172  optional, intent(in) :: mask !< An array used for weighting the i-mean
173 
174  ! Local variables
175  type(EFP_type), allocatable, dimension(:) :: asum, mask_sum
176  real :: mask_sum_r
177  integer :: is, ie, js, je, idg_off, jdg_off
178  integer :: i, j
179 
180  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
181  idg_off = g%idg_offset ; jdg_off = g%jdg_offset
182 
183  call reset_efp_overflow_error()
184 
185  allocate(asum(g%jsg:g%jeg))
186  if (present(mask)) then
187  allocate(mask_sum(g%jsg:g%jeg))
188 
189  do j=g%jsg,g%jeg
190  asum(j) = real_to_efp(0.0) ; mask_sum(j) = real_to_efp(0.0)
191  enddo
192 
193  do i=is,ie ; do j=js,je
194  asum(j+jdg_off) = asum(j+jdg_off) + real_to_efp(array(i,j)*mask(i,j))
195  mask_sum(j+jdg_off) = mask_sum(j+jdg_off) + real_to_efp(mask(i,j))
196  enddo ; enddo
197 
198  if (query_efp_overflow_error()) call mom_error(fatal, &
199  "global_i_mean overflow error occurred before sums across PEs.")
200 
201  call efp_list_sum_across_pes(asum(g%jsg:g%jeg), g%jeg-g%jsg+1)
202  call efp_list_sum_across_pes(mask_sum(g%jsg:g%jeg), g%jeg-g%jsg+1)
203 
204  if (query_efp_overflow_error()) call mom_error(fatal, &
205  "global_i_mean overflow error occurred during sums across PEs.")
206 
207  do j=js,je
208  mask_sum_r = efp_to_real(mask_sum(j+jdg_off))
209  if (mask_sum_r == 0.0 ) then ; i_mean(j) = 0.0 ; else
210  i_mean(j) = efp_to_real(asum(j+jdg_off)) / mask_sum_r
211  endif
212  enddo
213 
214  deallocate(mask_sum)
215  else
216  do j=g%jsg,g%jeg ; asum(j) = real_to_efp(0.0) ; enddo
217 
218  do i=is,ie ; do j=js,je
219  asum(j+jdg_off) = asum(j+jdg_off) + real_to_efp(array(i,j))
220  enddo ; enddo
221 
222  if (query_efp_overflow_error()) call mom_error(fatal, &
223  "global_i_mean overflow error occurred before sum across PEs.")
224 
225  call efp_list_sum_across_pes(asum(g%jsg:g%jeg), g%jeg-g%jsg+1)
226 
227  if (query_efp_overflow_error()) call mom_error(fatal, &
228  "global_i_mean overflow error occurred during sum across PEs.")
229 
230  do j=js,je
231  i_mean(j) = efp_to_real(asum(j+jdg_off)) / real(g%ieg-g%isg+1)
232  enddo
233  endif
234 
235  deallocate(asum)
236 

◆ 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 
)

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

Definition at line 242 of file MOM_spatial_means.F90.

242  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
243  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: array !< The variable being averaged
244  real, dimension(SZI_(G)), intent(out) :: j_mean !< Global mean of array along its j-axis
245  real, dimension(SZI_(G),SZJ_(G)), &
246  optional, intent(in) :: mask !< An array used for weighting the j-mean
247 
248  ! Local variables
249  type(EFP_type), allocatable, dimension(:) :: asum, mask_sum
250  real :: mask_sum_r
251  integer :: is, ie, js, je, idg_off, jdg_off
252  integer :: i, j
253 
254  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
255  idg_off = g%idg_offset ; jdg_off = g%jdg_offset
256 
257  call reset_efp_overflow_error()
258 
259  allocate(asum(g%isg:g%ieg))
260  if (present(mask)) then
261  allocate (mask_sum(g%isg:g%ieg))
262 
263  do i=g%isg,g%ieg
264  asum(i) = real_to_efp(0.0) ; mask_sum(i) = real_to_efp(0.0)
265  enddo
266 
267  do i=is,ie ; do j=js,je
268  asum(i+idg_off) = asum(i+idg_off) + real_to_efp(array(i,j)*mask(i,j))
269  mask_sum(i+idg_off) = mask_sum(i+idg_off) + real_to_efp(mask(i,j))
270  enddo ; enddo
271 
272  if (query_efp_overflow_error()) call mom_error(fatal, &
273  "global_j_mean overflow error occurred before sums across PEs.")
274 
275  call efp_list_sum_across_pes(asum(g%isg:g%ieg), g%ieg-g%isg+1)
276  call efp_list_sum_across_pes(mask_sum(g%isg:g%ieg), g%ieg-g%isg+1)
277 
278  if (query_efp_overflow_error()) call mom_error(fatal, &
279  "global_j_mean overflow error occurred during sums across PEs.")
280 
281  do i=is,ie
282  mask_sum_r = efp_to_real(mask_sum(i+idg_off))
283  if (mask_sum_r == 0.0 ) then ; j_mean(i) = 0.0 ; else
284  j_mean(i) = efp_to_real(asum(i+idg_off)) / mask_sum_r
285  endif
286  enddo
287 
288  deallocate(mask_sum)
289  else
290  do i=g%isg,g%ieg ; asum(i) = real_to_efp(0.0) ; enddo
291 
292  do i=is,ie ; do j=js,je
293  asum(i+idg_off) = asum(i+idg_off) + real_to_efp(array(i,j))
294  enddo ; enddo
295 
296  if (query_efp_overflow_error()) call mom_error(fatal, &
297  "global_j_mean overflow error occurred before sum across PEs.")
298 
299  call efp_list_sum_across_pes(asum(g%isg:g%ieg), g%ieg-g%isg+1)
300 
301  if (query_efp_overflow_error()) call mom_error(fatal, &
302  "global_j_mean overflow error occurred during sum across PEs.")
303 
304  do i=is,ie
305  j_mean(i) = efp_to_real(asum(i+idg_off)) / real(g%jeg-g%jsg+1)
306  enddo
307  endif
308 
309  deallocate(asum)
310 

◆ 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 
)

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]

Definition at line 65 of file MOM_spatial_means.F90.

65  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
66  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
67  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: var !< The variable to average
68  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
69  real, dimension(SZK_(GV)) :: global_layer_mean
70 
71  real, dimension(SZI_(G), SZJ_(G), SZK_(GV)) :: tmpForSumming, weight
72  real, dimension(SZK_(GV)) :: scalarij, weightij
73  real, dimension(SZK_(GV)) :: global_temp_scalar, global_weight_scalar
74  integer :: i, j, k, is, ie, js, je, nz
75  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
76 
77  tmpforsumming(:,:,:) = 0. ; weight(:,:,:) = 0.
78 
79  do k=1,nz ; do j=js,je ; do i=is,ie
80  weight(i,j,k) = (gv%H_to_m * h(i,j,k)) * (g%areaT(i,j) * g%mask2dT(i,j))
81  tmpforsumming(i,j,k) = var(i,j,k) * weight(i,j,k)
82  enddo ; enddo ; enddo
83 
84  global_temp_scalar = reproducing_sum(tmpforsumming,sums=scalarij)
85  global_weight_scalar = reproducing_sum(weight,sums=weightij)
86 
87  do k=1, nz
88  global_layer_mean(k) = scalarij(k) / weightij(k)
89  enddo
90 

◆ 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 
)

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.
Returns
The mass-weighted integral of var (or 1) in kg times the units of var

Definition at line 123 of file MOM_spatial_means.F90.

123  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
124  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
125  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
126  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
127  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
128  optional, intent(in) :: var !< The variable being integrated
129  logical, optional, intent(in) :: on_PE_only !< If present and true, the sum is only
130  !! done on the local PE, and it is _not_ order invariant.
131  real :: global_mass_integral !< The mass-weighted integral of var (or 1) in
132  !! kg times the units of var
133 
134  real, dimension(SZI_(G), SZJ_(G)) :: tmpForSumming
135  logical :: global_sum
136  integer :: i, j, k, is, ie, js, je, nz
137  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
138 
139  tmpforsumming(:,:) = 0.0
140 
141  if (present(var)) then
142  do k=1,nz ; do j=js,je ; do i=is,ie
143  tmpforsumming(i,j) = tmpforsumming(i,j) + var(i,j,k) * &
144  ((gv%H_to_kg_m2 * h(i,j,k)) * (g%areaT(i,j) * g%mask2dT(i,j)))
145  enddo ; enddo ; enddo
146  else
147  do k=1,nz ; do j=js,je ; do i=is,ie
148  tmpforsumming(i,j) = tmpforsumming(i,j) + &
149  ((gv%H_to_kg_m2 * h(i,j,k)) * (g%areaT(i,j) * g%mask2dT(i,j)))
150  enddo ; enddo ; enddo
151  endif
152  global_sum = .true. ; if (present(on_pe_only)) global_sum = .not.on_pe_only
153  if (global_sum) then
154  global_mass_integral = reproducing_sum(tmpforsumming)
155  else
156  global_mass_integral = 0.0
157  do j=js,je ; do i=is,ie
158  global_mass_integral = global_mass_integral + tmpforsumming(i,j)
159  enddo ; enddo
160  endif
161 

◆ 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 
)

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]
Returns
The thickness-weighted average of var

Definition at line 95 of file MOM_spatial_means.F90.

95  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
96  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
97  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
98  intent(in) :: var !< The variable being averaged
99  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
100  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
101  real :: global_volume_mean !< The thickness-weighted average of var
102 
103  real :: weight_here
104  real, dimension(SZI_(G), SZJ_(G)) :: tmpForSumming, sum_weight
105  integer :: i, j, k, is, ie, js, je, nz
106  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
107 
108  tmpforsumming(:,:) = 0. ; sum_weight(:,:) = 0.
109 
110  do k=1,nz ; do j=js,je ; do i=is,ie
111  weight_here = (gv%H_to_m * h(i,j,k)) * (g%areaT(i,j) * g%mask2dT(i,j))
112  tmpforsumming(i,j) = tmpforsumming(i,j) + var(i,j,k) * weight_here
113  sum_weight(i,j) = sum_weight(i,j) + weight_here
114  enddo ; enddo ; enddo
115  global_volume_mean = (reproducing_sum(tmpforsumming)) / &
116  (reproducing_sum(sum_weight))
117 
mom_coms::reproducing_sum
Find an accurate and order-invariant sum of distributed 2d or 3d fields.
Definition: MOM_coms.F90:54